CMS 3D CMS Logo

L1GlobalTriggerFDL Class Reference

Description: Final Decision Logic board. More...

#include <L1Trigger/GlobalTrigger/interface/L1GlobalTriggerFDL.h>

List of all members.

Public Member Functions

void fillDaqFdlBlock (const boost::uint16_t &activeBoardsGtDaq, const std::vector< L1GtBoard > &boardMaps, std::auto_ptr< L1GlobalTriggerReadoutRecord > &gtDaqReadoutRecord)
 fill the FDL block in the L1 GT DAQ record for iBxInEvent
void fillEvmFdlBlock (const boost::uint16_t &activeBoardsGtEvm, const std::vector< L1GtBoard > &boardMaps, std::auto_ptr< L1GlobalTriggerEvmReadoutRecord > &gtEvmReadoutRecord)
 fill the FDL block in the L1 GT EVM record for iBxInEvent
L1GtFdlWordgtFdlWord () const
 return the GtFdlWord
 L1GlobalTriggerFDL ()
 constructor
void reset ()
 clear FDL
void run (edm::Event &iEvent, const std::vector< int > &prescaleFactorsAlgoTrig, const std::vector< int > &prescaleFactorsTechTrig, const std::vector< unsigned int > &triggerMaskAlgoTrig, const std::vector< unsigned int > &triggerMaskTechTrig, const std::vector< unsigned int > &triggerMaskVetoAlgoTrig, const std::vector< unsigned int > &triggerMaskVetoTechTrig, const std::vector< L1GtBoard > &boardMaps, const int totalBxInEvent, const int iBxInEvent, const unsigned int numberPhysTriggers, const unsigned int numberTechnicalTriggers, const unsigned int numberDaqPartitions, const L1GlobalTriggerGTL *ptrGTL, const L1GlobalTriggerPSB *ptrPSB, const int pfAlgoSetIndex, const int pfTechSetIndex)
 run the FDL
virtual ~L1GlobalTriggerFDL ()
 destructor

Private Attributes

bool m_firstEv
 logical switches for the first event the first event in the luminosity segment and the first event in the run
bool m_firstEvLumiSegment
bool m_firstEvRun
L1GtFdlWordm_gtFdlWord
std::vector< std::vector< int > > m_prescaleCounterAlgoTrig
 prescale counters: NumberPhysTriggers counters per bunch cross in event
std::vector< std::vector< int > > m_prescaleCounterTechTrig
 prescale counters: technical trigger counters per bunch cross in event


Detailed Description

Description: Final Decision Logic board.

Implementation: <TODO: enter implementation details>

Author:
: M. Fierro - HEPHY Vienna - ORCA version

: Vasile Mihai Ghete - HEPHY Vienna - CMSSW version

$Date$ $Revision$

Definition at line 43 of file L1GlobalTriggerFDL.h.


Constructor & Destructor Documentation

L1GlobalTriggerFDL::L1GlobalTriggerFDL (  ) 

constructor

Definition at line 48 of file L1GlobalTriggerFDL.cc.

References m_firstEv, m_firstEvLumiSegment, m_firstEvRun, and m_gtFdlWord.

00049 {
00050 
00051     // create empty FDL word
00052     m_gtFdlWord = new L1GtFdlWord();
00053 
00054     // logical switches
00055     m_firstEv = true;
00056     m_firstEvLumiSegment = true;
00057     m_firstEvRun = true;
00058 
00059     // can not reserve memory here for prescale counters - no access to EventSetup
00060 
00061 }

L1GlobalTriggerFDL::~L1GlobalTriggerFDL (  )  [virtual]

destructor

Definition at line 64 of file L1GlobalTriggerFDL.cc.

References m_gtFdlWord, and reset().

00065 {
00066 
00067     reset();
00068     delete m_gtFdlWord;
00069 
00070 }


Member Function Documentation

void L1GlobalTriggerFDL::fillDaqFdlBlock ( const boost::uint16_t &  activeBoardsGtDaq,
const std::vector< L1GtBoard > &  boardMaps,
std::auto_ptr< L1GlobalTriggerReadoutRecord > &  gtDaqReadoutRecord 
)

