CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
Level1TriggerScalers.cc
Go to the documentation of this file.
1 /*
2  * File: DataFormats/Scalers/src/Level1TriggerScalers.cc (W.Badgett)
3  */
4 
7 
8 #include <iostream>
9 #include <ctime>
10 #include <cstdio>
11 
13  : version_(0),
14  trigType_(0),
15  eventID_(0),
16  sourceID_(0),
17  bunchNumber_(0),
18  collectionTime_(0, 0),
19  lumiSegmentNr_(0),
20  lumiSegmentOrbits_(0),
21  orbitNr_(0),
22  gtResets_(0),
23  bunchCrossingErrors_(0),
24  gtTriggers_(0),
25  gtEvents_(0),
26  gtTriggersRate_((float)0.0),
27  gtEventsRate_((float)0.0),
28  prescaleIndexAlgo_(0),
29  prescaleIndexTech_(0),
30  collectionTimeLumiSeg_(0, 0),
31  lumiSegmentNrLumiSeg_(0),
32  triggersPhysicsGeneratedFDL_(0),
33  triggersPhysicsLost_(0),
34  triggersPhysicsLostBeamActive_(0),
35  triggersPhysicsLostBeamInactive_(0),
36  l1AsPhysics_(0),
37  l1AsRandom_(0),
38  l1AsTest_(0),
39  l1AsCalibration_(0),
40  deadtime_(0),
41  deadtimeBeamActive_(0),
42  deadtimeBeamActiveTriggerRules_(0),
43  deadtimeBeamActiveCalibration_(0),
44  deadtimeBeamActivePrivateOrbit_(0),
45  deadtimeBeamActivePartitionController_(0),
46  deadtimeBeamActiveTimeSlot_(0),
47  gtAlgoCounts_(nLevel1Triggers),
48  gtTechCounts_(nLevel1TestTriggers),
49  lastOrbitCounter0_(0),
50  lastTestEnable_(0),
51  lastResync_(0),
52  lastStart_(0),
53  lastEventCounter0_(0),
54  lastHardReset_(0),
55  spare0_(0ULL),
56  spare1_(0ULL),
57  spare2_(0ULL) {}
58 
59 Level1TriggerScalers::Level1TriggerScalers(const unsigned char* rawData) {
61 
62  struct ScalersEventRecordRaw_v5 const* raw = reinterpret_cast<struct ScalersEventRecordRaw_v5 const*>(rawData);
63 
64  trigType_ = (raw->header >> 56) & 0xFULL;
65  eventID_ = (raw->header >> 32) & 0x00FFFFFFULL;
66  sourceID_ = (raw->header >> 8) & 0x00000FFFULL;
67  bunchNumber_ = (raw->header >> 20) & 0xFFFULL;
68 
69  version_ = raw->version;
70  if (version_ >= 3) {
71  collectionTime_.set_tv_sec(static_cast<long>(raw->trig.collectionTime_sec));
73 
76  orbitNr_ = raw->trig.orbitNr;
77  gtResets_ = raw->trig.gtResets;
80  gtEvents_ = raw->trig.gtEvents;
85 
88 
94 
97  l1AsTest_ = raw->trig.l1AsTest;
99  deadtime_ = raw->trig.deadtime;
106 
107  for (int i = 0; i < ScalersRaw::N_L1_TRIGGERS_v1; i++) {
108  gtAlgoCounts_.push_back(raw->trig.gtAlgoCounts[i]);
109  }
110 
111  for (int i = 0; i < ScalersRaw::N_L1_TEST_TRIGGERS_v1; i++) {
112  gtTechCounts_.push_back(raw->trig.gtTechCounts[i]);
113  }
114 
115  if (version_ >= 5) {
118  lastResync_ = raw->lastResync;
119  lastStart_ = raw->lastStart;
122  spare0_ = raw->spare[0];
123  spare1_ = raw->spare[1];
124  spare2_ = raw->spare[2];
125  } else {
126  lastOrbitCounter0_ = 0UL;
127  lastTestEnable_ = 0UL;
128  lastResync_ = 0UL;
129  lastStart_ = 0UL;
130  lastEventCounter0_ = 0UL;
131  lastHardReset_ = 0UL;
132  spare0_ = 0ULL;
133  spare1_ = 0ULL;
134  spare2_ = 0ULL;
135  }
136  }
137 }
138 
140 
141 double Level1TriggerScalers::rateLS(unsigned int counts) { return (rateLS(counts, firstShortLSRun)); }
142 
143 double Level1TriggerScalers::rateLS(unsigned long long counts) { return (rateLS(counts, firstShortLSRun)); }
144 
145 double Level1TriggerScalers::rateLS(unsigned int counts, int runNumber) {
146  unsigned long long counts64 = (unsigned long long)counts;
147  return (rateLS(counts64, runNumber));
148 }
149 
150 double Level1TriggerScalers::rateLS(unsigned long long counts, int runNumber) {
151  double rate;
152  if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
153  rate = ((double)counts) / 23.31040958083832;
154  } else {
155  rate = ((double)counts) / 93.24163832335329;
156  }
157  return (rate);
158 }
159 
160 double Level1TriggerScalers::percentLS(unsigned long long counts) { return (percentLS(counts, firstShortLSRun)); }
161 
162 double Level1TriggerScalers::percentLS(unsigned long long counts, int runNumber) {
163  double percent;
164  if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
165  percent = ((double)counts) / 9342812.16;
166  } else {
167  percent = ((double)counts) / 37371248.64;
168  }
169  if (percent > 100.0000) {
170  percent = 100.0;
171  }
172  return (percent);
173 }
174 
175 double Level1TriggerScalers::percentLSActive(unsigned long long counts) {
176  return (percentLSActive(counts, firstShortLSRun));
177 }
178 
179 double Level1TriggerScalers::percentLSActive(unsigned long long counts, int runNumber) {
180  double percent;
181  if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
182  percent = ((double)counts) / 7361003.52;
183  } else {
184  percent = ((double)counts) / 29444014.08;
185  }
186  if (percent > 100.0000) {
187  percent = 100.0;
188  }
189  return (percent);
190 }
191 
193 std::ostream& operator<<(std::ostream& s, Level1TriggerScalers const& c) {
194  s << "Level1TriggerScalers Version:" << c.version() << " SourceID: " << c.sourceID() << std::endl;
195  constexpr size_t kLineBufferSize = 164;
196  char line[kLineBufferSize];
197  char zeitHeaven[128];
198  struct tm* horaHeaven;
199 
200  snprintf(line,
201  kLineBufferSize,
202  " TrigType: %d EventID: %d BunchNumber: %d",
203  c.trigType(),
204  c.eventID(),
205  c.bunchNumber());
206  s << line << std::endl;
207 
208  struct timespec secondsToHeaven = c.collectionTime();
209  horaHeaven = gmtime(&secondsToHeaven.tv_sec);
210  strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
211  snprintf(line, kLineBufferSize, " CollectionTime: %s.%9.9d", zeitHeaven, (int)secondsToHeaven.tv_nsec);
212  s << line << std::endl;
213 
214  snprintf(line,
215  kLineBufferSize,
216  " LumiSegmentNr: %10u LumiSegmentOrbits: %10u",
217  c.lumiSegmentNr(),
218  c.lumiSegmentOrbits());
219  s << line << std::endl;
220 
221  snprintf(line,
222  kLineBufferSize,
223  " LumiSegmentNrLumiSeg: %10u OrbitNr: %10u ",
225  c.orbitNr());
226  s << line << std::endl;
227 
228  snprintf(line,
229  kLineBufferSize,
230  " GtResets: %10u BunchCrossingErrors: %10u",
231  c.gtResets(),
232  c.bunchCrossingErrors());
233  s << line << std::endl;
234 
235  snprintf(line,
236  kLineBufferSize,
237  " PrescaleIndexAlgo: %10d PrescaleIndexTech: %10d",
238  c.prescaleIndexAlgo(),
239  c.prescaleIndexTech());
240  s << line << std::endl;
241 
242  snprintf(
243  line, kLineBufferSize, " GtTriggers: %20llu %22.3f Hz", c.gtTriggers(), c.gtTriggersRate());
244  s << line << std::endl;
245 
246  snprintf(line, kLineBufferSize, " GtEvents: %20llu %22.3f Hz", c.gtEvents(), c.gtEventsRate());
247  s << line << std::endl;
248 
249  secondsToHeaven = c.collectionTimeLumiSeg();
250  horaHeaven = gmtime(&secondsToHeaven.tv_sec);
251  strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
252  snprintf(line, kLineBufferSize, " CollectionTimeLumiSeg: %s.%9.9d", zeitHeaven, (int)secondsToHeaven.tv_nsec);
253  s << line << std::endl;
254 
255  snprintf(line,
256  kLineBufferSize,
257  " TriggersPhysicsGeneratedFDL: %20llu %22.3f Hz",
260  s << line << std::endl;
261 
262  snprintf(line,
263  kLineBufferSize,
264  " TriggersPhysicsLost: %20llu %22.3f Hz",
267  s << line << std::endl;
268 
269  snprintf(line,
270  kLineBufferSize,
271  " TriggersPhysicsLostBeamActive: %20llu %22.3f Hz",
274  s << line << std::endl;
275 
276  snprintf(line,
277  kLineBufferSize,
278  " TriggersPhysicsLostBeamInactive: %20llu %22.3f Hz",
281  s << line << std::endl;
282 
283  snprintf(line,
284  kLineBufferSize,
285  " L1AsPhysics: %20llu %22.3f Hz",
286  c.l1AsPhysics(),
288  s << line << std::endl;
289 
290  snprintf(line,
291  kLineBufferSize,
292  " L1AsRandom: %20llu %22.3f Hz",
293  c.l1AsRandom(),
295  s << line << std::endl;
296 
297  snprintf(line,
298  kLineBufferSize,
299  " L1AsTest: %20llu %22.3f Hz",
300  c.l1AsTest(),
302  s << line << std::endl;
303 
304  snprintf(line,
305  kLineBufferSize,
306  " L1AsCalibration: %20llu %22.3f Hz",
307  c.l1AsCalibration(),
309  s << line << std::endl;
310 
311  snprintf(line,
312  kLineBufferSize,
313  " Deadtime: %20llu %17.3f%%",
314  c.deadtime(),
316  s << line << std::endl;
317 
318  snprintf(line,
319  kLineBufferSize,
320  " DeadtimeBeamActive: %20llu %17.3f%%",
321  c.deadtimeBeamActive(),
323  s << line << std::endl;
324 
325  snprintf(line,
326  kLineBufferSize,
327  " DeadtimeBeamActiveTriggerRules: %20llu %17.3f%%",
330  s << line << std::endl;
331 
332  snprintf(line,
333  kLineBufferSize,
334  " DeadtimeBeamActiveCalibration: %20llu %17.3f%%",
337  s << line << std::endl;
338 
339  snprintf(line,
340  kLineBufferSize,
341  " DeadtimeBeamActivePrivateOrbit: %20llu %17.3f%%",
344  s << line << std::endl;
345 
346  snprintf(line,
347  kLineBufferSize,
348  " DeadtimeBeamActivePartitionController:%20llu %17.3f%%",
351  s << line << std::endl;
352 
353  snprintf(line,
354  kLineBufferSize,
355  " DeadtimeBeamActiveTimeSlot: %20llu %17.3f%%",
358  s << line << std::endl;
359 
360  s << "Physics GtAlgoCounts" << std::endl;
361  const std::vector<unsigned int> gtAlgoCounts = c.gtAlgoCounts();
362  int length = gtAlgoCounts.size() / 4;
363  for (int i = 0; i < length; i++) {
364  snprintf(line,
365  kLineBufferSize,
366  " %3.3d: %10u %3.3d: %10u %3.3d: %10u %3.3d: %10u",
367  i,
368  gtAlgoCounts[i],
369  (i + length),
370  gtAlgoCounts[i + length],
371  (i + (length * 2)),
372  gtAlgoCounts[i + (length * 2)],
373  (i + (length * 3)),
374  gtAlgoCounts[i + (length * 3)]);
375  s << line << std::endl;
376  }
377 
378  s << "Test GtTechCounts" << std::endl;
379  const std::vector<unsigned int> gtTechCounts = c.gtTechCounts();
380  length = gtTechCounts.size() / 4;
381  for (int i = 0; i < length; i++) {
382  snprintf(line,
383  kLineBufferSize,
384  " %3.3d: %10u %3.3d: %10u %3.3d: %10u %3.3d: %10u",
385  i,
386  gtTechCounts[i],
387  (i + length),
388  gtTechCounts[i + length],
389  (i + (length * 2)),
390  gtTechCounts[i + (length * 2)],
391  (i + (length * 3)),
392  gtTechCounts[i + (length * 3)]);
393  s << line << std::endl;
394  }
395 
396  if (c.version() >= 5) {
397  snprintf(line, kLineBufferSize, " LastOrbitCounter0: %10u 0x%8.8X", c.lastOrbitCounter0(), c.lastOrbitCounter0());
398  s << line << std::endl;
399 
400  snprintf(line, kLineBufferSize, " LastTestEnable: %10u 0x%8.8X", c.lastTestEnable(), c.lastTestEnable());
401  s << line << std::endl;
402 
403  snprintf(line, kLineBufferSize, " LastResync: %10u 0x%8.8X", c.lastResync(), c.lastResync());
404  s << line << std::endl;
405 
406  snprintf(line, kLineBufferSize, " LastStart: %10u 0x%8.8X", c.lastStart(), c.lastStart());
407  s << line << std::endl;
408 
409  snprintf(line, kLineBufferSize, " LastEventCounter0: %10u 0x%8.8X", c.lastEventCounter0(), c.lastEventCounter0());
410  s << line << std::endl;
411 
412  snprintf(line, kLineBufferSize, " LastHardReset: %10u 0x%8.8X", c.lastHardReset(), c.lastHardReset());
413  s << line << std::endl;
414  }
415 
416  return s;
417 }
unsigned long long deadtimeBeamActivePartitionController_
unsigned int gtTechCounts[ScalersRaw::N_L1_TEST_TRIGGERS_v1]
Definition: ScalersRaw.h:113
unsigned long long l1AsTest
Definition: ScalersRaw.h:102
unsigned long long spare1_
unsigned int orbitNr
Definition: ScalersRaw.h:83
static double rateLS(unsigned long long counts)
const edm::EventSetup & c
unsigned long long deadtimeBeamActiveTriggerRules_
unsigned long long gtEvents() const
unsigned long long deadtimeBeamActive() const
unsigned long long deadtimeBeamActivePrivateOrbit() const
unsigned long long triggersPhysicsGeneratedFDL() const
unsigned int bunchNumber() const
unsigned int lastTestEnable() const
void set_tv_nsec(long value)
Definition: TimeSpec.h:19
unsigned int lastResync
Definition: ScalersRaw.h:230
unsigned long long spare[ScalersRaw::N_SPARE_v5]
Definition: ScalersRaw.h:234
unsigned long long l1AsRandom() const
unsigned long long triggersPhysicsLostBeamInactive_
unsigned long long deadtimeBeamActiveTimeSlot_
unsigned int lastResync() const
std::vector< unsigned int > gtTechCounts() const
unsigned long long l1AsRandom
Definition: ScalersRaw.h:101
unsigned int trigType() const
unsigned long long triggersPhysicsLost() const
unsigned long long triggersPhysicsLostBeamActive() const
unsigned int collectionTime_sec
Definition: ScalersRaw.h:79
std::ostream & operator<<(std::ostream &out, const ALILine &li)
Definition: ALILine.cc:167
unsigned int lastStart() const
unsigned int collectionTime_nsec
Definition: ScalersRaw.h:80
unsigned long long gtTriggers_
unsigned int gtResets() const
unsigned long long l1AsCalibration() const
std::vector< unsigned int > gtTechCounts_
std::vector< unsigned int > gtAlgoCounts_
unsigned long long l1AsCalibration
Definition: ScalersRaw.h:103
unsigned long long deadtimeBeamActivePartitionController
Definition: ScalersRaw.h:109
unsigned long long deadtimeBeamActiveTimeSlot
Definition: ScalersRaw.h:110
unsigned long long triggersPhysicsLostBeamInactive
Definition: ScalersRaw.h:99
unsigned long long deadtime_
unsigned long long deadtimeBeamActiveTriggerRules
Definition: ScalersRaw.h:106
unsigned long long l1AsTest_
unsigned int bunchCrossingErrors() const
unsigned int orbitNr() const
unsigned int lastEventCounter0() const
struct TriggerScalersRaw_v3 trig
Definition: ScalersRaw.h:224
std::vector< unsigned int > gtAlgoCounts() const
unsigned int collectionTimeLumiSeg_nsec
Definition: ScalersRaw.h:94
unsigned long long triggersPhysicsLost_
unsigned int lumiSegmentOrbits
Definition: ScalersRaw.h:82
unsigned long long l1AsPhysics() const
unsigned int lumiSegmentOrbits() const
static double percentLSActive(unsigned long long counts)
unsigned long long l1AsCalibration_
unsigned int lastHardReset
Definition: ScalersRaw.h:233
unsigned long long deadtimeBeamActiveTriggerRules() const
unsigned long long deadtimeBeamActivePrivateOrbit_
unsigned int lastHardReset() const
unsigned int lastStart
Definition: ScalersRaw.h:231
unsigned int lumiSegmentNrLumiSeg_
unsigned long long l1AsTest() const
unsigned int lastOrbitCounter0() const
unsigned int lumiSegmentNr() const
unsigned int lastEventCounter0
Definition: ScalersRaw.h:232
unsigned long long deadtimeBeamActive_
unsigned long long deadtimeBeamActivePrivateOrbit
Definition: ScalersRaw.h:108
unsigned long long deadtime() const
unsigned int gtAlgoCounts[ScalersRaw::N_L1_TRIGGERS_v1]
Definition: ScalersRaw.h:112
unsigned int lastTestEnable
Definition: ScalersRaw.h:229
unsigned long long gtTriggers() const
static double percentLS(unsigned long long counts)
unsigned int eventID() const
float gtTriggersRate() const
unsigned int sourceID() const
struct timespec collectionTimeLumiSeg() const
unsigned long long deadtimeBeamActiveTimeSlot() const
unsigned long long triggersPhysicsLostBeamActive
Definition: ScalersRaw.h:98
unsigned long long l1AsPhysics_
unsigned long long gtEvents
Definition: ScalersRaw.h:87
unsigned long long deadtimeBeamActiveCalibration() const
double rate(double x)
Definition: Constants.cc:3
unsigned long long deadtime
Definition: ScalersRaw.h:104
unsigned long long gtTriggers
Definition: ScalersRaw.h:86
unsigned long long triggersPhysicsLostBeamActive_
unsigned int collectionTimeLumiSeg_sec
Definition: ScalersRaw.h:93
unsigned int lastOrbitCounter0
Definition: ScalersRaw.h:228
unsigned long long deadtimeBeamActiveCalibration
Definition: ScalersRaw.h:107
unsigned int gtResets
Definition: ScalersRaw.h:84
unsigned long long triggersPhysicsLostBeamInactive() const
unsigned long long spare2_
void set_tv_sec(long value)
Definition: TimeSpec.h:18
unsigned long long header
Definition: ScalersRaw.h:222
unsigned int lumiSegmentNrLumiSeg
Definition: ScalersRaw.h:95
unsigned long long spare0_
unsigned long long deadtimeBeamActiveCalibration_
unsigned long long deadtimeBeamActivePartitionController() const
unsigned int lumiSegmentNrLumiSeg() const
unsigned long long deadtimeBeamActive
Definition: ScalersRaw.h:105
unsigned long long triggersPhysicsGeneratedFDL
Definition: ScalersRaw.h:96
unsigned long long triggersPhysicsLost
Definition: ScalersRaw.h:97
unsigned long long triggersPhysicsGeneratedFDL_
unsigned int bunchCrossingErrors
Definition: ScalersRaw.h:85
unsigned int lumiSegmentNr
Definition: ScalersRaw.h:81
struct timespec collectionTime() const
unsigned long long l1AsPhysics
Definition: ScalersRaw.h:100
unsigned long long gtEvents_
unsigned long long l1AsRandom_