CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BxTiming.cc
Go to the documentation of this file.
2 #include <cstdio>
3 
5 
6  verbose_ = iConfig.getUntrackedParameter<int>("VerboseFlag",0);
7  if(verbose())
8  std::cout << "BxTiming::BxTiming()...\n" << std::flush;
9 
10  fedRef_ = iConfig.getUntrackedParameter<int>("ReferenceFedId",813);
12  ("FedSource",edm::InputTag("source"));
13  fedSource_token_ = consumes<FEDRawDataCollection>(iConfig.getUntrackedParameter<edm::InputTag>
14  ("FedSource",edm::InputTag("source")));
16  ("GtSource",edm::InputTag("gtUnpack"));
17  gtSource_token_ = consumes<L1GlobalTriggerReadoutRecord>(iConfig.getUntrackedParameter<edm::InputTag>
18  ("GtSource",edm::InputTag("gtUnpack")));
20  ("HistFile","");
22  ("HistFolder", "L1T/BXSynch");
23 
24  runInFF_ = iConfig.getUntrackedParameter<bool> ("RunInFilterFarm", true);
25  if(runInFF_) histFolder_ = "L1T/BXSynch_EvF";
26  if(verbose())
27  std::cout << "Filter farm run setting?" << runInFF_
28  << "\n" << std::flush;
29 
30  listGtBits_ = iConfig.getUntrackedParameter<std::vector<int> > ("GtBitList", std::vector<int>(1,0));
31  if(listGtBits_.size()==1 && listGtBits_.at(0)==-1) {
32  int ngtbits = 128;
33  listGtBits_.reserve(ngtbits);
34  for(int i=0; i<ngtbits; i++)
35  listGtBits_[i]=i;
36  }
37 
38  if(verbose()) {
39  std::cout << "BxTiming: gt bits set for timing dqm:";
40  std::cout << "nbits:" << listGtBits_.size() << " list: " ;
41  for(size_t i=0; i!=listGtBits_.size(); i++)
42  std::cout << listGtBits_.at(i) << " " ;
43  std::cout << "\n" << std::flush;
44  }
45 
47 
48  dbe = NULL;
49  if (iConfig.getUntrackedParameter<bool>("DQMStore", false)) {
51  dbe->setVerbose(0);
52  }
53 
54  if(dbe!=NULL)
55  dbe->setCurrentFolder(histFolder_);
56 
57  nEvt_ = 0;
58 
59  if(verbose())
60  std::cout << "BxTiming::BxTiming constructor...done.\n" << std::flush;
61 }
62 
64 
65 void
67 
68  if(verbose())
69  std::cout << "BxTiming::beginJob() start\n" << std::flush;
70 
71  if(verbose())
72  std::cout << "BxTiming::beginJob() end.\n" << std::flush;
73 }
74 
75 
76 
77 void
78 BxTiming::beginRun(edm::Run const& iRun, edm::EventSetup const& iSetup)
79 {
80  if(dbe) {
82  }
83 
85  for(int i=0; i<nfed_;i++) {
86  nBxDiff[i][0]=0; nBxDiff[i][1]=nbig_; nBxDiff[i][2]=-1*nbig_;
87  nBxOccy[i][0]=0; nBxOccy[i][1]=nbig_; nBxOccy[i][2]=-1*nbig_;
88  }
89 
90  std::string lbl("");
91  std::string SysLabel[NSYS] = {
92  "PreShower", "ECAL", "HCAL", "GCT", "CSCTPG", "CSCTF", "DTTPG", "DTTF", "RPC", "GT"
93  };
94 
95  typedef std::pair<int, int> FEDRange;
96 
97  std::pair<int,int> fedRange[NSYS] = {
108  };
109  for(int i=0; i<NSYS; i++) fedRange_[i]=fedRange[i];
110 
111 
112  int fedRefSys=-1;
113  for(int i=0; i<NSYS; i++)
115  {fedRefSys=i; break;}
116  std::string refName("");
117  std::string spreadLabel[nspr_] = {"Spread","Min", "Max"};
118  if(fedRefSys>=0)
119  refName+=SysLabel[fedRefSys];
120  else
121  refName+=fedRef_;
122 
124 
125  const int dbx = nbig_;
126 
127  if(dbe) {
128 
130 
131  hBxDiffAllFed = dbe->bookProfile("BxDiffAllFed", "BxDiffAllFed",
132  nfed_ + 1, -0.5, nfed_+0.5,
133  2*dbx+1, -1*dbx-0.5,dbx+0.5
134  );
135 
136  for(int i=0; i<nspr_; i++) {
137  lbl.clear();lbl+="BxDiffAllFed";lbl+=spreadLabel[i];
138  hBxDiffAllFedSpread[i] = dbe->book1D(lbl.data(),lbl.data(), nfed_ + 1, -0.5, nfed_+0.5);
139  lbl.clear();lbl+="BxOccyAllFed";lbl+=spreadLabel[i];
140  hBxOccyAllFedSpread[i] = dbe->book1D(lbl.data(),lbl.data(), nfed_ + 1, -0.5, nfed_+0.5);
141  }
142 
143  lbl.clear();lbl+="BxOccyAllFed";
144  hBxOccyAllFed = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
145 
146  }
147 
148  // following histos defined only when not runing in the ff
149  if(dbe && !runInFF_) {
150 
152 
153  for(int i=0; i<NSYS; i++) {
154  lbl.clear();lbl+=SysLabel[i];lbl+="FedBxDiff";
155  int nfeds = fedRange_[i].second - fedRange_[i].first + 1;
156  nfeds = (nfeds>0)? nfeds:1;
157  hBxDiffSysFed[i] = dbe->bookProfile(lbl.data(),lbl.data(), nfeds,
158  fedRange_[i].first-0.5, fedRange_[i].second+0.5,
159  2*dbx+1,-1*dbx-0.5,dbx+0.5);
160 
161  }
162 
163 
164  lbl.clear();lbl+="BxOccyAllFed";
166  dbe->setCurrentFolder(histFolder_+"/SingleFed");
167  for(int i=0; i<nfed_; i++) {
168  lbl.clear(); lbl+="BxOccyOneFed";
169  char *ii = new char[1000]; std::sprintf(ii,"%d",i);lbl+=ii;
170  hBxOccyOneFed[i] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
171  delete [] ii;
172  }
173 
175  for(int i=0; i<nttype_; i++) {
176  lbl.clear();lbl+="BxOccyGtTrigType";
177  char *ii = new char[10]; std::sprintf(ii,"%d",i+1);lbl+=ii;
178  hBxOccyGtTrigType[i] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
179  delete [] ii;
180  }
181 
182  dbe->setCurrentFolder(histFolder_+"/SingleBit");
183  for(int i=0; i<NSYS; i++) {
184  hBxOccyTrigBit[i] = new MonitorElement*[listGtBits_.size()];
185  for(size_t j=0; j<listGtBits_.size(); j++) {
186  lbl.clear();lbl+=SysLabel[i];lbl+="BxOccyGtBit";
187  char *ii = new char[1000]; std::sprintf(ii,"%d",listGtBits_.at(j)); lbl+=ii;
188  hBxOccyTrigBit[i][j] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
189  delete [] ii;
190  }
191  }
192 
193  }
194 
196  hBxDiffAllFed->setAxisTitle("FED ID",1);
197  lbl.clear(); lbl+="BX(fed)-BX("; lbl+=refName; lbl+=")";
198  hBxDiffAllFed->setAxisTitle(lbl,2);
199  for(int i=0; i<nspr_; i++) {
200  lbl.clear(); lbl+="BX(fed)-BX("; lbl+=refName; lbl+=") "+spreadLabel[i];
201  hBxDiffAllFedSpread[i]->setAxisTitle("FED ID",1);
203  lbl.clear(); lbl+="Bx FED occupancy"; lbl+=" "; lbl+=spreadLabel[i];
204  hBxOccyAllFedSpread[i]->setAxisTitle("FED ID",1);
206  }
207 
208  hBxOccyAllFed->setAxisTitle("bx",1);
209  lbl.clear(); lbl+="Combined FED occupancy";
210  hBxOccyAllFed->setAxisTitle(lbl,2);
211 
212  // skip next if running in filter farm
213  if(runInFF_)
214  return;
215 
216  for(int i=0; i<NSYS; i++) {
217  lbl.clear(); lbl+=SysLabel[i]; lbl+=" FED ID";
218  hBxDiffSysFed[i]->setAxisTitle(lbl,1);
219  lbl.clear(); lbl+="BX("; lbl+=SysLabel[i]; lbl+=")-BX(";lbl+=refName; lbl+=")";
220  hBxDiffSysFed[i]->setAxisTitle(lbl,2);
221  }
222 
223  for(int i=0; i<nfed_; i++) {
224  hBxOccyOneFed[i] ->setAxisTitle("bx",1);
225  lbl.clear(); lbl+=" FED "; char *ii = new char[1000]; std::sprintf(ii,"%d",i);lbl+=ii; lbl+=" occupancy";
226  hBxOccyOneFed[i] ->setAxisTitle(lbl,2);
227  delete [] ii;
228  }
229  for(int i=0; i<nttype_; i++) {
230  hBxOccyGtTrigType[i]->setAxisTitle("bx",1);
231  lbl.clear(); lbl+="GT occupancy for trigger type "; char *ii = new char[10]; std::sprintf(ii,"%d",i+1);lbl+=ii;
233  delete [] ii;
234  }
235 
236  for(int i=0; i<NSYS; i++) {
237  for(size_t j=0; j<listGtBits_.size(); j++) {
238  hBxOccyTrigBit[i][j]->setAxisTitle("bx",1);
239  lbl.clear();lbl+=SysLabel[i];lbl+=" Bx occupancy for Trigger bit ";
240  char *ii = new char[10]; std::sprintf(ii,"%d",listGtBits_.at(j)); lbl+=ii;
241  hBxOccyTrigBit[i][j]->setAxisTitle(lbl,2);
242  delete [] ii;
243  }
244  }
245 }
246 
247 void
249 
250  if(verbose())
251  std::cout << "BxTiming::endJob Nevents: " << nEvt_ << "\n" << std::flush;
252 
253  if(histFile_.size()!=0 && dbe)
254  dbe->save(histFile_);
255 
256  if(verbose())
257  std::cout << "BxTiming::endJob() end.\n" << std::flush;
258 }
259 
260 
261 // ------------ method called to for each event ------------
262 void
264 
265  if(verbose())
266  std::cout << "BxTiming::analyze() start\n" << std::flush;
267 
268  nEvt_++;
269 
272  iEvent.getByToken(fedSource_token_, rawdata);
273 
274  if (!rawdata.isValid()) {
275 
276  if (verbose())
277  std::cout
278  << "BxTiming::analyze() | FEDRawDataCollection with input tag "
279  << fedSource_ << " not found.";
280 
281  return;
282  }
283 
284  // get the GT bits
286  iEvent.getByToken(gtSource_token_, gtdata);
287  std::vector<bool> gtbits;
288  int ngtbits = 128;
289  gtbits.reserve(ngtbits); for(int i=0; i<ngtbits; i++) gtbits[i]=false;
290  if(gtdata.isValid())
291  gtbits = gtdata->decisionWord();
292 
293  if(gtbits.size()==0) {
294  gtbits.push_back(true); // gtdata->decision();
295  if(verbose())
296  std::cout << "BxTiming::analyze() | unexpected empty decision bits!";
297  }
298 
299  if(verbose()) {
300  std::cout << "BxTiming::analyze() gt data valid:" << (int)(gtdata.isValid()?0:1)
301  << " decision word size:" << (int)(gtbits.size()) << " bits: ";
302  for(size_t i=0; i!=gtbits.size(); i++) {
303  int ii = gtbits.at(i)? 1:0;
304  std::cout << ii;
305  }
306  std::cout << ".\n" << std::flush;
307  }
308 
309 
310  // get reference bx
311  int bxRef = FEDHeader(rawdata->FEDData(fedRef_).data()).bxID();
312 
313  // triggerType
314  //trigger types: physics (1), calibration (2), random (3), traced physics (5), test (6)
315  int ttype = FEDHeader(rawdata->FEDData(812).data()).triggerType();
316 
317  // loop over feds
318  for (int i = 0; i<FEDNumbering::MAXFEDID+1; i++){
319  const FEDRawData& data = rawdata->FEDData(i);
320  size_t size=data.size();
321 
322  if(!size) continue;
323  FEDHeader header(data.data());
324  //int lvl1id = header.lvl1ID();//Level-1 event number generated by the TTC system
325  int bx = header.bxID(); // The bunch crossing number
326 
327  int bxDiff = calcBxDiff(bx,bxRef); // deviation from reference bx
328 
329  //min
330  if(nBxDiff[i][1]>bxDiff) nBxDiff[i][1] = bxDiff;
331  if(nBxOccy[i][1]>bx ) nBxOccy[i][1] = bx;
332  //max
333  if(nBxDiff[i][2]<bxDiff) nBxDiff[i][2] = bxDiff;
334  if(nBxOccy[i][2]<bx ) nBxOccy[i][2] = bx;
335 
336  if(verbose())
337  std::cout << " fed:" << i
338  << " bx:" << bx
339  << " bxRef:" << bxRef
340  << " diff:" << bxDiff
341  << " nBxDiff"<<" del:"<<nBxDiff[i][0]<<" min:"<<nBxDiff[i][1]<<" max:"<<nBxDiff[i][2]
342  << " nBxOccy"<<" del:"<<nBxOccy[i][0]<<" min:"<<nBxOccy[i][1]<<" max:"<<nBxOccy[i][2]
343  << "\n" << std::flush;
344 
345  hBxDiffAllFed->Fill(i,bxDiff);
346 
347  //if(ttype==1) //skip if not a physics trigger
348  hBxOccyAllFed->Fill(bx);
349 
350 
351  // done if running in filter farm
352  if(runInFF_)
353  continue;
354 
355  for(int j=0; j<NSYS; j++)
356  if(i>=fedRange_[j].first && i<=fedRange_[j].second)
357  hBxDiffSysFed[j]->Fill(i,bxDiff);
358 
359  for(size_t k=0; k!=listGtBits_.size(); k++) {
360  if((int)gtbits.size() <= listGtBits_.at(k)) {
361  if(verbose())
362  std::cout << "BxTiming analyze | problem with vector size!\n" << std::endl;
363  continue;
364  }
365  else if(!gtbits.at(listGtBits_.at(k)))
366  continue;
367  for(int j=0; j<NSYS; j++) {
368  if(i>=fedRange_[j].first && i<=fedRange_[j].second) {
369  hBxOccyTrigBit[j][k]->Fill(bx);
370  }
371  }
372  }
373 
374  if(i>=fedRange_[GLT].first && i<=fedRange_[GLT].second) //GT fed
375  if(ttype<nttype_)
376  hBxOccyGtTrigType[ttype-1]->Fill(bx);
377 
378  if(ttype!=1) continue; //skip if not a physics trigger
379  //hBxOccyAllFed->Fill(bx);
380  hBxOccyOneFed[i]->Fill(bx);
381 
382  }
383 
384  for(int i=0; i<nfed_;i++) {
385  nBxDiff[i][0]=nBxDiff[i][2]-nBxDiff[i][1];
386  nBxOccy[i][0]=nBxOccy[i][2]-nBxOccy[i][1];
387  if(nBxDiff[i][0]<0 || nBxOccy[i][0]<0) continue;
388  for(int j=0; j<nspr_; j++) {
389  hBxDiffAllFedSpread[j]->setBinContent(i,nBxDiff[i][j]);
390  hBxOccyAllFedSpread[j]->setBinContent(i,nBxOccy[i][j]);
391  }
392  if(verbose())
393  std::cout << "BxTiming fed:" << i
394  << " Bx-Bx(" << fedRef_ << ")::"
395  << " del:" << nBxDiff[i][0]
396  << " min:" << nBxDiff[i][1]
397  << " max:" << nBxDiff[i][2]
398  << " Occy: "
399  << " del:" << nBxOccy[i][0]
400  << " min:" << nBxOccy[i][1]
401  << " max:" << nBxOccy[i][2]
402  <<"\n" << std::flush;
403 
404  }
405 
406 
407  if(verbose())
408  std::cout << "BxTiming::analyze() end.\n" << std::flush;
409 }
410 
411 //----------------------------------------------------------------------
412 
413 int
414 BxTiming::calcBxDiff(int bx1, int bx2)
415 {
416  int diff = bx1 - bx2;
417 
418  while (diff < -half_norb_)
419  diff += norb_;
420 
421  while (diff > half_norb_)
422  diff -= norb_;
423 
424  return diff;
425 }
426 
427 //----------------------------------------------------------------------
static const int nttype_
Definition: BxTiming.h:85
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
static const int norb_
Definition: BxTiming.h:81
DQMStore * dbe
Definition: BxTiming.h:75
MonitorElement * hBxDiffAllFed
histograms
Definition: BxTiming.h:101
void setBinContent(int binx, double content)
set content of bin (1-D)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:872
MonitorElement * hBxDiffAllFedSpread[nspr_]
Definition: BxTiming.h:106
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
std::vector< int > listGtBits_
Definition: BxTiming.h:87
int verbose()
Definition: BxTiming.h:55
std::string histFile_
Definition: BxTiming.h:69
#define NULL
Definition: scimark2.h:8
int ii
Definition: cuy.py:588
MonitorElement * hBxOccyAllFed
Definition: BxTiming.h:103
MonitorElement ** hBxOccyOneFed
Definition: BxTiming.h:104
size_t size() const
Lenght of the data buffer in bytes.
Definition: FEDRawData.h:47
static const int nbig_
Definition: BxTiming.h:84
edm::EDGetTokenT< L1GlobalTriggerReadoutRecord > gtSource_token_
Definition: BxTiming.h:51
void Fill(long long x)
int nBxOccy[1500][nspr_]
Definition: BxTiming.h:98
U second(std::pair< T, U > const &p)
virtual void beginJob(void)
Definition: BxTiming.cc:66
MonitorElement * hBxOccyGtTrigType[nttype_]
Definition: BxTiming.h:109
int iEvent
Definition: GenABIO.cc:243
int nEvt_
Definition: BxTiming.h:66
int verbose_
Definition: BxTiming.h:54
int nBxDiff[1500][nspr_]
Definition: BxTiming.h:97
int j
Definition: DBlmapReader.cc:9
virtual void beginRun(edm::Run const &iRun, edm::EventSetup const &iSetup)
Definition: BxTiming.cc:78
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2296
static const int half_norb_
Definition: BxTiming.h:82
MonitorElement * bookProfile(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const char *option="s")
Definition: DQMStore.cc:1186
void setVerbose(unsigned level)
Definition: DQMStore.cc:548
edm::InputTag fedSource_
Definition: BxTiming.h:48
bool first
Definition: L1TdeRCT.cc:79
bool isValid() const
Definition: HandleBase.h:76
virtual void analyze(const edm::Event &, const edm::EventSetup &)
Definition: BxTiming.cc:263
int k[5][pyjets_maxn]
bool runInFF_
Definition: BxTiming.h:78
std::string histFolder_
Definition: BxTiming.h:72
edm::InputTag gtSource_
Definition: BxTiming.h:50
edm::EDGetTokenT< FEDRawDataCollection > fedSource_token_
Definition: BxTiming.h:49
int bxID()
The bunch crossing number.
Definition: FEDHeader.cc:24
int fedRef_
Definition: BxTiming.h:93
MonitorElement * hBxDiffSysFed[NSYS]
Definition: BxTiming.h:102
MonitorElement ** hBxOccyTrigBit[NSYS]
Definition: BxTiming.h:110
int nfed_
Definition: BxTiming.h:92
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
const unsigned char * data() const
Return a const pointer to the beginning of the data buffer.
Definition: FEDRawData.cc:28
int calcBxDiff(int bx1, int bx2)
Definition: BxTiming.cc:414
tuple cout
Definition: gather_cfg.py:121
virtual void endJob()
Definition: BxTiming.cc:248
dictionary rawdata
Definition: lumiPlot.py:393
std::pair< int, int > fedRange_[NSYS]
Definition: BxTiming.h:91
MonitorElement * hBxOccyAllFedSpread[nspr_]
Definition: BxTiming.h:107
void setAxisTitle(const std::string &title, int axis=1)
set x-, y- or z-axis title (axis=1, 2, 3 respectively)
~BxTiming()
Definition: BxTiming.cc:63
BxTiming(const edm::ParameterSet &)
Definition: BxTiming.cc:4
tuple size
Write out results.
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:584
static const int nspr_
Definition: BxTiming.h:96
Definition: Run.h:41