fill the FDL block in the L1 GT DAQ record for iBxInEvent

Definition at line 423 of file L1GlobalTriggerFDL.cc.

References FDL, and m_gtFdlWord.

Referenced by L1GlobalTrigger::produce().

00427 {
00428 
00429     typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
00430     for (CItBoardMaps
00431             itBoard = boardMaps.begin();
00432             itBoard != boardMaps.end(); ++itBoard) {
00433 
00434         int iPosition = itBoard->gtPositionDaqRecord();
00435         if (iPosition > 0) {
00436 
00437             int iActiveBit = itBoard->gtBitDaqActiveBoards();
00438             bool activeBoard = false;
00439 
00440             if (iActiveBit >= 0) {
00441                 activeBoard = activeBoardsGtDaq & (1 << iActiveBit);
00442             }
00443 
00444             if (activeBoard && (itBoard->gtBoardType() == FDL)) {
00445 
00446                 gtDaqReadoutRecord->setGtFdlWord(*m_gtFdlWord);
00447 
00448 
00449             }
00450 
00451         }
00452 
00453     }
00454 
00455 
00456 }

void L1GlobalTriggerFDL::fillEvmFdlBlock ( const boost::uint16_t &  activeBoardsGtEvm,
const std::vector< L1GtBoard > &  boardMaps,
std::auto_ptr< L1GlobalTriggerEvmReadoutRecord > &  gtEvmReadoutRecord 
)

fill the FDL block in the L1 GT EVM record for iBxInEvent

Definition at line 459 of file L1GlobalTriggerFDL.cc.

References FDL, and m_gtFdlWord.

Referenced by L1GlobalTrigger::produce().

00463 {
00464 
00465     typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
00466     for (CItBoardMaps
00467             itBoard = boardMaps.begin();
00468             itBoard != boardMaps.end(); ++itBoard) {
00469 
00470         int iPosition = itBoard->gtPositionEvmRecord();
00471         if (iPosition > 0) {
00472 
00473             int iActiveBit = itBoard->gtBitEvmActiveBoards();
00474             bool activeBoard = false;
00475 
00476             if (iActiveBit >= 0) {
00477                 activeBoard = activeBoardsGtEvm & (1 << iActiveBit);
00478             }
00479 
00480             if (activeBoard && (itBoard->gtBoardType() == FDL)) {
00481 
00482                 gtEvmReadoutRecord->setGtFdlWord(*m_gtFdlWord);
00483 
00484 
00485             }
00486 
00487         }
00488 
00489     }
00490 
00491 }

L1GtFdlWord* L1GlobalTriggerFDL::gtFdlWord (  )  const [inline]

return the GtFdlWord

Definition at line 89 of file L1GlobalTriggerFDL.h.

References m_gtFdlWord.

00090     {
00091         return m_gtFdlWord;
00092     }

void L1GlobalTriggerFDL::reset ( void   ) 

clear FDL

Definition at line 495 of file L1GlobalTriggerFDL.cc.

References m_gtFdlWord, and L1GtFdlWord::reset().

Referenced by L1GlobalTrigger::produce(), and ~L1GlobalTriggerFDL().

00496 {
00497 
00498     m_gtFdlWord->reset();
00499 
00500     // do NOT reset the prescale counters
00501 
00502 }

void L1GlobalTriggerFDL::run ( edm::Event iEvent,
const std::vector< int > &  prescaleFactorsAlgoTrig,
const std::vector< int > &  prescaleFactorsTechTrig,
const std::vector< unsigned int > &  triggerMaskAlgoTrig,
const std::vector< unsigned int > &  triggerMaskTechTrig,
const std::vector< unsigned int > &  triggerMaskVetoAlgoTrig,
const std::vector< unsigned int > &  triggerMaskVetoTechTrig,
const std::vector< L1GtBoard > &  boardMaps,
const int  totalBxInEvent,
const int  iBxInEvent,
const unsigned int  numberPhysTriggers,
const unsigned int  numberTechnicalTriggers,
const unsigned int  numberDaqPartitions,
const L1GlobalTriggerGTL ptrGTL,
const L1GlobalTriggerPSB ptrPSB,
const int  pfAlgoSetIndex,
const int  pfTechSetIndex 
)

