Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:05:09

0001 /*
0002  *   File: DataFormats/Scalers/src/Level1TriggerScalers.cc   (W.Badgett)
0003  */
0004 
0005 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
0006 #include "DataFormats/Scalers/interface/ScalersRaw.h"
0007 
0008 #include <iostream>
0009 #include <ctime>
0010 #include <cstdio>
0011 
0012 Level1TriggerScalers::Level1TriggerScalers()
0013     : version_(0),
0014       trigType_(0),
0015       eventID_(0),
0016       sourceID_(0),
0017       bunchNumber_(0),
0018       collectionTime_(0, 0),
0019       lumiSegmentNr_(0),
0020       lumiSegmentOrbits_(0),
0021       orbitNr_(0),
0022       gtResets_(0),
0023       bunchCrossingErrors_(0),
0024       gtTriggers_(0),
0025       gtEvents_(0),
0026       gtTriggersRate_((float)0.0),
0027       gtEventsRate_((float)0.0),
0028       prescaleIndexAlgo_(0),
0029       prescaleIndexTech_(0),
0030       collectionTimeLumiSeg_(0, 0),
0031       lumiSegmentNrLumiSeg_(0),
0032       triggersPhysicsGeneratedFDL_(0),
0033       triggersPhysicsLost_(0),
0034       triggersPhysicsLostBeamActive_(0),
0035       triggersPhysicsLostBeamInactive_(0),
0036       l1AsPhysics_(0),
0037       l1AsRandom_(0),
0038       l1AsTest_(0),
0039       l1AsCalibration_(0),
0040       deadtime_(0),
0041       deadtimeBeamActive_(0),
0042       deadtimeBeamActiveTriggerRules_(0),
0043       deadtimeBeamActiveCalibration_(0),
0044       deadtimeBeamActivePrivateOrbit_(0),
0045       deadtimeBeamActivePartitionController_(0),
0046       deadtimeBeamActiveTimeSlot_(0),
0047       gtAlgoCounts_(nLevel1Triggers),
0048       gtTechCounts_(nLevel1TestTriggers),
0049       lastOrbitCounter0_(0),
0050       lastTestEnable_(0),
0051       lastResync_(0),
0052       lastStart_(0),
0053       lastEventCounter0_(0),
0054       lastHardReset_(0),
0055       spare0_(0ULL),
0056       spare1_(0ULL),
0057       spare2_(0ULL) {}
0058 
0059 Level1TriggerScalers::Level1TriggerScalers(const unsigned char* rawData) {
0060   Level1TriggerScalers();
0061 
0062   struct ScalersEventRecordRaw_v5 const* raw = reinterpret_cast<struct ScalersEventRecordRaw_v5 const*>(rawData);
0063 
0064   trigType_ = (raw->header >> 56) & 0xFULL;
0065   eventID_ = (raw->header >> 32) & 0x00FFFFFFULL;
0066   sourceID_ = (raw->header >> 8) & 0x00000FFFULL;
0067   bunchNumber_ = (raw->header >> 20) & 0xFFFULL;
0068 
0069   version_ = raw->version;
0070   if (version_ >= 3) {
0071     collectionTime_.set_tv_sec(static_cast<long>(raw->trig.collectionTime_sec));
0072     collectionTime_.set_tv_nsec(raw->trig.collectionTime_nsec);
0073 
0074     lumiSegmentNr_ = raw->trig.lumiSegmentNr;
0075     lumiSegmentOrbits_ = raw->trig.lumiSegmentOrbits;
0076     orbitNr_ = raw->trig.orbitNr;
0077     gtResets_ = raw->trig.gtResets;
0078     bunchCrossingErrors_ = raw->trig.bunchCrossingErrors;
0079     gtTriggers_ = raw->trig.gtTriggers;
0080     gtEvents_ = raw->trig.gtEvents;
0081     gtTriggersRate_ = raw->trig.gtTriggersRate;
0082     gtEventsRate_ = raw->trig.gtEventsRate;
0083     prescaleIndexAlgo_ = raw->trig.prescaleIndexAlgo;
0084     prescaleIndexTech_ = raw->trig.prescaleIndexTech;
0085 
0086     collectionTimeLumiSeg_.set_tv_sec(static_cast<long>(raw->trig.collectionTimeLumiSeg_sec));
0087     collectionTimeLumiSeg_.set_tv_nsec(raw->trig.collectionTimeLumiSeg_nsec);
0088 
0089     lumiSegmentNrLumiSeg_ = raw->trig.lumiSegmentNrLumiSeg;
0090     triggersPhysicsGeneratedFDL_ = raw->trig.triggersPhysicsGeneratedFDL;
0091     triggersPhysicsLost_ = raw->trig.triggersPhysicsLost;
0092     triggersPhysicsLostBeamActive_ = raw->trig.triggersPhysicsLostBeamActive;
0093     triggersPhysicsLostBeamInactive_ = raw->trig.triggersPhysicsLostBeamInactive;
0094 
0095     l1AsPhysics_ = raw->trig.l1AsPhysics;
0096     l1AsRandom_ = raw->trig.l1AsRandom;
0097     l1AsTest_ = raw->trig.l1AsTest;
0098     l1AsCalibration_ = raw->trig.l1AsCalibration;
0099     deadtime_ = raw->trig.deadtime;
0100     deadtimeBeamActive_ = raw->trig.deadtimeBeamActive;
0101     deadtimeBeamActiveTriggerRules_ = raw->trig.deadtimeBeamActiveTriggerRules;
0102     deadtimeBeamActiveCalibration_ = raw->trig.deadtimeBeamActiveCalibration;
0103     deadtimeBeamActivePrivateOrbit_ = raw->trig.deadtimeBeamActivePrivateOrbit;
0104     deadtimeBeamActivePartitionController_ = raw->trig.deadtimeBeamActivePartitionController;
0105     deadtimeBeamActiveTimeSlot_ = raw->trig.deadtimeBeamActiveTimeSlot;
0106 
0107     for (int i = 0; i < ScalersRaw::N_L1_TRIGGERS_v1; i++) {
0108       gtAlgoCounts_.push_back(raw->trig.gtAlgoCounts[i]);
0109     }
0110 
0111     for (int i = 0; i < ScalersRaw::N_L1_TEST_TRIGGERS_v1; i++) {
0112       gtTechCounts_.push_back(raw->trig.gtTechCounts[i]);
0113     }
0114 
0115     if (version_ >= 5) {
0116       lastOrbitCounter0_ = raw->lastOrbitCounter0;
0117       lastTestEnable_ = raw->lastTestEnable;
0118       lastResync_ = raw->lastResync;
0119       lastStart_ = raw->lastStart;
0120       lastEventCounter0_ = raw->lastEventCounter0;
0121       lastHardReset_ = raw->lastHardReset;
0122       spare0_ = raw->spare[0];
0123       spare1_ = raw->spare[1];
0124       spare2_ = raw->spare[2];
0125     } else {
0126       lastOrbitCounter0_ = 0UL;
0127       lastTestEnable_ = 0UL;
0128       lastResync_ = 0UL;
0129       lastStart_ = 0UL;
0130       lastEventCounter0_ = 0UL;
0131       lastHardReset_ = 0UL;
0132       spare0_ = 0ULL;
0133       spare1_ = 0ULL;
0134       spare2_ = 0ULL;
0135     }
0136   }
0137 }
0138 
0139 Level1TriggerScalers::~Level1TriggerScalers() {}
0140 
0141 double Level1TriggerScalers::rateLS(unsigned int counts) { return (rateLS(counts, firstShortLSRun)); }
0142 
0143 double Level1TriggerScalers::rateLS(unsigned long long counts) { return (rateLS(counts, firstShortLSRun)); }
0144 
0145 double Level1TriggerScalers::rateLS(unsigned int counts, int runNumber) {
0146   unsigned long long counts64 = (unsigned long long)counts;
0147   return (rateLS(counts64, runNumber));
0148 }
0149 
0150 double Level1TriggerScalers::rateLS(unsigned long long counts, int runNumber) {
0151   double rate;
0152   if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
0153     rate = ((double)counts) / 23.31040958083832;
0154   } else {
0155     rate = ((double)counts) / 93.24163832335329;
0156   }
0157   return (rate);
0158 }
0159 
0160 double Level1TriggerScalers::percentLS(unsigned long long counts) { return (percentLS(counts, firstShortLSRun)); }
0161 
0162 double Level1TriggerScalers::percentLS(unsigned long long counts, int runNumber) {
0163   double percent;
0164   if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
0165     percent = ((double)counts) / 9342812.16;
0166   } else {
0167     percent = ((double)counts) / 37371248.64;
0168   }
0169   if (percent > 100.0000) {
0170     percent = 100.0;
0171   }
0172   return (percent);
0173 }
0174 
0175 double Level1TriggerScalers::percentLSActive(unsigned long long counts) {
0176   return (percentLSActive(counts, firstShortLSRun));
0177 }
0178 
0179 double Level1TriggerScalers::percentLSActive(unsigned long long counts, int runNumber) {
0180   double percent;
0181   if ((runNumber >= firstShortLSRun) || (runNumber <= 1)) {
0182     percent = ((double)counts) / 7361003.52;
0183   } else {
0184     percent = ((double)counts) / 29444014.08;
0185   }
0186   if (percent > 100.0000) {
0187     percent = 100.0;
0188   }
0189   return (percent);
0190 }
0191 
0192 /// Pretty-print operator for Level1TriggerScalers
0193 std::ostream& operator<<(std::ostream& s, Level1TriggerScalers const& c) {
0194   s << "Level1TriggerScalers    Version:" << c.version() << "   SourceID: " << c.sourceID() << std::endl;
0195   constexpr size_t kLineBufferSize = 164;
0196   char line[kLineBufferSize];
0197   char zeitHeaven[128];
0198   struct tm* horaHeaven;
0199 
0200   snprintf(line,
0201            kLineBufferSize,
0202            " TrigType: %d   EventID: %d    BunchNumber: %d",
0203            c.trigType(),
0204            c.eventID(),
0205            c.bunchNumber());
0206   s << line << std::endl;
0207 
0208   struct timespec secondsToHeaven = c.collectionTime();
0209   horaHeaven = gmtime(&secondsToHeaven.tv_sec);
0210   strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
0211   snprintf(line, kLineBufferSize, " CollectionTime:        %s.%9.9d", zeitHeaven, (int)secondsToHeaven.tv_nsec);
0212   s << line << std::endl;
0213 
0214   snprintf(line,
0215            kLineBufferSize,
0216            " LumiSegmentNr:        %10u   LumiSegmentOrbits:     %10u",
0217            c.lumiSegmentNr(),
0218            c.lumiSegmentOrbits());
0219   s << line << std::endl;
0220 
0221   snprintf(line,
0222            kLineBufferSize,
0223            " LumiSegmentNrLumiSeg: %10u   OrbitNr:               %10u ",
0224            c.lumiSegmentNrLumiSeg(),
0225            c.orbitNr());
0226   s << line << std::endl;
0227 
0228   snprintf(line,
0229            kLineBufferSize,
0230            " GtResets:             %10u   BunchCrossingErrors:   %10u",
0231            c.gtResets(),
0232            c.bunchCrossingErrors());
0233   s << line << std::endl;
0234 
0235   snprintf(line,
0236            kLineBufferSize,
0237            " PrescaleIndexAlgo:    %10d   PrescaleIndexTech:     %10d",
0238            c.prescaleIndexAlgo(),
0239            c.prescaleIndexTech());
0240   s << line << std::endl;
0241 
0242   snprintf(
0243       line, kLineBufferSize, " GtTriggers:                      %20llu %22.3f Hz", c.gtTriggers(), c.gtTriggersRate());
0244   s << line << std::endl;
0245 
0246   snprintf(line, kLineBufferSize, " GtEvents:                        %20llu %22.3f Hz", c.gtEvents(), c.gtEventsRate());
0247   s << line << std::endl;
0248 
0249   secondsToHeaven = c.collectionTimeLumiSeg();
0250   horaHeaven = gmtime(&secondsToHeaven.tv_sec);
0251   strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
0252   snprintf(line, kLineBufferSize, " CollectionTimeLumiSeg: %s.%9.9d", zeitHeaven, (int)secondsToHeaven.tv_nsec);
0253   s << line << std::endl;
0254 
0255   snprintf(line,
0256            kLineBufferSize,
0257            " TriggersPhysicsGeneratedFDL:     %20llu %22.3f Hz",
0258            c.triggersPhysicsGeneratedFDL(),
0259            Level1TriggerScalers::rateLS(c.triggersPhysicsGeneratedFDL()));
0260   s << line << std::endl;
0261 
0262   snprintf(line,
0263            kLineBufferSize,
0264            " TriggersPhysicsLost:             %20llu %22.3f Hz",
0265            c.triggersPhysicsLost(),
0266            Level1TriggerScalers::rateLS(c.triggersPhysicsLost()));
0267   s << line << std::endl;
0268 
0269   snprintf(line,
0270            kLineBufferSize,
0271            " TriggersPhysicsLostBeamActive:   %20llu %22.3f Hz",
0272            c.triggersPhysicsLostBeamActive(),
0273            Level1TriggerScalers::rateLS(c.triggersPhysicsLostBeamActive()));
0274   s << line << std::endl;
0275 
0276   snprintf(line,
0277            kLineBufferSize,
0278            " TriggersPhysicsLostBeamInactive: %20llu %22.3f Hz",
0279            c.triggersPhysicsLostBeamInactive(),
0280            Level1TriggerScalers::rateLS(c.triggersPhysicsLostBeamInactive()));
0281   s << line << std::endl;
0282 
0283   snprintf(line,
0284            kLineBufferSize,
0285            " L1AsPhysics:                     %20llu %22.3f Hz",
0286            c.l1AsPhysics(),
0287            Level1TriggerScalers::rateLS(c.l1AsPhysics()));
0288   s << line << std::endl;
0289 
0290   snprintf(line,
0291            kLineBufferSize,
0292            " L1AsRandom:                      %20llu %22.3f Hz",
0293            c.l1AsRandom(),
0294            Level1TriggerScalers::rateLS(c.l1AsRandom()));
0295   s << line << std::endl;
0296 
0297   snprintf(line,
0298            kLineBufferSize,
0299            " L1AsTest:                        %20llu %22.3f Hz",
0300            c.l1AsTest(),
0301            Level1TriggerScalers::rateLS(c.l1AsTest()));
0302   s << line << std::endl;
0303 
0304   snprintf(line,
0305            kLineBufferSize,
0306            " L1AsCalibration:                 %20llu %22.3f Hz",
0307            c.l1AsCalibration(),
0308            Level1TriggerScalers::rateLS(c.l1AsCalibration()));
0309   s << line << std::endl;
0310 
0311   snprintf(line,
0312            kLineBufferSize,
0313            " Deadtime:                             %20llu %17.3f%%",
0314            c.deadtime(),
0315            Level1TriggerScalers::percentLS(c.deadtime()));
0316   s << line << std::endl;
0317 
0318   snprintf(line,
0319            kLineBufferSize,
0320            " DeadtimeBeamActive:                   %20llu %17.3f%%",
0321            c.deadtimeBeamActive(),
0322            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActive()));
0323   s << line << std::endl;
0324 
0325   snprintf(line,
0326            kLineBufferSize,
0327            " DeadtimeBeamActiveTriggerRules:       %20llu %17.3f%%",
0328            c.deadtimeBeamActiveTriggerRules(),
0329            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveTriggerRules()));
0330   s << line << std::endl;
0331 
0332   snprintf(line,
0333            kLineBufferSize,
0334            " DeadtimeBeamActiveCalibration:        %20llu %17.3f%%",
0335            c.deadtimeBeamActiveCalibration(),
0336            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveCalibration()));
0337   s << line << std::endl;
0338 
0339   snprintf(line,
0340            kLineBufferSize,
0341            " DeadtimeBeamActivePrivateOrbit:       %20llu %17.3f%%",
0342            c.deadtimeBeamActivePrivateOrbit(),
0343            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActivePrivateOrbit()));
0344   s << line << std::endl;
0345 
0346   snprintf(line,
0347            kLineBufferSize,
0348            " DeadtimeBeamActivePartitionController:%20llu %17.3f%%",
0349            c.deadtimeBeamActivePartitionController(),
0350            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActivePartitionController()));
0351   s << line << std::endl;
0352 
0353   snprintf(line,
0354            kLineBufferSize,
0355            " DeadtimeBeamActiveTimeSlot:           %20llu %17.3f%%",
0356            c.deadtimeBeamActiveTimeSlot(),
0357            Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveTimeSlot()));
0358   s << line << std::endl;
0359 
0360   s << "Physics GtAlgoCounts" << std::endl;
0361   const std::vector<unsigned int> gtAlgoCounts = c.gtAlgoCounts();
0362   int length = gtAlgoCounts.size() / 4;
0363   for (int i = 0; i < length; i++) {
0364     snprintf(line,
0365              kLineBufferSize,
0366              " %3.3d: %10u    %3.3d: %10u    %3.3d: %10u    %3.3d: %10u",
0367              i,
0368              gtAlgoCounts[i],
0369              (i + length),
0370              gtAlgoCounts[i + length],
0371              (i + (length * 2)),
0372              gtAlgoCounts[i + (length * 2)],
0373              (i + (length * 3)),
0374              gtAlgoCounts[i + (length * 3)]);
0375     s << line << std::endl;
0376   }
0377 
0378   s << "Test GtTechCounts" << std::endl;
0379   const std::vector<unsigned int> gtTechCounts = c.gtTechCounts();
0380   length = gtTechCounts.size() / 4;
0381   for (int i = 0; i < length; i++) {
0382     snprintf(line,
0383              kLineBufferSize,
0384              " %3.3d: %10u    %3.3d: %10u    %3.3d: %10u    %3.3d: %10u",
0385              i,
0386              gtTechCounts[i],
0387              (i + length),
0388              gtTechCounts[i + length],
0389              (i + (length * 2)),
0390              gtTechCounts[i + (length * 2)],
0391              (i + (length * 3)),
0392              gtTechCounts[i + (length * 3)]);
0393     s << line << std::endl;
0394   }
0395 
0396   if (c.version() >= 5) {
0397     snprintf(line, kLineBufferSize, " LastOrbitCounter0:  %10u  0x%8.8X", c.lastOrbitCounter0(), c.lastOrbitCounter0());
0398     s << line << std::endl;
0399 
0400     snprintf(line, kLineBufferSize, " LastTestEnable:     %10u  0x%8.8X", c.lastTestEnable(), c.lastTestEnable());
0401     s << line << std::endl;
0402 
0403     snprintf(line, kLineBufferSize, " LastResync:         %10u  0x%8.8X", c.lastResync(), c.lastResync());
0404     s << line << std::endl;
0405 
0406     snprintf(line, kLineBufferSize, " LastStart:          %10u  0x%8.8X", c.lastStart(), c.lastStart());
0407     s << line << std::endl;
0408 
0409     snprintf(line, kLineBufferSize, " LastEventCounter0:  %10u  0x%8.8X", c.lastEventCounter0(), c.lastEventCounter0());
0410     s << line << std::endl;
0411 
0412     snprintf(line, kLineBufferSize, " LastHardReset:      %10u  0x%8.8X", c.lastHardReset(), c.lastHardReset());
0413     s << line << std::endl;
0414   }
0415 
0416   return s;
0417 }