CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
APVCyclePhaseMonitor Class Reference

#include <DPGAnalysis/SiStripTools/plugins/APVCyclePhaseMonitor.cc>

Inheritance diagram for APVCyclePhaseMonitor:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

 APVCyclePhaseMonitor (const edm::ParameterSet &)
 
 ~APVCyclePhaseMonitor ()
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &) override
 
virtual void beginJob () override
 
virtual void beginRun (const edm::Run &, const edm::EventSetup &) override
 
virtual void endJob () override
 
virtual void endRun (const edm::Run &, const edm::EventSetup &) override
 

Private Attributes

const edm::InputTag _apvphasecollection
 
std::map< std::string, TH1F * > _hphases
 
std::map< std::string, TProfile * > _hphasevsorbit
 
std::map< std::string, TH1F ** > _hselectedphases
 
std::map< std::string, TH1F ** > _hselectedphasessize
 
std::map< std::string, TH1F ** > _hselectedphasesvector
 
std::map< std::string,
TProfile ** > 
_hselectedphasevectorvsorbit
 
std::map< std::string,
TProfile ** > 
_hselectedphasevsorbit
 
unsigned int _nevents
 
std::vector< std::string > _selectedparts
 
std::vector< std::string > _selectedvectorparts
 
const unsigned int m_LSfrac
 
const unsigned int m_maxLS
 
RunHistogramManager m_rhm
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 51 of file APVCyclePhaseMonitor.cc.

Constructor & Destructor Documentation

APVCyclePhaseMonitor::APVCyclePhaseMonitor ( const edm::ParameterSet iConfig)
explicit

Definition at line 93 of file APVCyclePhaseMonitor.cc.

References _apvphasecollection, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _selectedparts, _selectedvectorparts, m_LSfrac, m_maxLS, m_rhm, RunHistogramManager::makeTH1F(), and RunHistogramManager::makeTProfile().