run the FDL

Definition at line 75 of file L1GlobalTriggerFDL.cc.

References edm::Event::bunchCrossing(), lat::endl(), FDL, L1GlobalTriggerGTL::getAlgorithmOR(), L1GlobalTriggerPSB::getGtTechnicalTriggers(), edm::Event::id(), LogDebug, edm::Event::luminosityBlock(), m_firstEv, m_firstEvLumiSegment, m_gtFdlWord, m_prescaleCounterAlgoTrig, m_prescaleCounterTechTrig, L1GlobalTriggerReadoutSetup::NumberPhysTriggers, edm::Event::orbitNumber(), L1GtFdlWord::setBoardId(), L1GtFdlWord::setBxInEvent(), L1GtFdlWord::setBxNr(), L1GtFdlWord::setEventNr(), L1GtFdlWord::setFinalOR(), L1GtFdlWord::setGtDecisionWord(), L1GtFdlWord::setGtPrescaleFactorIndexAlgo(), L1GtFdlWord::setGtPrescaleFactorIndexTech(), L1GtFdlWord::setGtTechnicalTriggerWord(), L1GtFdlWord::setLocalBxNr(), L1GtFdlWord::setLumiSegmentNr(), and L1GtFdlWord::setOrbitNr().

Referenced by L1GlobalTrigger::produce().

