CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/DataFormats/Scalers/src/Level1TriggerScalers.cc

Go to the documentation of this file.
00001 /*
00002  *   File: DataFormats/Scalers/src/Level1TriggerScalers.cc   (W.Badgett)
00003  */
00004 
00005 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
00006 #include "DataFormats/Scalers/interface/ScalersRaw.h"
00007 
00008 #include <iostream>
00009 #include <time.h>
00010 #include <cstdio>
00011 
00012 Level1TriggerScalers::Level1TriggerScalers():
00013   version_(0),
00014   trigType_(0),
00015   eventID_(0),
00016   sourceID_(0),
00017   bunchNumber_(0),
00018   collectionTime_(0,0),
00019   lumiSegmentNr_(0),
00020   lumiSegmentOrbits_(0),
00021   orbitNr_(0),
00022   gtResets_(0),
00023   bunchCrossingErrors_(0),
00024   gtTriggers_(0),
00025   gtEvents_(0),
00026   gtTriggersRate_((float)0.0),
00027   gtEventsRate_((float)0.0),
00028   prescaleIndexAlgo_(0),
00029   prescaleIndexTech_(0),
00030   collectionTimeLumiSeg_(0,0),
00031   lumiSegmentNrLumiSeg_(0),
00032   triggersPhysicsGeneratedFDL_(0),
00033   triggersPhysicsLost_(0),
00034   triggersPhysicsLostBeamActive_(0),
00035   triggersPhysicsLostBeamInactive_(0),
00036   l1AsPhysics_(0),
00037   l1AsRandom_(0),
00038   l1AsTest_(0),
00039   l1AsCalibration_(0),
00040   deadtime_(0),
00041   deadtimeBeamActive_(0),
00042   deadtimeBeamActiveTriggerRules_(0),
00043   deadtimeBeamActiveCalibration_(0),
00044   deadtimeBeamActivePrivateOrbit_(0),
00045   deadtimeBeamActivePartitionController_(0),
00046   deadtimeBeamActiveTimeSlot_(0),
00047   gtAlgoCounts_(nLevel1Triggers),
00048   gtTechCounts_(nLevel1TestTriggers),
00049   lastOrbitCounter0_(0),
00050   lastTestEnable_(0),
00051   lastResync_(0),
00052   lastStart_(0),
00053   lastEventCounter0_(0),
00054   lastHardReset_(0),
00055   spare0_(0ULL),
00056   spare1_(0ULL),
00057   spare2_(0ULL)
00058 { 
00059 }
00060 
00061 Level1TriggerScalers::Level1TriggerScalers(const unsigned char * rawData)
00062 { 
00063   Level1TriggerScalers();
00064 
00065   struct ScalersEventRecordRaw_v5 * raw 
00066     = (struct ScalersEventRecordRaw_v5 *)rawData;
00067 
00068   trigType_     = ( raw->header >> 56 ) &        0xFULL;
00069   eventID_      = ( raw->header >> 32 ) & 0x00FFFFFFULL;
00070   sourceID_     = ( raw->header >>  8 ) & 0x00000FFFULL;
00071   bunchNumber_  = ( raw->header >> 20 ) &      0xFFFULL;
00072 
00073   version_      = raw->version;
00074   if ( version_ >= 3 )
00075   {
00076     collectionTime_.set_tv_sec( static_cast<long>(
00077       raw->trig.collectionTime_sec));
00078     collectionTime_.set_tv_nsec( 
00079       raw->trig.collectionTime_nsec);
00080 
00081     lumiSegmentNr_        = raw->trig.lumiSegmentNr;
00082     lumiSegmentOrbits_    = raw->trig.lumiSegmentOrbits;
00083     orbitNr_              = raw->trig.orbitNr;
00084     gtResets_             = raw->trig.gtResets;
00085     bunchCrossingErrors_  = raw->trig.bunchCrossingErrors;
00086     gtTriggers_           = raw->trig.gtTriggers;
00087     gtEvents_             = raw->trig.gtEvents;
00088     gtTriggersRate_       = raw->trig.gtTriggersRate;
00089     gtEventsRate_         = raw->trig.gtEventsRate;
00090     prescaleIndexAlgo_    = raw->trig.prescaleIndexAlgo;
00091     prescaleIndexTech_    = raw->trig.prescaleIndexTech;
00092 
00093     collectionTimeLumiSeg_.set_tv_sec( static_cast<long>(
00094       raw->trig.collectionTimeLumiSeg_sec));
00095     collectionTimeLumiSeg_.set_tv_nsec( 
00096       raw->trig.collectionTimeLumiSeg_nsec);
00097 
00098     lumiSegmentNrLumiSeg_           = raw->trig.lumiSegmentNrLumiSeg;
00099     triggersPhysicsGeneratedFDL_    = raw->trig.triggersPhysicsGeneratedFDL;
00100     triggersPhysicsLost_            = raw->trig.triggersPhysicsLost;
00101     triggersPhysicsLostBeamActive_  = raw->trig.triggersPhysicsLostBeamActive;
00102     triggersPhysicsLostBeamInactive_ = 
00103       raw->trig.triggersPhysicsLostBeamInactive;
00104 
00105     l1AsPhysics_                    = raw->trig.l1AsPhysics;
00106     l1AsRandom_                     = raw->trig.l1AsRandom;
00107     l1AsTest_                       = raw->trig.l1AsTest;
00108     l1AsCalibration_                = raw->trig.l1AsCalibration;
00109     deadtime_                       = raw->trig.deadtime;
00110     deadtimeBeamActive_             = raw->trig.deadtimeBeamActive;
00111     deadtimeBeamActiveTriggerRules_ = raw->trig.deadtimeBeamActiveTriggerRules;
00112     deadtimeBeamActiveCalibration_  = raw->trig.deadtimeBeamActiveCalibration;
00113     deadtimeBeamActivePrivateOrbit_ = raw->trig.deadtimeBeamActivePrivateOrbit;
00114     deadtimeBeamActivePartitionController_ = 
00115       raw->trig.deadtimeBeamActivePartitionController;
00116     deadtimeBeamActiveTimeSlot_ = raw->trig.deadtimeBeamActiveTimeSlot;
00117 
00118     for ( int i=0; i<ScalersRaw::N_L1_TRIGGERS_v1; i++)
00119     { gtAlgoCounts_.push_back( raw->trig.gtAlgoCounts[i]);}
00120 
00121     for ( int i=0; i<ScalersRaw::N_L1_TEST_TRIGGERS_v1; i++)
00122     { gtTechCounts_.push_back( raw->trig.gtTechCounts[i]);}
00123 
00124     if ( version_ >= 5 )
00125     {
00126       lastOrbitCounter0_ = raw->lastOrbitCounter0;
00127       lastTestEnable_    = raw->lastTestEnable;
00128       lastResync_        = raw->lastResync;
00129       lastStart_         = raw->lastStart;
00130       lastEventCounter0_ = raw->lastEventCounter0;
00131       lastHardReset_     = raw->lastHardReset;
00132       spare0_            = raw->spare[0];
00133       spare1_            = raw->spare[1];
00134       spare2_            = raw->spare[2];
00135     }
00136     else
00137     {
00138       lastOrbitCounter0_ = 0UL;
00139       lastTestEnable_    = 0UL;
00140       lastResync_        = 0UL;
00141       lastStart_         = 0UL;
00142       lastEventCounter0_ = 0UL;
00143       lastHardReset_     = 0UL;
00144       spare0_            = 0ULL;
00145       spare1_            = 0ULL;
00146       spare2_            = 0ULL;
00147     }
00148   }
00149 }
00150 
00151 Level1TriggerScalers::~Level1TriggerScalers() { } 
00152 
00153 double Level1TriggerScalers::rateLS(unsigned int counts)
00154 { return(rateLS(counts,firstShortLSRun));}
00155 
00156 double Level1TriggerScalers::rateLS(unsigned long long counts)
00157 { return(rateLS(counts,firstShortLSRun));}
00158 
00159 double Level1TriggerScalers::rateLS(unsigned int counts,
00160                                     int runNumber)
00161 { 
00162   unsigned long long counts64 = (unsigned long long)counts;
00163   return(rateLS(counts64,runNumber));
00164 }
00165 
00166 double Level1TriggerScalers::rateLS(unsigned long long counts,
00167                                     int runNumber)
00168 { 
00169   double rate;
00170   if (( runNumber >= firstShortLSRun ) || ( runNumber <= 1 ))
00171   {
00172     rate = ((double)counts) / 23.31040958083832;
00173   }
00174   else
00175   {
00176     rate = ((double)counts) / 93.24163832335329;
00177   }
00178   return(rate);
00179 }
00180 
00181 double Level1TriggerScalers::percentLS(unsigned long long counts)
00182 { return(percentLS(counts,firstShortLSRun));}
00183 
00184 double Level1TriggerScalers::percentLS(unsigned long long counts,
00185                                        int runNumber)
00186 { 
00187   double percent;
00188   if (( runNumber >= firstShortLSRun ) || ( runNumber <= 1 ))
00189   {
00190     percent = ((double)counts) /  9342812.16;
00191   }
00192   else
00193   {
00194     percent = ((double)counts) / 37371248.64;
00195   }
00196   if ( percent > 100.0000 ) { percent = 100.0;}
00197   return(percent);
00198 }
00199 
00200 double Level1TriggerScalers::percentLSActive(unsigned long long counts)
00201 { return(percentLSActive(counts,firstShortLSRun));}
00202 
00203 double Level1TriggerScalers::percentLSActive(unsigned long long counts,
00204                                        int runNumber)
00205 { 
00206   double percent;
00207   if (( runNumber >= firstShortLSRun ) || ( runNumber <= 1 ))
00208   {
00209     percent = ((double)counts) /  7361003.52;
00210   }
00211   else
00212   {
00213     percent = ((double)counts) / 29444014.08;
00214   }
00215   if ( percent > 100.0000 ) { percent = 100.0;}
00216   return(percent);
00217 }
00218 
00220 std::ostream& operator<<(std::ostream& s,Level1TriggerScalers const &c) 
00221 {
00222   s << "Level1TriggerScalers    Version:" << c.version() <<
00223     "   SourceID: " << c.sourceID() << std::endl;
00224   char line[128];
00225   char zeitHeaven[128];
00226   struct tm * horaHeaven;
00227 
00228   sprintf(line, " TrigType: %d   EventID: %d    BunchNumber: %d", 
00229           c.trigType(), c.eventID(), c.bunchNumber());
00230   s << line << std::endl;
00231 
00232   struct timespec secondsToHeaven = c.collectionTime();
00233   horaHeaven = gmtime(&secondsToHeaven.tv_sec);
00234   strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
00235   sprintf(line, " CollectionTime:        %s.%9.9d" , 
00236           zeitHeaven, (int)secondsToHeaven.tv_nsec);
00237   s << line << std::endl;
00238 
00239   sprintf(line,
00240           " LumiSegmentNr:        %10u   LumiSegmentOrbits:     %10u",
00241           c.lumiSegmentNr(), c.lumiSegmentOrbits());
00242   s << line << std::endl;
00243 
00244   sprintf(line,
00245           " LumiSegmentNrLumiSeg: %10u   OrbitNr:               %10u ",
00246           c.lumiSegmentNrLumiSeg(),  c.orbitNr());
00247   s << line << std::endl;
00248 
00249   sprintf(line,
00250           " GtResets:             %10u   BunchCrossingErrors:   %10u",
00251           c.gtResets(), c.bunchCrossingErrors());
00252   s << line << std::endl;
00253 
00254   sprintf(line,
00255           " PrescaleIndexAlgo:    %10d   PrescaleIndexTech:     %10d",
00256           c.prescaleIndexAlgo(), c.prescaleIndexTech());
00257   s << line << std::endl;
00258 
00259   sprintf(line, " GtTriggers:                      %20llu %22.3f Hz", 
00260           c.gtTriggers(), c.gtTriggersRate());
00261   s << line << std::endl;
00262 
00263   sprintf(line, " GtEvents:                        %20llu %22.3f Hz", 
00264           c.gtEvents(), c.gtEventsRate());
00265   s << line << std::endl;
00266 
00267   secondsToHeaven = c.collectionTimeLumiSeg();
00268   horaHeaven = gmtime(&secondsToHeaven.tv_sec);
00269   strftime(zeitHeaven, sizeof(zeitHeaven), "%Y.%m.%d %H:%M:%S", horaHeaven);
00270   sprintf(line, " CollectionTimeLumiSeg: %s.%9.9d" , 
00271           zeitHeaven, (int)secondsToHeaven.tv_nsec);
00272   s << line << std::endl;
00273 
00274 
00275   sprintf(line, " TriggersPhysicsGeneratedFDL:     %20llu %22.3f Hz",
00276           c.triggersPhysicsGeneratedFDL(),
00277           Level1TriggerScalers::rateLS(c.triggersPhysicsGeneratedFDL()));
00278   s << line << std::endl;
00279 
00280   sprintf(line, " TriggersPhysicsLost:             %20llu %22.3f Hz",
00281           c.triggersPhysicsLost(),
00282           Level1TriggerScalers::rateLS(c.triggersPhysicsLost()));
00283   s << line << std::endl;
00284 
00285   sprintf(line, " TriggersPhysicsLostBeamActive:   %20llu %22.3f Hz",
00286           c.triggersPhysicsLostBeamActive(),
00287           Level1TriggerScalers::rateLS(c.triggersPhysicsLostBeamActive()));
00288   s << line << std::endl;
00289 
00290   sprintf(line, " TriggersPhysicsLostBeamInactive: %20llu %22.3f Hz",
00291           c.triggersPhysicsLostBeamInactive(),
00292           Level1TriggerScalers::rateLS(c.triggersPhysicsLostBeamInactive()));
00293   s << line << std::endl;
00294 
00295   sprintf(line, " L1AsPhysics:                     %20llu %22.3f Hz",
00296           c.l1AsPhysics(),
00297           Level1TriggerScalers::rateLS(c.l1AsPhysics()));
00298   s << line << std::endl;
00299 
00300   sprintf(line, " L1AsRandom:                      %20llu %22.3f Hz",
00301           c.l1AsRandom(),
00302           Level1TriggerScalers::rateLS(c.l1AsRandom()));
00303   s << line << std::endl;
00304 
00305   sprintf(line, " L1AsTest:                        %20llu %22.3f Hz",
00306           c.l1AsTest(),
00307           Level1TriggerScalers::rateLS(c.l1AsTest()));
00308   s << line << std::endl;
00309 
00310   sprintf(line, " L1AsCalibration:                 %20llu %22.3f Hz",
00311           c.l1AsCalibration(),
00312           Level1TriggerScalers::rateLS(c.l1AsCalibration()));
00313   s << line << std::endl;
00314 
00315   sprintf(line, " Deadtime:                             %20llu %17.3f%%",
00316           c.deadtime(),
00317           Level1TriggerScalers::percentLS(c.deadtime()));
00318   s << line << std::endl;
00319 
00320   sprintf(line, " DeadtimeBeamActive:                   %20llu %17.3f%%",
00321           c.deadtimeBeamActive(),
00322           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActive()));
00323   s << line << std::endl;
00324 
00325   sprintf(line, " DeadtimeBeamActiveTriggerRules:       %20llu %17.3f%%",
00326           c.deadtimeBeamActiveTriggerRules(),
00327           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveTriggerRules()));
00328   s << line << std::endl;
00329 
00330   sprintf(line, " DeadtimeBeamActiveCalibration:        %20llu %17.3f%%",
00331           c.deadtimeBeamActiveCalibration(),
00332           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveCalibration()));
00333   s << line << std::endl;
00334 
00335   sprintf(line, " DeadtimeBeamActivePrivateOrbit:       %20llu %17.3f%%",
00336           c.deadtimeBeamActivePrivateOrbit(),
00337           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActivePrivateOrbit()));
00338   s << line << std::endl;
00339 
00340   sprintf(line, " DeadtimeBeamActivePartitionController:%20llu %17.3f%%",
00341           c.deadtimeBeamActivePartitionController(),
00342           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActivePartitionController()));
00343   s << line << std::endl;
00344 
00345   sprintf(line, " DeadtimeBeamActiveTimeSlot:           %20llu %17.3f%%",
00346           c.deadtimeBeamActiveTimeSlot(),
00347           Level1TriggerScalers::percentLSActive(c.deadtimeBeamActiveTimeSlot()));
00348   s << line << std::endl;
00349 
00350   s << "Physics GtAlgoCounts" << std::endl;
00351   const std::vector<unsigned int> gtAlgoCounts = c.gtAlgoCounts();
00352   int length = gtAlgoCounts.size() / 4;
00353   for ( int i=0; i<length; i++)
00354   {
00355     sprintf(line," %3.3d: %10u    %3.3d: %10u    %3.3d: %10u    %3.3d: %10u",
00356             i,              gtAlgoCounts[i], 
00357             (i+length),     gtAlgoCounts[i+length], 
00358             (i+(length*2)), gtAlgoCounts[i+(length*2)], 
00359             (i+(length*3)), gtAlgoCounts[i+(length*3)]);
00360     s << line << std::endl;
00361   }
00362 
00363   s << "Test GtTechCounts" << std::endl;
00364   const std::vector<unsigned int> gtTechCounts = c.gtTechCounts();
00365   length = gtTechCounts.size() / 4;
00366   for ( int i=0; i<length; i++)
00367   {
00368     sprintf(line," %3.3d: %10u    %3.3d: %10u    %3.3d: %10u    %3.3d: %10u",
00369             i,              gtTechCounts[i], 
00370             (i+length),     gtTechCounts[i+length], 
00371             (i+(length*2)), gtTechCounts[i+(length*2)], 
00372             (i+(length*3)), gtTechCounts[i+(length*3)]);
00373     s << line << std::endl;
00374   }
00375 
00376   if ( c.version() >= 5 )
00377   {
00378     sprintf(line," LastOrbitCounter0:  %10u  0x%8.8X", c.lastOrbitCounter0(),
00379             c.lastOrbitCounter0()); 
00380     s << line << std::endl;
00381     
00382     sprintf(line," LastTestEnable:     %10u  0x%8.8X", c.lastTestEnable(),
00383             c.lastTestEnable()); 
00384     s << line << std::endl;
00385     
00386     sprintf(line," LastResync:         %10u  0x%8.8X", c.lastResync(),
00387             c.lastResync()); 
00388     s << line << std::endl;
00389     
00390     sprintf(line," LastStart:          %10u  0x%8.8X", c.lastStart(),
00391             c.lastStart()); 
00392     s << line << std::endl;
00393     
00394     sprintf(line," LastEventCounter0:  %10u  0x%8.8X", c.lastEventCounter0(),
00395             c.lastEventCounter0()); 
00396     s << line << std::endl;
00397     
00398     sprintf(line," LastHardReset:      %10u  0x%8.8X", c.lastHardReset(),
00399             c.lastHardReset()); 
00400     s << line << std::endl;
00401   }
00402 
00403   return s;
00404 }