93  :
94  _apvphasecollection(iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection")),
95  _selectedparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedPartitions",std::vector<std::string>())),
96  _selectedvectorparts(iConfig.getUntrackedParameter<std::vector<std::string> >("selectedVectorPartitions",std::vector<std::string>())),
97  m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin",125)),
98  m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction",16)),
101  _nevents(0)
102 {
103  //now do what ever initialization is needed
104 
105  edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << _apvphasecollection << " used";
106 
107  for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
108 
109  char hname[300];
110 
111  sprintf(hname,"selected_phase_%s",part->c_str());
112  edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
113  _hselectedphases[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
114 
115  sprintf(hname,"selected_phasevsorbit_%s",part->c_str());
116  edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
117  _hselectedphasevsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
118  }
119 
120  for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
121  part!=_selectedvectorparts.end();++part) {
122 
123  char hname[300];
124 
125  sprintf(hname,"selected_phase_vector_%s",part->c_str());
126  edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
127  _hselectedphasesvector[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
128 
129  sprintf(hname,"selected_phase_size_%s",part->c_str());
130  edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
131  _hselectedphasessize[*part] = m_rhm.makeTH1F(hname,hname,10,-0.5,9.5);
132 
133  sprintf(hname,"selected_phasevectorvsorbit_%s",part->c_str());
134  edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
135  _hselectedphasevectorvsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
136  }
137 
138 
139 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::map< std::string, TH1F ** > _hselectedphasesvector
std::vector< std::string > _selectedparts
std::map< std::string, TProfile ** > _hselectedphasevectorvsorbit
std::map< std::string, TProfile ** > _hselectedphasevsorbit
std::map< std::string, TH1F * > _hphases
TH1F ** makeTH1F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
TProfile ** makeTProfile(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
std::map< std::string, TH1F ** > _hselectedphasessize
const unsigned int m_maxLS
std::vector< std::string > _selectedvectorparts
part
Definition: HCALResponse.h:20
std::map< std::string, TH1F ** > _hselectedphases
const edm::InputTag _apvphasecollection
std::map< std::string, TProfile * > _hphasevsorbit
RunHistogramManager m_rhm
const unsigned int m_LSfrac
APVCyclePhaseMonitor::~APVCyclePhaseMonitor ( )

Definition at line 142 of file APVCyclePhaseMonitor.cc.

143 {
144 
145  // do anything here that needs to be done at desctruction time
146  // (e.g. close files, deallocate resources etc.)
147 
148 }

Member Function Documentation

void APVCyclePhaseMonitor::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overrideprivatevirtual

Implements edm::EDAnalyzer.

Definition at line 157 of file APVCyclePhaseMonitor.cc.

References _apvphasecollection, _hphases, _hphasevsorbit, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _nevents, _selectedparts, _selectedvectorparts, HcalObjRepresent::Fill(), edm::Event::getByLabel(), m_LSfrac, m_maxLS, TFileDirectory::mkdir(), edm::EventBase::orbitNumber(), and edm::Event::run().

158 {
159  using namespace edm;
160 
161  _nevents++;
162 
164  iEvent.getByLabel(_apvphasecollection,apvphases);
165 
166  // improve the matchin between default and actual partitions
167 
169 
170  for(std::map<std::string,int>::const_iterator phase = apvphases->get().begin(); phase != apvphases->get().end(); ++phase) {
171 
172  if(_hphases.find(phase->first)==_hphases.end()) {
173  char dirname[300];
174  sprintf(dirname,"run_%d",iEvent.run());
175  TFileDirectory subrun = tfserv->mkdir(dirname);
176 
177  char hname[300];
178 
179  sprintf(hname,"phase_%s",phase->first.c_str());
180  edm::LogInfo("TH1FBeingBooked") << "TH1F " << hname << " being booked" ;
181  _hphases[phase->first] = subrun.make<TH1F>(hname,hname,70,-0.5,69.5);
182  _hphases[phase->first]->GetXaxis()->SetTitle("BX mod 70"); _hphases[phase->first]->GetYaxis()->SetTitle("Events");
183 
184  sprintf(hname,"phasevsorbit_%s",phase->first.c_str());
185  edm::LogInfo("TProfileBeingBooked") << "TProfile " << hname << " being booked" ;
186  _hphasevsorbit[phase->first] = subrun.make<TProfile>(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
187  _hphasevsorbit[phase->first]->SetBit(TH1::kCanRebin);
188  _hphasevsorbit[phase->first]->GetXaxis()->SetTitle("time [orbit#]"); _hphasevsorbit[phase->first]->GetYaxis()->SetTitle("Phase");
189 
190  }
191  _hphases[phase->first]->Fill(phase->second);
192  _hphasevsorbit[phase->first]->Fill(iEvent.orbitNumber(),phase->second);
193  }
194 
195  // selected partitions
196 
197  for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
199  (*_hselectedphases[*part])->Fill(apvphases->getPhase(*part));
200  }
202  (*_hselectedphasevsorbit[*part])->Fill(iEvent.orbitNumber(),apvphases->getPhase(*part));
203  }
204  }
205 
206  for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
207  part!=_selectedvectorparts.end();++part) {
208 
209  const std::vector<int> phases = apvphases->getPhases(*part);
210 
212  (*_hselectedphasessize[*part])->Fill(phases.size());
213  }
214 
215  for(std::vector<int>::const_iterator phase=phases.begin();phase!=phases.end();++phase) {
217  (*_hselectedphasesvector[*part])->Fill(*phase);
218  }
221  (*_hselectedphasevectorvsorbit[*part])->Fill(iEvent.orbitNumber(),*phase);
222  }
223  }
224 
225  }
226 }
std::map< std::string, TH1F ** > _hselectedphasesvector
std::vector< std::string > _selectedparts
std::map< std::string, TProfile ** > _hselectedphasevectorvsorbit
std::map< std::string, TProfile ** > _hselectedphasevsorbit
std::map< std::string, TH1F * > _hphases
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
RunNumber_t run() const
Definition: Event.h:88
std::map< std::string, TH1F ** > _hselectedphasessize
int orbitNumber() const
Definition: EventBase.h:63
const unsigned int m_maxLS
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
std::vector< std::string > _selectedvectorparts
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
part
Definition: HCALResponse.h:20
std::map< std::string, TH1F ** > _hselectedphases
const edm::InputTag _apvphasecollection
std::map< std::string, TProfile * > _hphasevsorbit
const unsigned int m_LSfrac
void APVCyclePhaseMonitor::beginJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 277 of file APVCyclePhaseMonitor.cc.

