#include #include #include #include #include #include #include #include #include #include #include #include #include #include using namespace std; namespace RAT { DQTriggerProc::DQTriggerProc() : DataQualityProc( "dqtriggerproc") { } DQTriggerProc::~DQTriggerProc() { } void DQTriggerProc::BeginOfRun( DS::Run& run ) { DataQualityProc::BeginOfRun( run ); // Set criteria as class members and save to RATDB via JSON object fDQChecks = DB::Get()->GetLink( "DQCHECKS", "neutrino" ); fMinN100LRate = fDQChecks->GetD( "min_nhit100l_rate" ); fCriteria[ "min_nhit100l_rate" ] = fMinN100LRate; fMinEsumHRate = fDQChecks->GetD( "min_esum_hi_rate" ); fCriteria[ "min_esum_hi_rate" ] = fMinEsumHRate; fMissingGTIDLimit = fDQChecks->GetI("max_num_missing_gtids"); fCriteria[ "max_num_missing_gtids" ] = fMissingGTIDLimit; fBitFlipGTIDLimit = fDQChecks->GetI("max_num_bitflip_gtids"); fCriteria[ "max_num_bitflip_gtids" ] = fBitFlipGTIDLimit; firstEventFlag = true; passedSecondEvent = false; fNTrigTypeBits = 27; // Both vectors should contain one entry per bit fTrigTypeCounts.resize( fNTrigTypeBits ); fTrigTypeRates.resize( fNTrigTypeBits ); //Set the orphan count to 0 fNOrphans = 0; fBitFlipGTIDCount = 0; fMissingGTIDCount = 0; } Processor::Result DQTriggerProc::DSEvent(DS::Run&, DS::Entry& ds) { for( size_t iEV = 0; iEV < ds.GetEVCount(); iEV++ ) Event(ds, ds.GetEV(iEV)); return OK; } Processor::Result DQTriggerProc::Event( DS::Entry&, DS::EV& ev ) { int trigType = ev.GetTrigType(); vector::iterator iTrigTypeCount = fTrigTypeCounts.begin(); int bit = 0; while( iTrigTypeCount != fTrigTypeCounts.end() ) { int query = 1 << bit; if( ( trigType & query ) == query ) { *iTrigTypeCount = *iTrigTypeCount + 1; // increment count } iTrigTypeCount++; // move iterator to next bit bit++; // next bit } int gtid = ev.GetGTID(); if(gtid == 0 && fPreviousGTID!=0xfefffe){ fNOrphans++; } if(firstEventFlag){ fPreviousGTID = gtid; firstEventFlag = false; } //If the clock is at 0xfefffe the next GTID will be 0 due to clock rollover else if(gtid==0 && fPreviousGTID==0xfefffe){ } //If orphan dont do any checks else if(gtid==0){ //If we are passing the second event if(!firstEventFlag && !passedSecondEvent){ passedSecondEvent = true; } } int gtidOffset = gtid-fPreviousGTID; //GTIDS should always increase if negative possible bit flip on sign bit flag as bitflip gtid if(gtid < 0){ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } //If clock rollover or orphan dont do any bitflip testing else if(gtid==0){ } else if(fPreviousPreviousGTID==0xfefffe && fPreviousGTID==0){ } else{ //If we have a negative offset check if previous event was a bitflip if(gtidOffset < 0){ //If there are no previous bitflip GTIDS this must be a bitflip GTID if (fBitFlipGTIDs.size() == 0){ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } //If previous event was a bit flip make sure that there is an offset of two from the event before the bitflip event else if(fPreviousGTID == fBitFlipGTIDs.back()){ if(passedSecondEvent){ int previousPreviousOffset = gtid - fPreviousPreviousGTID; //If the offset is not equal to two it is likely this event was a bitflip as well if(previousPreviousOffset != 1){ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } } //If it is the second event we cant do the check and we should flag the gtid else{ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } } //If the previous event is not a bitflip we have a bitflip from 1 to 0 in the GTID else{ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } } else if(gtidOffset > 1){ //If previous event is an orphan check the GTID before the orphan if(fPreviousGTID == 0){ gtidOffset = gtid - fPreviousPreviousGTID; } //If offset is two missing and missing GTID does not end with 0xFFFF //(These are expected to be missing) add in unexpected missing GTID if(gtidOffset == 2 && (((fPreviousGTID+1) & 0xFFFF) != 0xFFFF)){ fMissingGTIDCount++; fMissingGTIDs.push_back(fPreviousGTID+1); } //Check if the offset is equal to a bitflip of 0 to 1 on any given bit (check if offset is greater than 1 incase the previous event was an orphan) else if(gtidOffset>1){ int bitNum = 2; for(int i=0; i<29; i++){ //Shift bits in into to the right once to check next bit (Do this before to prevent going over INT_MAX) bitNum = bitNum << 1; //+1 for the iteration bitNum is the bit flipped if(gtidOffset == bitNum+1){ fBitFlipGTIDCount++; fBitFlipGTIDs.push_back(gtid); } } } } } //If we are passing the second event if(!firstEventFlag && !passedSecondEvent){ passedSecondEvent = true; } //Add in the previous gtids fPreviousPreviousGTID = fPreviousGTID; fPreviousGTID = gtid; } void DQTriggerProc::EndOfRun( DS::Run& run ) { // Book histogram fTrigTypeRateHist = new TH1D( "fTrigTypeRateHist", "Rate for each trigger channel", fNTrigTypeBits, 0, fNTrigTypeBits ); const DU::TrigBits trigBits = RAT::DU::Utility::Get()->GetTrigBits(); for( int nBin = 1; nBin <= fTrigTypeRateHist->GetNbinsX(); nBin++ ) fTrigTypeRateHist->GetXaxis()->SetBinLabel ( nBin, trigBits.GetBitName( nBin - 1 ).c_str() ); fTrigTypeRateHist->SetYTitle("Rate (Hz)"); fTrigTypeRateHist->SetStats(0); // Fill rate vector vector::iterator iTrigTypeCount = fTrigTypeCounts.begin(); vector::iterator iTrigTypeRate = fTrigTypeRates.begin(); int bit = 0; debug << "DQTriggerProc::EndOfRun: trigType rates are: " << newline; debug << " "; while( iTrigTypeCount != fTrigTypeCounts.end() && iTrigTypeRate != fTrigTypeRates.end() ) { debug << "( " << *iTrigTypeCount << ", "; *iTrigTypeRate = static_cast( *iTrigTypeCount ) / run.GetRunLength(); // Uses transient run lenght from DQRunProc debug << *iTrigTypeRate << " ), "; fTrigTypeRateHist->Fill( bit, *iTrigTypeRate ); iTrigTypeCount++; // move iterator to next bit iTrigTypeRate++; // move iterator to next bit bit++; // next bit } debug << newline; // Perform checks // Check 1: N100L rate check if( fTrigTypeRates[0] >= fMinN100LRate ) // passes check { detail << "DQTimeProc::EndOfRun: run " << GetRunID() << " passed N100L trigger rate check" << newline; fCheckResult = true; } else // does not pass check { detail << "DQTimeProc::EndOfRun: run " << GetRunID() << " failed run header check" << newline; detail << " --> N100L rate " << fTrigTypeRates[0] << " < " << fMinN100LRate << " minimum" << newline; fCheckResult = false; } AddToRATDB( "n100l_trigger_rate", fTrigTypeRates[0] ); AddToRATDB( "orphans_count", fNOrphans); Update( "n100l_trigger_rate", fCheckResult ); AddToRATDB( "missing_gtids",fMissingGTIDs); if(fMissingGTIDCount > fMissingGTIDLimit){ Update("triggerProcMissingGTID",false); } else{ Update("triggerProcMissingGTID",true); } if(fBitFlipGTIDCount > fBitFlipGTIDLimit){ Update("triggerProcBitFlipGTID",false); } else{ Update("triggerProcBitFlipGTID",true); } AddToRATDB("bitflip_gtids",fBitFlipGTIDs); // Check 2: Esum hi rate check if( fTrigTypeRates[6] >= fMinEsumHRate ) // passes check { detail << "DQTimeProc::EndOfRun: run " << GetRunID() << " passed EsumH trigger rate check" << newline; fCheckResult = true; } else // does not pass check { detail << "DQTimeProc::EndOfRun: run " << GetRunID() << " failed run header check" << newline; detail << " --> EsumH rate " << fTrigTypeRates[6] << " < " << fMinEsumHRate << " minimum" << newline; fCheckResult = false; } AddToRATDB( "esumh_trigger_rate", fTrigTypeRates[6] ); Update( "esumh_trigger_rate", fCheckResult ); // Save histograms if( fOutputImages ) Plot1DHist( fTrigTypeRateHist, "", "Rate" ); WriteToRoot( fTrigTypeRateHist, "fTrigTypeRateHist" ); TCanvas *C1 = new TCanvas("C1"); string filename = fTrigTypeRateHist->GetName(); filename = filename+".png"; C1->cd(); fTrigTypeRateHist->SetXTitle("Trigger Channel"); fTrigTypeRateHist->SetYTitle("Rate (Hz)"); fTrigTypeRateHist->Draw(); C1->Print( filename.c_str() ); delete C1; C1 = NULL; // Delete objects delete fTrigTypeRateHist; fTrigTypeRateHist = NULL; DataQualityProc::EndOfRun( run ); } }// namespace RAT