00092 {
00093 
00094     // FIXME get rid of bitset in GTL in order to use only EventSetup 
00095     const unsigned int numberPhysTriggersSet =
00096         L1GlobalTriggerReadoutSetup::NumberPhysTriggers;
00097 
00098     // get gtlDecisionWord from GTL
00099     std::bitset<numberPhysTriggersSet> gtlDecisionWord = ptrGTL->getAlgorithmOR();
00100 
00101     // convert decision word from std::bitset to std::vector<bool>
00102     DecisionWord algoDecisionWord(numberPhysTriggers);
00103 
00104     for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
00105 
00106         bool bitValue = gtlDecisionWord.test( iBit );
00107         algoDecisionWord[ iBit ] = bitValue;
00108     }
00109 
00110     // prescale counters are reset at the beginning of the luminosity segment
00111 
00112     if (m_firstEv) {
00113 
00114         // prescale counters: numberPhysTriggers counters per bunch cross
00115         m_prescaleCounterAlgoTrig.reserve(numberPhysTriggers*totalBxInEvent);
00116 
00117         for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
00118 
00119             m_prescaleCounterAlgoTrig.push_back(prescaleFactorsAlgoTrig);
00120         }
00121 
00122         // prescale counters: numberTechnicalTriggers counters per bunch cross
00123         m_prescaleCounterTechTrig.reserve(numberTechnicalTriggers*totalBxInEvent);
00124 
00125         for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
00126 
00127             m_prescaleCounterTechTrig.push_back(prescaleFactorsTechTrig);
00128         }
00129 
00130         m_firstEv = false;
00131     }
00132 
00133     // TODO FIXME find the beginning of the luminosity segment
00134     if (m_firstEvLumiSegment) {
00135 
00136         m_prescaleCounterAlgoTrig.clear();
00137         for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
00138             m_prescaleCounterAlgoTrig.push_back(prescaleFactorsAlgoTrig);
00139         }
00140 
00141         m_prescaleCounterTechTrig.clear();
00142         for (int iBxInEvent = 0; iBxInEvent <= totalBxInEvent; ++iBxInEvent) {
00143             m_prescaleCounterTechTrig.push_back(prescaleFactorsTechTrig);
00144         }
00145 
00146         m_firstEvLumiSegment = false;
00147 
00148     }
00149 
00150 
00151     // prescale the algorithm, if necessary
00152 
00153     // iBxInEvent is ... -2 -1 0 1 2 ... while counters are 0 1 2 3 4 ...
00154     int inBxInEvent =  totalBxInEvent/2 + iBxInEvent;
00155 
00156     for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
00157 
00158         if (prescaleFactorsAlgoTrig.at(iBit) != 1) {
00159 
00160             bool bitValue = algoDecisionWord.at( iBit );
00161             if (bitValue) {
00162 
00163                 (m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit))--;
00164                 if (m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit) == 0) {
00165 
00166                     // bit already true in algoDecisionWord, just reset counter
00167                     m_prescaleCounterAlgoTrig.at(inBxInEvent).at(iBit) = 
00168                         prescaleFactorsAlgoTrig.at(iBit);
00169 
00170                     //LogTrace("L1GlobalTriggerFDL")
00171                     //<< "\nPrescaled algorithm: " << iBit << ". Reset counter to "
00172                     //<< prescaleFactorsAlgoTrig.at(iBit) << "\n"
00173                     //<< std::endl;
00174 
00175                 } else {
00176 
00177                     // change bit to false
00178                     algoDecisionWord[iBit] = false;;
00179 
00180                     //LogTrace("L1GlobalTriggerFDL")
00181                     //<< "\nPrescaled algorithm: " << iBit << ". Result set to false"
00182                     //<< std::endl;
00183 
00184                 }
00185             }
00186         }
00187     }
00188 
00189     // algo decision word written in the FDL readout before the trigger mask 
00190     // in order to allow multiple DAQ partitions
00191 
00192     //
00193     // technical triggers
00194     //
00195     
00196     std::vector<bool> techDecisionWord = *(ptrPSB->getGtTechnicalTriggers());
00197     
00198     // prescale the technical trigger, if necessary
00199 
00200     for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
00201 
00202         if (prescaleFactorsTechTrig.at(iBit) != 1) {
00203 
00204             bool bitValue = techDecisionWord.at( iBit );
00205             if (bitValue) {
00206 
00207                 (m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit))--;
00208                 if (m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit) == 0) {
00209 
00210                     // bit already true in techDecisionWord, just reset counter
00211                     m_prescaleCounterTechTrig.at(inBxInEvent).at(iBit) = 
00212                         prescaleFactorsTechTrig.at(iBit);
00213 
00214                     //LogTrace("L1GlobalTriggerFDL")
00215                     //<< "\nPrescaled algorithm: " << iBit << ". Reset counter to "
00216                     //<< prescaleFactorsTechTrig.at(iBit) << "\n"
00217                     //<< std::endl;
00218 
00219                 } else {
00220 
00221                     // change bit to false
00222                     techDecisionWord[iBit] = false;
00223 
00224                     //LogTrace("L1GlobalTriggerFDL")
00225                     //<< "\nPrescaled technical trigger: " << iBit << ". Result set to false"
00226                     //<< std::endl;
00227 
00228                 }
00229             }
00230         }
00231     }
00232 
00233     // technical trigger decision word written in the FDL readout before the trigger mask 
00234     // in order to allow multiple DAQ partitions
00235     
00236     //
00237     // compute the final decision word per DAQ partition
00238     //
00239 
00240     boost::uint16_t finalOrValue = 0;
00241 
00242     for (unsigned int iDaq = 0; iDaq < numberDaqPartitions; ++iDaq) {
00243 
00244         bool daqPartitionFinalOR = false;
00245 
00246         // starts with technical trigger veto mask to minimize computation
00247         // no algorithm trigger veto mask is implemented up to now in hardware,
00248         // therefore do not implement it here
00249         bool vetoTechTrig = false;
00250 
00251         for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
00252 
00253             int triggerMaskVetoTechTrigBit = 
00254                 triggerMaskVetoTechTrig[iBit] & (1 << iDaq);
00255             //LogTrace("L1GlobalTriggerFDL")
00256             //<< "\nTechnical trigger bit: " << iBit
00257             //<< " mask = " << triggerMaskVetoTechTrigBit 
00258             //<< " DAQ partition " << iDaq
00259             //<< std::endl;
00260 
00261             if (triggerMaskVetoTechTrigBit && techDecisionWord[iBit]) {
00262 
00263                 daqPartitionFinalOR = false;
00264                 vetoTechTrig = true;
00265 
00266                 //LogTrace("L1GlobalTriggerFDL")
00267                 //<< "\nVeto mask technical trigger: " << iBit 
00268                 // << ". FinalOR for DAQ partition " << iDaq << " set to false"
00269                 //<< std::endl;
00270 
00271                 break;
00272             }
00273 
00274         }
00275 
00276         // apply algorithm and technical trigger masks only if no veto from technical trigger
00277         if (!vetoTechTrig) {
00278 
00279             // algorithm trigger mask
00280             bool algoFinalOr = false;
00281 
00282             for (unsigned int iBit = 0; iBit < numberPhysTriggers; ++iBit) {
00283 
00284                 bool iBitDecision = false;
00285                 
00286                 int triggerMaskAlgoTrigBit = triggerMaskAlgoTrig[iBit] & (1 << iDaq);
00287                 //LogTrace("L1GlobalTriggerFDL")
00288                 //<< "\nAlgorithm trigger bit: " << iBit 
00289                 //<< " mask = " << triggerMaskAlgoTrigBit
00290                 //<< " DAQ partition " << iDaq
00291                 //<< std::endl;
00292 
00293                 if (triggerMaskAlgoTrigBit) {
00294                     iBitDecision = false;
00295 
00296                     //LogTrace("L1GlobalTriggerFDL")
00297                     //<< "\nMasked algorithm trigger: " << iBit << ". Result set to false"
00298                     //<< std::endl;
00299                 } else {
00300                     iBitDecision = algoDecisionWord[iBit];
00301                 }
00302 
00303                 algoFinalOr = algoFinalOr || iBitDecision;
00304 
00305             }
00306 
00307             // set the technical trigger mask: block the corresponding algorithm if bit value is 1
00308 
00309             bool techFinalOr = false;
00310 
00311             for (unsigned int iBit = 0; iBit < numberTechnicalTriggers; ++iBit) {
00312 
00313                 bool iBitDecision = false;
00314 
00315                 int triggerMaskTechTrigBit = triggerMaskTechTrig[iBit] & (1 << iDaq);
00316                 //LogTrace("L1GlobalTriggerFDL")
00317                 //<< "\nTechnical trigger bit: " << iBit 
00318                 //<< " mask = " << triggerMaskTechTrigBit
00319                 //<< std::endl;
00320 
00321                 if (triggerMaskTechTrigBit) {
00322 
00323                     iBitDecision = false;
00324 
00325                     //LogTrace("L1GlobalTriggerFDL")
00326                     //<< "\nMasked technical trigger: " << iBit << ". Result set to false"
00327                     //<< std::endl;
00328                 } else {
00329                     iBitDecision = techDecisionWord[iBit];
00330                 }
00331 
00332                 techFinalOr = techFinalOr || iBitDecision;
00333             }
00334             
00335             daqPartitionFinalOR = algoFinalOr || techFinalOr;
00336             
00337         } else {
00338             
00339             daqPartitionFinalOR = false; // vetoTechTrig 
00340         
00341         }
00342         
00343         // push it in finalOrValue
00344         boost::uint16_t daqPartitionFinalORValue = 
00345             static_cast<boost::uint16_t>(daqPartitionFinalOR);
00346             
00347         finalOrValue = finalOrValue | (daqPartitionFinalORValue << iDaq);
00348 
00349     }
00350     
00351     
00352     // fill everything we know in the L1GtFdlWord
00353 
00354     typedef std::vector<L1GtBoard>::const_iterator CItBoardMaps;
00355     for (CItBoardMaps
00356             itBoard = boardMaps.begin();
00357             itBoard != boardMaps.end(); ++itBoard) {
00358 
00359 
00360         if ((itBoard->gtBoardType() == FDL)) {
00361 
00362             m_gtFdlWord->setBoardId( itBoard->gtBoardId() );
00363 
00364             // BxInEvent
00365             m_gtFdlWord->setBxInEvent(iBxInEvent);
00366             
00367             // bunch crossing
00368             
00369             // fill in emulator the same bunch crossing (12 bits - hardwired number of bits...)
00370             // and the same local bunch crossing for all boards
00371             int bxCross = iEvent.bunchCrossing();
00372             boost::uint16_t bxCrossHw = 0;
00373             if ((bxCross & 0xFFF) == bxCross) {
00374                 bxCrossHw = static_cast<boost::uint16_t> (bxCross);
00375             }
00376             else {
00377                 bxCrossHw = 0; // Bx number too large, set to 0!
00378                 LogDebug("L1GlobalTrigger")
00379                     << "\nBunch cross number [hex] = "
00380                     << std::hex << bxCross
00381                     << "\n  larger than 12 bits. Set to 0! \n"
00382                     << std::dec << std::endl;
00383             }
00384 
00385             m_gtFdlWord->setBxNr(bxCrossHw);
00386 
00387             // set event number since last L1 reset generated in FDL
00388             m_gtFdlWord->setEventNr(
00389                 static_cast<boost::uint32_t>(iEvent.id().event()) );
00390 
00391             // technical trigger decision word
00392             m_gtFdlWord->setGtTechnicalTriggerWord(techDecisionWord);
00393 
00394             // algorithm trigger decision word
00395             m_gtFdlWord->setGtDecisionWord(algoDecisionWord);
00396             
00397             // index of prescale factor set - technical triggers and algo
00398             m_gtFdlWord->setGtPrescaleFactorIndexTech(static_cast<boost::uint16_t>(pfTechSetIndex));
00399             m_gtFdlWord->setGtPrescaleFactorIndexAlgo(static_cast<boost::uint16_t>(pfAlgoSetIndex));
00400 
00401             // NoAlgo bit FIXME
00402             
00403             // finalOR
00404             m_gtFdlWord->setFinalOR(finalOrValue);
00405 
00406             // orbit number
00407             m_gtFdlWord->setOrbitNr(static_cast<boost::uint32_t>(iEvent.orbitNumber()) );
00408 
00409             // luminosity segment number
00410             m_gtFdlWord->setLumiSegmentNr(static_cast<boost::uint16_t>(iEvent.luminosityBlock()));
00411 
00412             // local bunch crossing - set identical with absolute BxNr
00413             m_gtFdlWord->setLocalBxNr(bxCrossHw);
00414 
00415 
00416         }
00417 
00418     }
00419 
00420 }


