CMS 3D CMS Logo

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 () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

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

Private Attributes

edm::EDGetTokenT< APVCyclePhaseCollection_apvphasecollectionToken
 
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
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- 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 _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _selectedparts, _selectedvectorparts, edm::ParameterSet::getParameter(), m_LSfrac, m_maxLS, m_rhm, RunHistogramManager::makeTH1F(), and RunHistogramManager::makeTProfile().

93  :
94  _apvphasecollectionToken(consumes<APVCyclePhaseCollection>(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)),
102  _nevents(0)
103 {
104  //now do what ever initialization is needed
105 
106  edm::LogInfo("UsedAPVCyclePhaseCollection") << " APVCyclePhaseCollection " << iConfig.getParameter<edm::InputTag>("apvCyclePhaseCollection") << " used";
107 
108  for(std::vector<std::string>::const_iterator part=_selectedparts.begin();part!=_selectedparts.end();++part) {
109 
110  char hname[300];
111 
112  sprintf(hname,"selected_phase_%s",part->c_str());
113  edm::LogInfo("SelectedTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
114  _hselectedphases[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
115 
116  sprintf(hname,"selected_phasevsorbit_%s",part->c_str());
117  edm::LogInfo("SelectedTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
118  _hselectedphasevsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
119  }
120 
121  for(std::vector<std::string>::const_iterator part=_selectedvectorparts.begin();
122  part!=_selectedvectorparts.end();++part) {
123 
124  char hname[300];
125 
126  sprintf(hname,"selected_phase_vector_%s",part->c_str());
127  edm::LogInfo("SelectedVectTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
128  _hselectedphasesvector[*part] = m_rhm.makeTH1F(hname,hname,70,-0.5,69.5);
129 
130  sprintf(hname,"selected_phase_size_%s",part->c_str());
131  edm::LogInfo("SelectedVectSizeTH1FBeingBooked") << "TH1F " << hname << " being booked" ;
132  _hselectedphasessize[*part] = m_rhm.makeTH1F(hname,hname,10,-0.5,9.5);
133 
134  sprintf(hname,"selected_phasevectorvsorbit_%s",part->c_str());
135  edm::LogInfo("SelectedVectTProfileBeingBooked") << "TProfile " << hname << " being booked" ;
136  _hselectedphasevectorvsorbit[*part] = m_rhm.makeTProfile(hname,hname,m_LSfrac*m_maxLS,0,m_maxLS*262144);
137  }
138 
139 
140 }
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
edm::EDGetTokenT< APVCyclePhaseCollection > _apvphasecollectionToken
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)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
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
std::map< std::string, TProfile * > _hphasevsorbit
RunHistogramManager m_rhm
const unsigned int m_LSfrac
APVCyclePhaseMonitor::~APVCyclePhaseMonitor ( )
override

Definition at line 143 of file APVCyclePhaseMonitor.cc.

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

Member Function Documentation

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

Definition at line 158 of file APVCyclePhaseMonitor.cc.

References _apvphasecollectionToken, _hphases, _hphasevsorbit, _hselectedphases, _hselectedphasessize, _hselectedphasesvector, _hselectedphasevectorvsorbit, _hselectedphasevsorbit, _nevents, _selectedparts, _selectedvectorparts, compare_using_db::dirname, HcalObjRepresent::Fill(), APVCyclePhaseCollection::get(), edm::Event::getByToken(), APVCyclePhaseCollection::getPhase(), APVCyclePhaseCollection::getPhases(), m_LSfrac, m_maxLS, TFileService::mkdir(), edm::EventBase::orbitNumber(), hcalRecAlgoESProd_cfi::phase, phases, and edm::Event::run().

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

Reimplemented from edm::EDAnalyzer.

Definition at line 278 of file APVCyclePhaseMonitor.cc.

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

Definition at line 230 of file APVCyclePhaseMonitor.cc.

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

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

References _nevents, and DEFINE_FWK_MODULE.

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

Definition at line 271 of file APVCyclePhaseMonitor.cc.

272 {
273 }

Member Data Documentation

edm::EDGetTokenT<APVCyclePhaseCollection> APVCyclePhaseMonitor::_apvphasecollectionToken
private

Definition at line 66 of file APVCyclePhaseMonitor.cc.

Referenced by analyze().

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().