278 {
279 
280 }
void APVCyclePhaseMonitor::beginRun ( const edm::Run iRun,
const edm::EventSetup  
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 229 of file APVCyclePhaseMonitor.cc.

References _hphases, _hphasevsorbit, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, RunHistogramManager::beginRun(), estimatePileup::hist, and m_rhm.

230 {
231 
232  _hphases.clear();
233  _hphasevsorbit.clear();
234 
235  m_rhm.beginRun(iRun);
236 
237  for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphases.begin();hist!=_hselectedphases.end();++hist) {
238  if(*(hist->second)) {
239  (*(hist->second))->GetXaxis()->SetTitle("BX mod 70"); (*(hist->second))->GetYaxis()->SetTitle("Events");
240  }
241  }
242  for(std::map<std::string,TProfile**>::const_iterator prof=_hselectedphasevsorbit.begin();prof!=_hselectedphasevsorbit.end();++prof) {
243  if(*(prof->second)) {
244  (*(prof->second))->SetBit(TH1::kCanRebin);
245  (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
246  (*(prof->second))->GetYaxis()->SetTitle("Phase");
247  }
248  }
249  for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphasesvector.begin();hist!=_hselectedphasesvector.end();++hist) {
250  if(*(hist->second)) {
251  (*(hist->second))->GetXaxis()->SetTitle("BX mod 70"); (*(hist->second))->GetYaxis()->SetTitle("Events");
252  }
253  }
254  for(std::map<std::string,TH1F**>::const_iterator hist=_hselectedphasessize.begin();hist!=_hselectedphasessize.end();++hist) {
255  if(*(hist->second)) {
256  (*(hist->second))->GetXaxis()->SetTitle("Number of Phases"); (*(hist->second))->GetYaxis()->SetTitle("Events");
257  }
258  }
259  for(std::map<std::string,TProfile**>::const_iterator prof=_hselectedphasevectorvsorbit.begin();prof!=_hselectedphasevectorvsorbit.end();++prof) {
260  if(*(prof->second)) {
261  (*(prof->second))->SetBit(TH1::kCanRebin);
262  (*(prof->second))->GetXaxis()->SetTitle("time [orbit#]");
263  (*(prof->second))->GetYaxis()->SetTitle("Phase");
264  }
265  }
266 
267 }
std::map< std::string, TH1F ** > _hselectedphasesvector
std::map< std::string, TProfile ** > _hselectedphasevectorvsorbit
std::map< std::string, TProfile ** > _hselectedphasevsorbit
std::map< std::string, TH1F * > _hphases
std::map< std::string, TH1F ** > _hselectedphasessize
void beginRun(const edm::Run &iRun)
std::map< std::string, TH1F ** > _hselectedphases
std::map< std::string, TProfile * > _hphasevsorbit
RunHistogramManager m_rhm
void APVCyclePhaseMonitor::endJob ( void  )
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 284 of file APVCyclePhaseMonitor.cc.

References _nevents.

284  {
285 
286  edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
287 
288 }
void APVCyclePhaseMonitor::endRun ( const edm::Run iRun,
const edm::EventSetup  
)
overrideprivatevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 270 of file APVCyclePhaseMonitor.cc.

271 {
272 }

Member Data Documentation

const edm::InputTag APVCyclePhaseMonitor::_apvphasecollection
private

Definition at line 66 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

std::map<std::string,TH1F*> APVCyclePhaseMonitor::_hphases
private

Definition at line 72 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TProfile*> APVCyclePhaseMonitor::_hphasevsorbit
private

Definition at line 76 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and beginRun().

std::map<std::string,TH1F**> APVCyclePhaseMonitor::_hselectedphases
private

Definition at line 73 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), APVCyclePhaseMonitor(), and beginRun().

std::map<std::string,TH1F**> APVCyclePhaseMonitor::_hselectedphasessize
private

Definition at line 75 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), APVCyclePhaseMonitor(), and beginRun().

std::map<std::string,TH1F**> APVCyclePhaseMonitor::_hselectedphasesvector
private

Definition at line 74 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), APVCyclePhaseMonitor(), and beginRun().

std::map<std::string,TProfile**> APVCyclePhaseMonitor::_hselectedphasevectorvsorbit
private

Definition at line 78 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), APVCyclePhaseMonitor(), and beginRun().

std::map<std::string,TProfile**> APVCyclePhaseMonitor::_hselectedphasevsorbit
private

Definition at line 77 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), APVCyclePhaseMonitor(), and beginRun().

unsigned int APVCyclePhaseMonitor::_nevents
private

Definition at line 79 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and endJob().

std::vector<std::string> APVCyclePhaseMonitor::_selectedparts
private

Definition at line 67 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

std::vector<std::string> APVCyclePhaseMonitor::_selectedvectorparts
private

Definition at line 68 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

const unsigned int APVCyclePhaseMonitor::m_LSfrac
private

Definition at line 70 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

const unsigned int APVCyclePhaseMonitor::m_maxLS
private

Definition at line 69 of file APVCyclePhaseMonitor.cc.

Referenced by analyze(), and APVCyclePhaseMonitor().

RunHistogramManager APVCyclePhaseMonitor::m_rhm
private

Definition at line 71 of file APVCyclePhaseMonitor.cc.

Referenced by APVCyclePhaseMonitor(), and beginRun().