CMS 3D CMS Logo

APVCyclePhaseMonitor.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripTools
4 // Class: APVCyclePhaseMonitor
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Tue Jul 19 11:56:00 CEST 2009
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
24 #include "TH1F.h"
25 #include "TProfile.h"
26 
29 
33 
35 
37 
40 
42 
44 
46 
47 //
48 // class decleration
49 //
50 
52  public:
53  explicit APVCyclePhaseMonitor(const edm::ParameterSet&);
54  ~APVCyclePhaseMonitor() override;
55 
56 
57  private:
58  void beginJob() override ;
59  void beginRun(const edm::Run&, const edm::EventSetup&) override;
60  void endRun(const edm::Run&, const edm::EventSetup&) override;
61  void analyze(const edm::Event&, const edm::EventSetup&) override;
62  void endJob() override ;
63 
64  // ----------member data ---------------------------
65 
67  std::vector<std::string> _selectedparts;
68  std::vector<std::string> _selectedvectorparts;
69  const unsigned int m_maxLS;
70  const unsigned int m_LSfrac;
72  std::map<std::string,TH1F*> _hphases;
73  std::map<std::string,TH1F**> _hselectedphases;
74  std::map<std::string,TH1F**> _hselectedphasesvector;
75  std::map<std::string,TH1F**> _hselectedphasessize;
76  std::map<std::string,TProfile*> _hphasevsorbit;
77  std::map<std::string,TProfile**> _hselectedphasevsorbit;
78  std::map<std::string,TProfile**> _hselectedphasevectorvsorbit;
79  unsigned int _nevents;
80 };
81 
82 //
83 // constants, enums and typedefs
84 //
85 
86 //
87 // static data member definitions
88 //
89 
90 //
91 // constructors and destructor
92 //
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 }
141 
142 
144 {
145 
146  // do anything here that needs to be done at desctruction time
147  // (e.g. close files, deallocate resources etc.)
148 
149 }
150 
151 
152 //
153 // member functions
154 //
155 
156 // ------------ method called to for each event ------------
157 void
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 }
228 
229 void
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 }
269 
270 void
272 {
273 }
274 
275 
276 // ------------ method called once each job just before starting event loop ------------
277 void
279 {
280 
281 }
282 
283 // ------------ method called once each job just after ending the event loop ------------
284 void
286 
287  edm::LogInfo("EndOfJob") << _nevents << " analyzed events";
288 
289 }
290 
291 
292 //define this as a plug-in
T getParameter(std::string const &) const
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:517
std::map< std::string, TProfile ** > _hselectedphasevsorbit
std::map< std::string, TH1F * > _hphases
void endRun(const edm::Run &, const edm::EventSetup &) override
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
const std::map< std::string, int > & get() const
TH1F ** makeTH1F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax)
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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)
void beginRun(const edm::Run &, const edm::EventSetup &) override
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)
RunNumber_t run() const
Definition: Event.h:101
APVCyclePhaseMonitor(const edm::ParameterSet &)
std::map< std::string, TH1F ** > _hselectedphasessize
int orbitNumber() const
Definition: EventBase.h:65
void analyze(const edm::Event &, const edm::EventSetup &) override
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
void beginRun(const edm::Run &iRun)
std::map< std::string, TH1F ** > _hselectedphases
HLT enums.
std::map< std::string, TProfile * > _hphasevsorbit
RunHistogramManager m_rhm
const unsigned int m_LSfrac
Definition: Run.h:45