CMS 3D CMS Logo

EventTimeDistribution.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: SiStripTools
4 // Class: EventTimeDistribution
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Tue Jul 19 11:56:00 CEST 2009
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
23 #include <string>
24 #include <vector>
25 
26 #include "TH1F.h"
27 #include "TH2F.h"
28 
31 
35 
37 
39 
42 
44 
47 
49 //
50 // class decleration
51 //
52 
54 public:
56  ~EventTimeDistribution() override;
57 
58 private:
59  void beginJob() override;
60  void beginRun(const edm::Run&, const edm::EventSetup&) override;
61  void endRun(const edm::Run&, const edm::EventSetup&) override;
62  void analyze(const edm::Event&, const edm::EventSetup&) override;
63  void endJob() override;
64 
65  // ----------member data ---------------------------
66 
70  const bool _wantdbxvsbxincycle;
71  const bool _wantdbxvsbx;
72  const bool _wantbxincyclevsbx;
74  unsigned int _nevents;
75  const unsigned int m_maxLS;
76  const unsigned int m_LSfrac;
77  const bool m_ewhdepthHisto;
78 
80 
81  TH1F** _dbx;
82  std::vector<TH1F**> m_dbxhistos;
83  std::vector<std::pair<unsigned int, unsigned int> > m_dbxindices;
84  TH1F** _bx;
85  TH1F** _bxincycle;
86  TH1F** _orbit;
88  TH2F** _dbxvsbx;
91  TH1F** m_ewhdepth;
92 };
93 
94 //
95 // constants, enums and typedefs
96 //
97 
98 //
99 // static data member definitions
100 //
101 
102 //
103 // constructors and destructor
104 //
106  : _historyProductToken(consumes<EventWithHistory>(iConfig.getParameter<edm::InputTag>("historyProduct"))),
107  _apvphasecollToken(consumes<APVCyclePhaseCollection>(iConfig.getParameter<edm::InputTag>("apvPhaseCollection"))),
108  _phasepart(iConfig.getUntrackedParameter<std::string>("phasePartition", "None")),
109  _wantdbxvsbxincycle(iConfig.getUntrackedParameter<bool>("wantDBXvsBXincycle", false)),
110  _wantdbxvsbx(iConfig.getUntrackedParameter<bool>("wantDBXvsBX", false)),
111  _wantbxincyclevsbx(iConfig.getUntrackedParameter<bool>("wantBXincyclevsBX", false)),
112  _wantorbitvsbxincycle(iConfig.getUntrackedParameter<bool>("wantOrbitvsBXincycle", false)),
113  _nevents(0),
114  m_maxLS(iConfig.getUntrackedParameter<unsigned int>("maxLSBeforeRebin", 100)),
115  m_LSfrac(iConfig.getUntrackedParameter<unsigned int>("startingLSFraction", 4)),
116  m_ewhdepthHisto(iConfig.getUntrackedParameter<bool>("wantEWHDepthHisto", false)),
119  _dbxvsbx(nullptr),
123  //now do what ever initialization is needed
124 
125  std::vector<edm::ParameterSet> dbxhistoparams(iConfig.getUntrackedParameter<std::vector<edm::ParameterSet> >(
126  "dbxHistosParams", std::vector<edm::ParameterSet>()));
127 
128  for (std::vector<edm::ParameterSet>::const_iterator params = dbxhistoparams.begin(); params != dbxhistoparams.end();
129  ++params) {
130  m_dbxindices.push_back(std::pair<unsigned int, unsigned int>(params->getParameter<unsigned int>("firstEvent"),
131  params->getParameter<unsigned int>("secondEvent")));
132  char hname[300];
133  sprintf(hname,
134  "dbx_%d_%d",
135  params->getParameter<unsigned int>("firstEvent"),
136  params->getParameter<unsigned int>("secondEvent"));
137  char htitle[300];
138  sprintf(htitle,
139  "dbx(%d,%d)",
140  params->getParameter<unsigned int>("firstEvent"),
141  params->getParameter<unsigned int>("secondEvent"));
142 
143  m_dbxhistos.push_back(_rhm.makeTH1F(hname,
144  htitle,
145  params->getParameter<int>("nbins"),
146  params->getParameter<double>("min"),
147  params->getParameter<double>("max")));
148  LogDebug("DBXHistoPreBooking") << "Booked DBX histo named " << hname << " untitled " << htitle;
149  }
150 
151  _dbx = _rhm.makeTH1F("dbx", "dbx", 1000, -0.5, 999.5);
152  _bx = _rhm.makeTH1F("bx", "BX number", 3564, -0.5, 3563.5);
153  _bxincycle = _rhm.makeTH1F("bxcycle", "bxcycle", 70, -0.5, 69.5);
154  _orbit = _rhm.makeTH1F("orbit", "orbit", m_LSfrac * m_maxLS, 0, m_maxLS * 262144);
156  _dbxvsbxincycle = _rhm.makeTH2F("dbxvsbxincycle", "dbxvsbxincycle", 70, -0.5, 69.5, 1000, -0.5, 999.5);
157  if (_wantdbxvsbx)
158  _dbxvsbx = _rhm.makeTH2F("dbxvsbx", "dbxvsbx", 3564, -0.5, 3563.5, 1000, -0.5, 999.5);
159  if (_wantbxincyclevsbx)
160  _bxincyclevsbx = _rhm.makeTH2F("bxincyclevsbx", "bxincyclevsbx", 3564, -0.5, 3563.5, 70, -0.5, 69.5);
163  _rhm.makeTH2F("orbitvsbxincycle", "orbitvsbxincycle", 70, -0.5, 69.5, m_maxLS, 0, m_maxLS * 262144);
164  if (m_ewhdepthHisto)
165  m_ewhdepth = _rhm.makeTH1F("ewhdepth", "EventWithHistory Depth", 11, -0.5, 10.5);
166 
167  edm::LogInfo("UsedAPVCyclePhaseCollection")
168  << " APVCyclePhaseCollection " << iConfig.getParameter<edm::InputTag>("apvPhaseCollection") << " used";
169 }
170 
172  // do anything here that needs to be done at desctruction time
173  // (e.g. close files, deallocate resources etc.)
174 }
175 
176 //
177 // member functions
178 //
179 
180 // ------------ method called to for each event ------------
182  using namespace edm;
183 
184  _nevents++;
185 
187  iEvent.getByToken(_historyProductToken, he);
188 
189  // improve the matchin between default and actual partitions
190 
191  (*_dbx)->Fill(he->deltaBX());
192  std::vector<std::pair<unsigned int, unsigned int> >::const_iterator indices = m_dbxindices.begin();
193  for (std::vector<TH1F**>::const_iterator dbxhist = m_dbxhistos.begin(); dbxhist != m_dbxhistos.end();
194  ++dbxhist, ++indices) {
195  (*(*dbxhist))->Fill(he->deltaBX(indices->first, indices->second));
196  }
197 
198  (*_bx)->Fill(iEvent.bunchCrossing() % 3564);
199  (*_orbit)->Fill(iEvent.orbitNumber());
200  if (_dbxvsbx && *_dbxvsbx)
201  (*_dbxvsbx)->Fill(iEvent.bunchCrossing() % 3564, he->deltaBX());
202  if (m_ewhdepth && *m_ewhdepth)
203  (*m_ewhdepth)->Fill(he->depth());
204 
206  iEvent.getByToken(_apvphasecollToken, apvphase);
207 
208  long long tbx = he->absoluteBX();
209  if (apvphase.isValid() && !apvphase.failedToGet()) {
210  const int thephase = apvphase->getPhase(_phasepart);
213  tbx -= thephase;
214  (*_bxincycle)->Fill(tbx % 70);
216  (*_dbxvsbxincycle)->Fill(tbx % 70, he->deltaBX());
218  (*_bxincyclevsbx)->Fill(iEvent.bunchCrossing() % 3564, tbx % 70);
220  (*_orbitvsbxincycle)->Fill(tbx % 70, iEvent.orbitNumber());
221 
222  } else {
223  LogDebug("InvalidPhase") << "Invalid APVCyclePhase value : " << _phasepart << " " << thephase;
224  }
225  }
226 }
227 
229  _rhm.beginRun(iRun);
230 
231  if (*_dbx) {
232  (*_dbx)->GetXaxis()->SetTitle("#DeltaBX");
233  }
234 
235  LogDebug("NomberOfHistos") << m_dbxhistos.size();
236  for (std::vector<TH1F**>::const_iterator dbxhist = m_dbxhistos.begin(); dbxhist != m_dbxhistos.end(); ++dbxhist) {
237  LogDebug("HistoPointer") << *dbxhist;
238  if (*(*dbxhist)) {
239  (*(*dbxhist))->GetXaxis()->SetTitle("#DeltaBX");
240  }
241  }
242  LogDebug("LabelDone") << "all labels set";
243 
244  if (*_bx) {
245  (*_bx)->GetXaxis()->SetTitle("BX");
246  }
247 
248  if (*_bxincycle) {
249  (*_bxincycle)->GetXaxis()->SetTitle("Event BX mod(70)");
250  }
251 
252  if (*_orbit) {
253  (*_orbit)->SetCanExtend(TH1::kXaxis);
254  (*_orbit)->GetXaxis()->SetTitle("time [Orb#]");
255  }
256 
257  LogDebug("StdPlotsDone") << "all labels in std plots set";
258 
260  (*_dbxvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)");
261  (*_dbxvsbxincycle)->GetYaxis()->SetTitle("#DeltaBX");
262  }
263 
264  if (_dbxvsbx && *_dbxvsbx) {
265  (*_dbxvsbx)->GetXaxis()->SetTitle("BX");
266  (*_dbxvsbx)->GetYaxis()->SetTitle("#DeltaBX");
267  }
268 
269  if (_bxincyclevsbx && *_bxincyclevsbx) {
270  (*_bxincyclevsbx)->GetXaxis()->SetTitle("BX");
271  (*_bxincyclevsbx)->GetYaxis()->SetTitle("Event BX mod(70)");
272  }
273 
275  (*_orbitvsbxincycle)->SetCanExtend(TH1::kYaxis);
276  (*_orbitvsbxincycle)->GetXaxis()->SetTitle("Event BX mod(70)");
277  (*_orbitvsbxincycle)->GetYaxis()->SetTitle("time [Orb#]");
278  }
279 
280  if (m_ewhdepth && *m_ewhdepth) {
281  (*m_ewhdepth)->GetXaxis()->SetTitle("Depth");
282  }
283 }
284 
286 
287 // ------------ method called once each job just before starting event loop ------------
289 
290 // ------------ method called once each job just after ending the event loop ------------
291 void EventTimeDistribution::endJob() { edm::LogInfo("EndOfJob") << _nevents << " analyzed events"; }
292 
293 //define this as a plug-in
#define LogDebug(id)
unsigned int depth() const
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
long long deltaBX(const unsigned int ev2, const unsigned int ev1) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define nullptr
int bunchCrossing() const
Definition: EventBase.h:64
long long absoluteBX(const unsigned int ev1) const
edm::EDGetTokenT< APVCyclePhaseCollection > _apvphasecollToken
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< EventWithHistory > _historyProductToken
void endRun(const edm::Run &, const edm::EventSetup &) override
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
std::vector< std::pair< unsigned int, unsigned int > > m_dbxindices
const int getPhase(const std::string partition) const
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
const unsigned int m_LSfrac
int orbitNumber() const
Definition: EventBase.h:65
bool isValid() const
Definition: HandleBase.h:70
RunHistogramManager _rhm
const std::string _phasepart
EventTimeDistribution(const edm::ParameterSet &)
bool failedToGet() const
Definition: HandleBase.h:72
void beginRun(const edm::Run &iRun)
const unsigned int m_maxLS
HLT enums.
void analyze(const edm::Event &, const edm::EventSetup &) override
TH2F ** makeTH2F(const char *name, const char *title, const unsigned int nbinx, const double xmin, const double xmax, const unsigned int nbiny, const double ymin, const double ymax)
std::vector< TH1F ** > m_dbxhistos
Definition: Run.h:45