Member Data Documentation

bool L1GlobalTriggerFDL::m_firstEv [private]

logical switches for the first event the first event in the luminosity segment and the first event in the run

Definition at line 110 of file L1GlobalTriggerFDL.h.

Referenced by L1GlobalTriggerFDL(), and run().

bool L1GlobalTriggerFDL::m_firstEvLumiSegment [private]

Definition at line 111 of file L1GlobalTriggerFDL.h.

Referenced by L1GlobalTriggerFDL(), and run().

bool L1GlobalTriggerFDL::m_firstEvRun [private]

Definition at line 112 of file L1GlobalTriggerFDL.h.

Referenced by L1GlobalTriggerFDL().

L1GtFdlWord* L1GlobalTriggerFDL::m_gtFdlWord [private]

Definition at line 98 of file L1GlobalTriggerFDL.h.

Referenced by fillDaqFdlBlock(), fillEvmFdlBlock(), gtFdlWord(), L1GlobalTriggerFDL(), reset(), run(), and ~L1GlobalTriggerFDL().

std::vector<std::vector<int> > L1GlobalTriggerFDL::m_prescaleCounterAlgoTrig [private]

prescale counters: NumberPhysTriggers counters per bunch cross in event

Definition at line 101 of file L1GlobalTriggerFDL.h.

Referenced by run().

std::vector<std::vector<int> > L1GlobalTriggerFDL::m_prescaleCounterTechTrig [private]

prescale counters: technical trigger counters per bunch cross in event

Definition at line 104 of file L1GlobalTriggerFDL.h.

Referenced by run().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:26:40 2009 for CMSSW by  doxygen 1.5.4