Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:20:35

0001 /**
0002  * \class CorrThreeBodyCondition
0003  *
0004  * \orig author: Elisa Fontanesi - Boston University
0005  *               CorrCondition and CorrWithOverlapRemovalCondition classes used as a starting point
0006  *
0007  * Description: L1 Global Trigger three-body correlation conditions:                                                                                     
0008  *              evaluation of a three-body correlation condition (= three-muon invariant mass)
0009  *
0010  */
0011 
0012 // this class header
0013 #include "L1Trigger/L1TGlobal/interface/CorrCondition.h"
0014 #include "L1Trigger/L1TGlobal/interface/CorrThreeBodyCondition.h"
0015 
0016 // system include files
0017 #include <iostream>
0018 #include <iomanip>
0019 
0020 #include <string>
0021 #include <vector>
0022 #include <algorithm>
0023 
0024 // user include files
0025 //   base classes
0026 #include "L1Trigger/L1TGlobal/interface/CorrelationTemplate.h"
0027 #include "L1Trigger/L1TGlobal/interface/CorrelationThreeBodyTemplate.h"
0028 #include "L1Trigger/L1TGlobal/interface/ConditionEvaluation.h"
0029 
0030 #include "L1Trigger/L1TGlobal/interface/MuCondition.h"
0031 #include "L1Trigger/L1TGlobal/interface/MuonTemplate.h"
0032 #include "L1Trigger/L1TGlobal/interface/GlobalScales.h"
0033 #include "L1Trigger/L1TGlobal/interface/GlobalBoard.h"
0034 
0035 #include "DataFormats/L1Trigger/interface/L1Candidate.h"
0036 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0037 #include "FWCore/MessageLogger/interface/MessageDrop.h"
0038 
0039 // constructors
0040 //     default
0041 l1t::CorrThreeBodyCondition::CorrThreeBodyCondition() : ConditionEvaluation() {}
0042 
0043 //     from base template condition (from event setup usually)
0044 l1t::CorrThreeBodyCondition::CorrThreeBodyCondition(const GlobalCondition* corrTemplate,
0045                                                     const GlobalCondition* cond0Condition,
0046                                                     const GlobalCondition* cond1Condition,
0047                                                     const GlobalCondition* cond2Condition,
0048                                                     const GlobalBoard* ptrGTB)
0049     : ConditionEvaluation(),
0050       m_gtCorrelationThreeBodyTemplate(static_cast<const CorrelationThreeBodyTemplate*>(corrTemplate)),
0051       m_gtCond0(cond0Condition),
0052       m_gtCond1(cond1Condition),
0053       m_gtCond2(cond2Condition),
0054       m_uGtB(ptrGTB) {}
0055 
0056 // copy constructor
0057 void l1t::CorrThreeBodyCondition::copy(const l1t::CorrThreeBodyCondition& cp) {
0058   m_gtCorrelationThreeBodyTemplate = cp.gtCorrelationThreeBodyTemplate();
0059   m_uGtB = cp.getuGtB();
0060 
0061   m_condMaxNumberObjects = cp.condMaxNumberObjects();
0062   m_condLastResult = cp.condLastResult();
0063   m_combinationsInCond = cp.getCombinationsInCond();
0064 
0065   m_verbosity = cp.m_verbosity;
0066 }
0067 
0068 l1t::CorrThreeBodyCondition::CorrThreeBodyCondition(const l1t::CorrThreeBodyCondition& cp) : ConditionEvaluation() {
0069   copy(cp);
0070 }
0071 
0072 // destructor
0073 l1t::CorrThreeBodyCondition::~CorrThreeBodyCondition() {
0074   // empty
0075 }
0076 
0077 // equal operator
0078 l1t::CorrThreeBodyCondition& l1t::CorrThreeBodyCondition::operator=(const l1t::CorrThreeBodyCondition& cp) {
0079   copy(cp);
0080   return *this;
0081 }
0082 
0083 ///   set the pointer to uGT GlobalBoard
0084 void l1t::CorrThreeBodyCondition::setuGtB(const GlobalBoard* ptrGTB) { m_uGtB = ptrGTB; }
0085 
0086 void l1t::CorrThreeBodyCondition::setScales(const GlobalScales* sc) { m_gtScales = sc; }
0087 
0088 // try all object permutations and check spatial correlations, if required
0089 const bool l1t::CorrThreeBodyCondition::evaluateCondition(const int bxEval) const {
0090   if (m_verbosity) {
0091     std::ostringstream myCout;
0092     m_gtCorrelationThreeBodyTemplate->print(myCout);
0093     LogDebug("L1TGlobal") << "Three-body Correlation Condition Evaluation..." << std::endl;
0094   }
0095 
0096   bool condResult = false;
0097 
0098   // number of objects in the condition (three) and their type
0099   int nObjInCond = 3;
0100   std::vector<GlobalObject> cndObjTypeVec(nObjInCond);
0101 
0102   // evaluate first the three subconditions (Type1s)
0103   const GtConditionCategory cond0Categ = m_gtCorrelationThreeBodyTemplate->cond0Category();
0104   const GtConditionCategory cond1Categ = m_gtCorrelationThreeBodyTemplate->cond1Category();
0105   const GtConditionCategory cond2Categ = m_gtCorrelationThreeBodyTemplate->cond2Category();
0106 
0107   const MuonTemplate* corrMuon = nullptr;
0108 
0109   CombinationsInCond cond0Comb;
0110   CombinationsInCond cond1Comb;
0111   CombinationsInCond cond2Comb;
0112 
0113   int cond0bx(0);
0114   int cond1bx(0);
0115   int cond2bx(0);
0116 
0117   // FIRST OBJECT
0118   bool reqObjResult = false;
0119 
0120   if (cond0Categ == CondMuon) {
0121     LogDebug("L1TGlobal") << "\n --------------------- First muon checks ---------------------" << std::endl;
0122     corrMuon = static_cast<const MuonTemplate*>(m_gtCond0);
0123     MuCondition muCondition(corrMuon, m_uGtB, 0, 0);
0124 
0125     muCondition.evaluateConditionStoreResult(bxEval);
0126     reqObjResult = muCondition.condLastResult();
0127 
0128     cond0Comb = (muCondition.getCombinationsInCond());
0129     cond0bx = bxEval + (corrMuon->condRelativeBx());
0130     cndObjTypeVec[0] = (corrMuon->objectType())[0];
0131 
0132     if (m_verbosity) {
0133       std::ostringstream myCout;
0134       muCondition.print(myCout);
0135       LogDebug("L1TGlobal") << myCout.str() << std::endl;
0136     }
0137   } else {
0138     // Interested only in three-muon correlations
0139     LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 0" << std::endl;
0140     return false;
0141   }
0142 
0143   // return if first subcondition is false
0144   if (!reqObjResult) {
0145     LogDebug("L1TGlobal") << "\n  First subcondition false, second subcondition not evaluated and not printed."
0146                           << std::endl;
0147     return false;
0148   }
0149 
0150   // SECOND OBJECT
0151   reqObjResult = false;
0152 
0153   if (cond1Categ == CondMuon) {
0154     LogDebug("L1TGlobal") << "\n --------------------- Second muon checks ---------------------" << std::endl;
0155     corrMuon = static_cast<const MuonTemplate*>(m_gtCond1);
0156     MuCondition muCondition(corrMuon, m_uGtB, 0, 0);
0157 
0158     muCondition.evaluateConditionStoreResult(bxEval);
0159     reqObjResult = muCondition.condLastResult();
0160 
0161     cond1Comb = (muCondition.getCombinationsInCond());
0162     cond1bx = bxEval + (corrMuon->condRelativeBx());
0163     cndObjTypeVec[1] = (corrMuon->objectType())[0];
0164 
0165     if (m_verbosity) {
0166       std::ostringstream myCout;
0167       muCondition.print(myCout);
0168       LogDebug("L1TGlobal") << myCout.str() << std::endl;
0169     }
0170   }
0171 
0172   else {
0173     // Interested only in three-muon correlations
0174     LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 1" << std::endl;
0175     return false;
0176   }
0177 
0178   // return if second subcondition is false
0179   if (!reqObjResult) {
0180     LogDebug("L1TGlobal") << "\n  Second subcondition false, third subcondition not evaluated and not printed."
0181                           << std::endl;
0182     return false;
0183   }
0184 
0185   // THIRD OBJECT
0186   reqObjResult = false;
0187 
0188   if (cond2Categ == CondMuon) {
0189     LogDebug("L1TGlobal") << "\n --------------------- Third muon checks ---------------------" << std::endl;
0190     corrMuon = static_cast<const MuonTemplate*>(m_gtCond2);
0191     MuCondition muCondition(corrMuon, m_uGtB, 0, 0);
0192 
0193     muCondition.evaluateConditionStoreResult(bxEval);
0194     reqObjResult = muCondition.condLastResult();
0195 
0196     cond2Comb = (muCondition.getCombinationsInCond());
0197     cond2bx = bxEval + (corrMuon->condRelativeBx());
0198     cndObjTypeVec[2] = (corrMuon->objectType())[0];
0199 
0200     if (m_verbosity) {
0201       std::ostringstream myCout;
0202       muCondition.print(myCout);
0203       LogDebug("L1TGlobal") << myCout.str() << std::endl;
0204     }
0205   }
0206 
0207   else {
0208     // Interested only in three-muon correlations
0209     LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 2" << std::endl;
0210     return false;
0211   }
0212 
0213   // Return if third subcondition is false
0214   if (!reqObjResult) {
0215     return false;
0216   } else {
0217     LogDebug("L1TGlobal")
0218         << "\n"
0219         << "Found three objects satisfying subconditions: evaluate three-body correlation requirements.\n"
0220         << std::endl;
0221   }
0222 
0223   // Since we have three good legs get the correlation parameters
0224   CorrelationThreeBodyTemplate::CorrelationThreeBodyParameter corrPar =
0225       *(m_gtCorrelationThreeBodyTemplate->correlationThreeBodyParameter());
0226 
0227   // Vector to store the indices of the objects involved in the condition evaluation
0228   SingleCombInCond objectsInComb;
0229   objectsInComb.reserve(nObjInCond);
0230 
0231   // Clear the m_combinationsInCond vector:
0232   // it will store the set of objects satisfying the condition evaluated as true
0233   (combinationsInCond()).clear();
0234 
0235   // Pointers to objects
0236   const BXVector<const l1t::Muon*>* candMuVec = nullptr;
0237 
0238   // Make the conversions of the indices, depending on the combination of objects involved
0239   int phiIndex0 = 0;
0240   double phi0Phy = 0.;
0241   int phiIndex1 = 0;
0242   double phi1Phy = 0.;
0243   int phiIndex2 = 0;
0244   double phi2Phy = 0.;
0245 
0246   int etaIndex0 = 0;
0247   double eta0Phy = 0.;
0248   int etaBin0 = 0;
0249   int etaIndex1 = 0;
0250   double eta1Phy = 0.;
0251   int etaBin1 = 0;
0252   int etaIndex2 = 0;
0253   double eta2Phy = 0.;
0254   int etaBin2 = 0;
0255 
0256   int etIndex0 = 0;
0257   int etBin0 = 0;
0258   double et0Phy = 0.;
0259   int etIndex1 = 0;
0260   int etBin1 = 0;
0261   double et1Phy = 0.;
0262   int etIndex2 = 0;
0263   int etBin2 = 0;
0264   double et2Phy = 0.;
0265 
0266   // Charges to take into account the charge correlation
0267   int chrg0 = -1;
0268   int chrg1 = -1;
0269   int chrg2 = -1;
0270 
0271   // Determine the number of phi bins to get cutoff at pi
0272   int phiBound = 0;
0273   if (cond0Categ == CondMuon || cond1Categ == CondMuon || cond2Categ == CondMuon) {
0274     const GlobalScales::ScaleParameters& par = m_gtScales->getMUScales();
0275     phiBound = (int)((par.phiMax - par.phiMin) / par.phiStep) / 2;
0276   } else {
0277     //Assumes all objects are on same phi scale
0278     const GlobalScales::ScaleParameters& par = m_gtScales->getEGScales();
0279     phiBound = (int)((par.phiMax - par.phiMin) / par.phiStep) / 2;
0280   }
0281   LogDebug("L1TGlobal") << "Phi Bound = " << phiBound << std::endl;
0282 
0283   // Keep track of objects for LUTS
0284   std::string lutObj0 = "NULL";
0285   std::string lutObj1 = "NULL";
0286   std::string lutObj2 = "NULL";
0287 
0288   LogTrace("L1TGlobal") << "  Number of objects satisfying the subcondition 0: " << (cond0Comb.size()) << std::endl;
0289   LogTrace("L1TGlobal") << "  Number of objects satisfying the subcondition 1: " << (cond1Comb.size()) << std::endl;
0290   LogTrace("L1TGlobal") << "  Number of objects satisfying the subcondition 2: " << (cond2Comb.size()) << std::endl;
0291 
0292   ////////////////////////////////
0293   // LOOP OVER ALL COMBINATIONS //
0294   ////////////////////////////////
0295   unsigned int preShift = 0;
0296 
0297   // *** Looking for a set of three objects
0298   for (std::vector<SingleCombInCond>::const_iterator it0Comb = cond0Comb.begin(); it0Comb != cond0Comb.end();
0299        it0Comb++) {
0300     // Type1s: there is 1 object only, no need for a loop, index 0 should be OK in (*it0Comb)[0]
0301     // ... but add protection to not crash
0302     LogDebug("L1TGlobal") << "Looking at first subcondition" << std::endl;
0303     int obj0Index = -1;
0304 
0305     if (!(*it0Comb).empty()) {
0306       obj0Index = (*it0Comb)[0];
0307     } else {
0308       LogTrace("L1TGlobal") << "\n  SingleCombInCond (*it0Comb).size() " << ((*it0Comb).size()) << std::endl;
0309       return false;
0310     }
0311 
0312     // FIRST OBJECT: Collect the information on the first leg of the correlation
0313     if (cond0Categ == CondMuon) {
0314       lutObj0 = "MU";
0315       candMuVec = m_uGtB->getCandL1Mu();
0316       phiIndex0 = (candMuVec->at(cond0bx, obj0Index))->hwPhiAtVtx();
0317       etaIndex0 = (candMuVec->at(cond0bx, obj0Index))->hwEtaAtVtx();
0318       etIndex0 = (candMuVec->at(cond0bx, obj0Index))->hwPt();
0319       chrg0 = (candMuVec->at(cond0bx, obj0Index))->hwCharge();
0320 
0321       etaBin0 = etaIndex0;
0322       if (etaBin0 < 0)
0323         etaBin0 = m_gtScales->getMUScales().etaBins.size() + etaBin0;
0324 
0325       etBin0 = etIndex0;
0326       int ssize = m_gtScales->getMUScales().etBins.size();
0327       if (etBin0 >= ssize) {
0328         etBin0 = ssize - 1;
0329         LogTrace("L1TGlobal") << "MU0 hw et" << etBin0 << " out of scale range.  Setting to maximum.";
0330       }
0331 
0332       // Determine Floating Pt numbers for floating point calculation
0333       std::pair<double, double> binEdges = m_gtScales->getMUScales().phiBins.at(phiIndex0);
0334       phi0Phy = 0.5 * (binEdges.second + binEdges.first);
0335       binEdges = m_gtScales->getMUScales().etaBins.at(etaBin0);
0336       eta0Phy = 0.5 * (binEdges.second + binEdges.first);
0337       binEdges = m_gtScales->getMUScales().etBins.at(etBin0);
0338       et0Phy = 0.5 * (binEdges.second + binEdges.first);
0339       LogDebug("L1TGlobal") << "Found all quantities for MU0" << std::endl;
0340     } else {
0341       // Interested only in three-muon correlations
0342       LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 0" << std::endl;
0343       return false;
0344     }
0345 
0346     // SECOND OBJECT: Now loop over the second leg to get its information
0347     for (std::vector<SingleCombInCond>::const_iterator it1Comb = cond1Comb.begin(); it1Comb != cond1Comb.end();
0348          it1Comb++) {
0349       LogDebug("L1TGlobal") << "Looking at second subcondition" << std::endl;
0350       int obj1Index = -1;
0351 
0352       if (!(*it1Comb).empty()) {
0353         obj1Index = (*it1Comb)[0];
0354       } else {
0355         LogTrace("L1TGlobal") << "\n  SingleCombInCond (*it1Comb).size() " << ((*it1Comb).size()) << std::endl;
0356         return false;
0357       }
0358 
0359       // If we are dealing with the same object type avoid the two legs either being the same object
0360       if (cndObjTypeVec[0] == cndObjTypeVec[1] && obj0Index == obj1Index && cond0bx == cond1bx) {
0361         LogDebug("L1TGlobal") << "Corr Condition looking at same leg...skip" << std::endl;
0362         continue;
0363       }
0364 
0365       if (cond1Categ == CondMuon) {
0366         lutObj1 = "MU";
0367         candMuVec = m_uGtB->getCandL1Mu();
0368         phiIndex1 = (candMuVec->at(cond1bx, obj1Index))->hwPhiAtVtx();
0369         etaIndex1 = (candMuVec->at(cond1bx, obj1Index))->hwEtaAtVtx();
0370         etIndex1 = (candMuVec->at(cond1bx, obj1Index))->hwPt();
0371         chrg1 = (candMuVec->at(cond1bx, obj1Index))->hwCharge();
0372 
0373         etaBin1 = etaIndex1;
0374         if (etaBin1 < 0)
0375           etaBin1 = m_gtScales->getMUScales().etaBins.size() + etaBin1;
0376 
0377         etBin1 = etIndex1;
0378         int ssize = m_gtScales->getMUScales().etBins.size();
0379         if (etBin1 >= ssize) {
0380           LogTrace("L1TGlobal") << "MU1 hw et" << etBin1 << " out of scale range.  Setting to maximum.";
0381           etBin1 = ssize - 1;
0382         }
0383 
0384         // Determine Floating Pt numbers for floating point calculation
0385         std::pair<double, double> binEdges = m_gtScales->getMUScales().phiBins.at(phiIndex1);
0386         phi1Phy = 0.5 * (binEdges.second + binEdges.first);
0387         binEdges = m_gtScales->getMUScales().etaBins.at(etaBin1);
0388         eta1Phy = 0.5 * (binEdges.second + binEdges.first);
0389         binEdges = m_gtScales->getMUScales().etBins.at(etBin1);
0390         et1Phy = 0.5 * (binEdges.second + binEdges.first);
0391         LogDebug("L1TGlobal") << "Found all quantities for MU1" << std::endl;
0392       } else {
0393         // Interested only in three-muon correlations
0394         LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 1" << std::endl;
0395         return false;
0396       }
0397 
0398       // THIRD OBJECT: Finally loop over the third leg to get its information
0399       for (std::vector<SingleCombInCond>::const_iterator it2Comb = cond2Comb.begin(); it2Comb != cond2Comb.end();
0400            it2Comb++) {
0401         LogDebug("L1TGlobal") << "Looking at the third object for the three-body condition" << std::endl;
0402         int obj2Index = -1;
0403 
0404         if (!(*it2Comb).empty()) {
0405           obj2Index = (*it2Comb)[0];
0406         } else {
0407           LogTrace("L1TGlobal") << "\n  SingleCombInCond (*it2Comb).size() " << ((*it2Comb).size()) << std::endl;
0408           return false;
0409         }
0410 
0411         // If we are dealing with the same object type avoid the two legs
0412         // either being the same object
0413         if ((cndObjTypeVec[0] == cndObjTypeVec[2] && obj0Index == obj2Index && cond0bx == cond2bx) ||
0414             (cndObjTypeVec[1] == cndObjTypeVec[2] && obj1Index == obj2Index && cond1bx == cond2bx)) {
0415           LogDebug("L1TGlobal") << "Corr Condition looking at same leg...skip" << std::endl;
0416           continue;
0417         }
0418 
0419         if (cond2Categ == CondMuon) {
0420           lutObj2 = "MU";
0421           candMuVec = m_uGtB->getCandL1Mu();
0422           phiIndex2 = (candMuVec->at(cond2bx, obj2Index))->hwPhiAtVtx();
0423           etaIndex2 = (candMuVec->at(cond2bx, obj2Index))->hwEtaAtVtx();
0424           etIndex2 = (candMuVec->at(cond2bx, obj2Index))->hwPt();
0425           chrg2 = (candMuVec->at(cond2bx, obj2Index))->hwCharge();
0426 
0427           etaBin2 = etaIndex2;
0428           if (etaBin2 < 0)
0429             etaBin2 = m_gtScales->getMUScales().etaBins.size() + etaBin2;
0430 
0431           etBin2 = etIndex2;
0432           int ssize = m_gtScales->getMUScales().etBins.size();
0433           if (etBin2 >= ssize) {
0434             LogTrace("L1TGlobal") << "MU2 hw et" << etBin2 << " out of scale range.  Setting to maximum.";
0435             etBin2 = ssize - 1;
0436           }
0437 
0438           // Determine Floating Pt numbers for floating point calculation
0439           std::pair<double, double> binEdges = m_gtScales->getMUScales().phiBins.at(phiIndex2);
0440           phi2Phy = 0.5 * (binEdges.second + binEdges.first);
0441           binEdges = m_gtScales->getMUScales().etaBins.at(etaBin2);
0442           eta2Phy = 0.5 * (binEdges.second + binEdges.first);
0443           binEdges = m_gtScales->getMUScales().etBins.at(etBin2);
0444           et2Phy = 0.5 * (binEdges.second + binEdges.first);
0445           LogDebug("L1TGlobal") << "Found all quantities for MU2" << std::endl;
0446         }
0447 
0448         else {
0449           // Interested only in three-muon correlations
0450           LogDebug("L1TGlobal") << "CondMuon not satisfied for Leg 2" << std::endl;
0451           return false;
0452         };
0453 
0454         if (m_verbosity) {
0455           LogDebug("L1TGlobal") << "\n >>>>>> THREE-MUON EVENT!" << std::endl;
0456           LogDebug("L1TGlobal") << ">>>>>> Object involved in the three-body correlation condition are ["
0457                                 << l1t::GlobalObjectEnumToString(cndObjTypeVec[0]) << ", "
0458                                 << l1t::GlobalObjectEnumToString(cndObjTypeVec[1]) << ", "
0459                                 << l1t::GlobalObjectEnumToString(cndObjTypeVec[2]) << "] with collection indices ["
0460                                 << obj0Index << ", " << obj1Index << obj2Index << "] "
0461                                 << " having: \n"
0462                                 << "     Et  values  = [" << etIndex0 << ", " << etIndex1 << ", " << etIndex2 << "]\n"
0463                                 << "     phi indices = [" << phiIndex0 << ", " << phiIndex1 << ", " << phiIndex2
0464                                 << "]\n"
0465                                 << "     eta indices = [" << etaIndex0 << ", " << etaIndex1 << ", " << etaIndex2
0466                                 << "]\n"
0467                                 << "     charge values = [" << chrg0 << ", " << chrg1 << ", " << chrg2 << "]\n";
0468         }
0469 
0470         // Now perform the desired correlation on these three objects:
0471         // reqResult will be set true in case all checks were successful for a given combination of three muons
0472         bool reqResult = false;
0473         bool chrgCorrel = true;
0474 
0475         // Check the three-muon charge correlation, if requested.
0476         // NOTE that the charge can be 1 (positive) or 0 (negative), so [SS] SUM(chrg) == 3 OR 0, [OS] SUM(chrg) == 1 OR 2
0477         if (cond0Categ == CondMuon && cond1Categ == CondMuon && cond2Categ == CondMuon) {
0478           // Check for opp-sign requirement:
0479           if (corrPar.chargeCorrelation == 4 && ((chrg0 + chrg1 + chrg2) == 3 || (chrg0 + chrg1 + chrg2) == 0)) {
0480             chrgCorrel = false;
0481           }
0482           // Check for same-sign
0483           if (corrPar.chargeCorrelation == 2 && ((chrg0 + chrg1 + chrg2) == 1 || (chrg0 + chrg1 + chrg2) == 2)) {
0484             chrgCorrel = false;
0485           }
0486           // Ignore the charge correlation requirement
0487           if (corrPar.chargeCorrelation == 1) {
0488             chrgCorrel = true;
0489           }
0490         }
0491 
0492         // Clear the vector containing indices of the objects of the combination involved in the condition evaluation
0493         objectsInComb.clear();
0494         objectsInComb.push_back(obj0Index);
0495         objectsInComb.push_back(obj1Index);
0496         objectsInComb.push_back(obj2Index);
0497 
0498         // Delta eta and phi calculations needed to evaluate the three-body invariant mass
0499         double deltaPhiPhy_01 = fabs(phi1Phy - phi0Phy);
0500         if (deltaPhiPhy_01 > M_PI)
0501           deltaPhiPhy_01 = 2. * M_PI - deltaPhiPhy_01;
0502         double deltaEtaPhy_01 = fabs(eta1Phy - eta0Phy);
0503 
0504         double deltaPhiPhy_02 = fabs(phi2Phy - phi0Phy);
0505         if (deltaPhiPhy_02 > M_PI)
0506           deltaPhiPhy_02 = 2. * M_PI - deltaPhiPhy_02;
0507         double deltaEtaPhy_02 = fabs(eta2Phy - eta0Phy);
0508 
0509         double deltaPhiPhy_12 = fabs(phi2Phy - phi1Phy);
0510         if (deltaPhiPhy_12 > M_PI)
0511           deltaPhiPhy_12 = 2. * M_PI - deltaPhiPhy_12;
0512         double deltaEtaPhy_12 = fabs(eta2Phy - eta1Phy);
0513 
0514         // Determine the integer based delta eta and delta phi
0515         int deltaPhiFW_01 = abs(phiIndex0 - phiIndex1);
0516         if (deltaPhiFW_01 >= phiBound)
0517           deltaPhiFW_01 = 2 * phiBound - deltaPhiFW_01;
0518         std::string lutName_01 = lutObj0;
0519         lutName_01 += "-";
0520         lutName_01 += lutObj1;
0521         long long deltaPhiLUT_01 = m_gtScales->getLUT_DeltaPhi(lutName_01, deltaPhiFW_01);
0522         unsigned int precDeltaPhiLUT_01 = m_gtScales->getPrec_DeltaPhi(lutName_01);
0523 
0524         int deltaEtaFW_01 = abs(etaIndex0 - etaIndex1);
0525         long long deltaEtaLUT_01 = 0;
0526         unsigned int precDeltaEtaLUT_01 = 0;
0527         deltaEtaLUT_01 = m_gtScales->getLUT_DeltaEta(lutName_01, deltaEtaFW_01);
0528         precDeltaEtaLUT_01 = m_gtScales->getPrec_DeltaEta(lutName_01);
0529         ///
0530         int deltaPhiFW_02 = abs(phiIndex0 - phiIndex2);
0531         if (deltaPhiFW_02 >= phiBound)
0532           deltaPhiFW_02 = 2 * phiBound - deltaPhiFW_02;
0533         std::string lutName_02 = lutObj0;
0534         lutName_02 += "-";
0535         lutName_02 += lutObj2;
0536         long long deltaPhiLUT_02 = m_gtScales->getLUT_DeltaPhi(lutName_02, deltaPhiFW_02);
0537         unsigned int precDeltaPhiLUT_02 = m_gtScales->getPrec_DeltaPhi(lutName_02);
0538 
0539         int deltaEtaFW_02 = abs(etaIndex0 - etaIndex2);
0540         long long deltaEtaLUT_02 = 0;
0541         unsigned int precDeltaEtaLUT_02 = 0;
0542         deltaEtaLUT_02 = m_gtScales->getLUT_DeltaEta(lutName_02, deltaEtaFW_02);
0543         precDeltaEtaLUT_02 = m_gtScales->getPrec_DeltaEta(lutName_02);
0544         ///
0545         int deltaPhiFW_12 = abs(phiIndex1 - phiIndex2);
0546         if (deltaPhiFW_12 >= phiBound)
0547           deltaPhiFW_12 = 2 * phiBound - deltaPhiFW_12;
0548         std::string lutName_12 = lutObj1;
0549         lutName_12 += "-";
0550         lutName_12 += lutObj2;
0551         long long deltaPhiLUT_12 = m_gtScales->getLUT_DeltaPhi(lutName_12, deltaPhiFW_12);
0552         unsigned int precDeltaPhiLUT_12 = m_gtScales->getPrec_DeltaPhi(lutName_12);
0553 
0554         int deltaEtaFW_12 = abs(etaIndex1 - etaIndex2);
0555         long long deltaEtaLUT_12 = 0;
0556         unsigned int precDeltaEtaLUT_12 = 0;
0557         deltaEtaLUT_12 = m_gtScales->getLUT_DeltaEta(lutName_12, deltaEtaFW_12);
0558         precDeltaEtaLUT_12 = m_gtScales->getPrec_DeltaEta(lutName_12);
0559         ///
0560 
0561         LogDebug("L1TGlobal") << "### Obj0 phiFW = " << phiIndex0 << " Obj1 phiFW = " << phiIndex1 << "\n"
0562                               << "    DeltaPhiFW = " << deltaPhiFW_01 << "    LUT Name 01= " << lutName_01
0563                               << " Prec = " << precDeltaPhiLUT_01 << "\n"
0564                               << "    LUT Name 02= " << lutName_02 << " Prec = " << precDeltaPhiLUT_02 << "\n"
0565                               << "    LUT Name 12= " << lutName_12 << " Prec = " << precDeltaPhiLUT_12 << "\n"
0566                               << "    DeltaPhiLUT_01 = " << deltaPhiLUT_01 << "\n"
0567                               << "    DeltaPhiLUT_02 = " << deltaPhiLUT_02 << "\n"
0568                               << "    DeltaPhiLUT_12 = " << deltaPhiLUT_12 << "\n"
0569                               << "### Obj0 etaFW = " << etaIndex0 << " Obj1 etaFW = " << etaIndex1 << "\n"
0570                               << "    DeltaEtaFW = " << deltaEtaFW_01 << "    LUT Name 01 = " << lutName_01
0571                               << " Prec 01 = " << precDeltaEtaLUT_01 << "\n"
0572                               << "    LUT Name 02 = " << lutName_02 << " Prec 02 = " << precDeltaEtaLUT_02 << "\n"
0573                               << "    LUT Name 12 = " << lutName_12 << " Prec 12 = " << precDeltaEtaLUT_12 << "\n"
0574                               << "    DeltaEtaLUT_01 = " << deltaEtaLUT_01 << "    DeltaEtaLUT_02 = " << deltaEtaLUT_02
0575                               << "    DeltaEtaLUT_12 = " << deltaEtaLUT_12 << std::endl;
0576 
0577         if (corrPar.corrCutType & 0x9) {
0578           // Invariant mass calculation based for each pair on
0579           // M = sqrt(2*p1*p2(cosh(eta1-eta2) - cos(phi1 - phi2)))
0580           // NOTE: we calculate (1/2)M^2
0581           ///
0582           double cosDeltaPhiPhy_01 = cos(deltaPhiPhy_01);
0583           double coshDeltaEtaPhy_01 = cosh(deltaEtaPhy_01);
0584           double massSqPhy_01 = et0Phy * et1Phy * (coshDeltaEtaPhy_01 - cosDeltaPhiPhy_01);
0585 
0586           long long cosDeltaPhiLUT_01 = m_gtScales->getLUT_DeltaPhi_Cos(lutName_01, deltaPhiFW_01);
0587           unsigned int precCosLUT_01 = m_gtScales->getPrec_DeltaPhi_Cos(lutName_01);
0588 
0589           long long coshDeltaEtaLUT_01;
0590           coshDeltaEtaLUT_01 = m_gtScales->getLUT_DeltaEta_Cosh(lutName_01, deltaEtaFW_01);
0591           unsigned int precCoshLUT_01 = m_gtScales->getPrec_DeltaEta_Cosh(lutName_01);
0592           if (precCoshLUT_01 - precCosLUT_01 != 0)
0593             LogDebug("L1TGlobal") << "Warning: Cos and Cosh LUTs on different precision" << std::endl;
0594 
0595           double cosDeltaPhiPhy_02 = cos(deltaPhiPhy_02);
0596           double coshDeltaEtaPhy_02 = cosh(deltaEtaPhy_02);
0597           if (corrPar.corrCutType & 0x10)
0598             coshDeltaEtaPhy_02 = 1.;
0599           double massSqPhy_02 = et0Phy * et2Phy * (coshDeltaEtaPhy_02 - cosDeltaPhiPhy_02);
0600           long long cosDeltaPhiLUT_02 = m_gtScales->getLUT_DeltaPhi_Cos(lutName_02, deltaPhiFW_02);
0601           unsigned int precCosLUT_02 = m_gtScales->getPrec_DeltaPhi_Cos(lutName_02);
0602           long long coshDeltaEtaLUT_02;
0603           if (corrPar.corrCutType & 0x10) {
0604             coshDeltaEtaLUT_02 = 1 * pow(10, precCosLUT_02);
0605           } else {
0606             coshDeltaEtaLUT_02 = m_gtScales->getLUT_DeltaEta_Cosh(lutName_02, deltaEtaFW_02);
0607             unsigned int precCoshLUT_02 = m_gtScales->getPrec_DeltaEta_Cosh(lutName_02);
0608             if (precCoshLUT_02 - precCosLUT_02 != 0)
0609               LogDebug("L1TGlobal") << "Warning: Cos and Cosh LUTs on different precision" << std::endl;
0610           }
0611 
0612           double cosDeltaPhiPhy_12 = cos(deltaPhiPhy_12);
0613           double coshDeltaEtaPhy_12 = cosh(deltaEtaPhy_12);
0614           if (corrPar.corrCutType & 0x10)
0615             coshDeltaEtaPhy_12 = 1.;
0616           double massSqPhy_12 = et1Phy * et2Phy * (coshDeltaEtaPhy_12 - cosDeltaPhiPhy_12);
0617           long long cosDeltaPhiLUT_12 = m_gtScales->getLUT_DeltaPhi_Cos(lutName_12, deltaPhiFW_12);
0618           unsigned int precCosLUT_12 = m_gtScales->getPrec_DeltaPhi_Cos(lutName_12);
0619           long long coshDeltaEtaLUT_12;
0620           if (corrPar.corrCutType & 0x10) {
0621             coshDeltaEtaLUT_12 = 1 * pow(10, precCosLUT_12);
0622           } else {
0623             coshDeltaEtaLUT_12 = m_gtScales->getLUT_DeltaEta_Cosh(lutName_12, deltaEtaFW_12);
0624             unsigned int precCoshLUT_12 = m_gtScales->getPrec_DeltaEta_Cosh(lutName_12);
0625             if (precCoshLUT_12 - precCosLUT_12 != 0)
0626               LogDebug("L1TGlobal") << "Warning: Cos and Cosh LUTs on different precision" << std::endl;
0627           }
0628 
0629           std::string lutName = lutObj0;
0630           lutName += "-ET";
0631           long long ptObj0 = m_gtScales->getLUT_Pt("Mass_" + lutName, etIndex0);
0632           unsigned int precPtLUTObj0 = m_gtScales->getPrec_Pt("Mass_" + lutName);
0633 
0634           lutName = lutObj1;
0635           lutName += "-ET";
0636           long long ptObj1 = m_gtScales->getLUT_Pt("Mass_" + lutName, etIndex1);
0637           unsigned int precPtLUTObj1 = m_gtScales->getPrec_Pt("Mass_" + lutName);
0638 
0639           lutName = lutObj2;
0640           lutName += "-ET";
0641           long long ptObj2 = m_gtScales->getLUT_Pt("Mass_" + lutName, etIndex2);
0642           unsigned int precPtLUTObj2 = m_gtScales->getPrec_Pt("Mass_" + lutName);
0643 
0644           // Pt and Angles are at different precision
0645           long long massSq_01 = ptObj0 * ptObj1 * (coshDeltaEtaLUT_01 - cosDeltaPhiLUT_01);
0646           long long massSq_02 = ptObj0 * ptObj2 * (coshDeltaEtaLUT_02 - cosDeltaPhiLUT_02);
0647           long long massSq_12 = ptObj1 * ptObj2 * (coshDeltaEtaLUT_12 - cosDeltaPhiLUT_12);
0648 
0649           // Note: There is an assumption here that Cos and Cosh have the same precision
0650           //       unsigned int preShift_01 = precPtLUTObj0 + precPtLUTObj1 + precCosLUT - corrPar.precMassCut;
0651           unsigned int preShift_01 = precPtLUTObj0 + precPtLUTObj1 + precCosLUT_01 - corrPar.precMassCut;
0652           unsigned int preShift_02 = precPtLUTObj0 + precPtLUTObj2 + precCosLUT_02 - corrPar.precMassCut;
0653           unsigned int preShift_12 = precPtLUTObj1 + precPtLUTObj2 + precCosLUT_12 - corrPar.precMassCut;
0654 
0655           LogDebug("L1TGlobal") << "####################################\n";
0656           LogDebug("L1TGlobal") << "    Testing the dimuon invariant mass between the FIRST PAIR 0-1 (" << lutObj0
0657                                 << "," << lutObj1 << ") \n"
0658                                 << "    massSq/2     = " << massSq_01 << "\n"
0659                                 << "    Precision Shift = " << preShift_01 << "\n"
0660                                 << "    massSq   (shift)= " << (massSq_01 / pow(10, preShift_01 + corrPar.precMassCut))
0661                                 << "\n"
0662                                 << "    massSqPhy/2  = " << massSqPhy_01
0663                                 << "  sqrt(|massSq|) = " << sqrt(fabs(2. * massSqPhy_01)) << std::endl;
0664 
0665           LogDebug("L1TGlobal") << "####################################\n";
0666           LogDebug("L1TGlobal") << "    Testing the dimuon invariant mass between the SECOND PAIR 0-2 (" << lutObj0
0667                                 << "," << lutObj2 << ") \n"
0668                                 << "    massSq/2     = " << massSq_02 << "\n"
0669                                 << "    Precision Shift = " << preShift_02 << "\n"
0670                                 << "    massSq   (shift)= " << (massSq_02 / pow(10, preShift_02 + corrPar.precMassCut))
0671                                 << "\n"
0672                                 << "    massSqPhy/2  = " << massSqPhy_02
0673                                 << "  sqrt(|massSq|) = " << sqrt(fabs(2. * massSqPhy_02)) << std::endl;
0674 
0675           LogDebug("L1TGlobal") << "####################################\n";
0676           LogDebug("L1TGlobal") << "    Testing the dimuon invariant mass between the THIRD PAIR 1-2 (" << lutObj1
0677                                 << "," << lutObj2 << ") \n"
0678                                 << "    massSq/2     = " << massSq_12 << "\n"
0679                                 << "    Precision Shift = " << preShift_12 << "\n"
0680                                 << "    massSq   (shift)= " << (massSq_12 / pow(10, preShift_12 + corrPar.precMassCut))
0681                                 << "\n"
0682                                 << "    massSqPhy/2  = " << massSqPhy_12
0683                                 << "  sqrt(|massSq|) = " << sqrt(fabs(2. * massSqPhy_12)) << std::endl;
0684 
0685           LogDebug("L1TGlobal") << "\n ########### THREE-BODY INVARIANT MASS #########################\n";
0686           long long massSq = 0;
0687 
0688           if (preShift_01 == preShift_02 && preShift_01 == preShift_12 && preShift_02 == preShift_12) {
0689             LogDebug("L1TGlobal") << "Check the preshift value: " << preShift_01 << " = " << preShift_02 << " = "
0690                                   << preShift_12 << std::endl;
0691             preShift = preShift_01;
0692           } else {
0693             LogDebug("L1TGlobal")
0694                 << "Preshift values considered for the sum of the dimuon invariant masses are different!" << std::endl;
0695           }
0696 
0697           if ((massSq_01 != massSq_02) && (massSq_01 != massSq_12) && (massSq_02 != massSq_12)) {
0698             massSq = massSq_01 + massSq_02 + massSq_12;
0699             LogDebug("L1TGlobal") << "massSq = " << massSq << std::endl;
0700           } else {
0701             LogDebug("L1TGlobal") << "Same pair of muons considered, three-body invariant mass do not computed"
0702                                   << std::endl;
0703           }
0704 
0705           if (massSq >= 0 && massSq >= (long long)(corrPar.minMassCutValue * pow(10, preShift)) &&
0706               massSq <= (long long)(corrPar.maxMassCutValue * pow(10, preShift))) {
0707             LogDebug("L1TGlobal") << "    Passed Invariant Mass Cut ["
0708                                   << (long long)(corrPar.minMassCutValue * pow(10, preShift)) << ","
0709                                   << (long long)(corrPar.maxMassCutValue * pow(10, preShift)) << "]" << std::endl;
0710             reqResult = true;
0711           } else {
0712             LogDebug("L1TGlobal") << "    Failed Invariant Mass Cut ["
0713                                   << (long long)(corrPar.minMassCutValue * pow(10, preShift)) << ","
0714                                   << (long long)(corrPar.maxMassCutValue * pow(10, preShift)) << "]" << std::endl;
0715             reqResult = false;
0716           }
0717         }
0718 
0719         if (reqResult && chrgCorrel) {
0720           condResult = true;
0721           (combinationsInCond()).push_back(objectsInComb);
0722         }
0723 
0724       }  //end loop over third leg
0725     }    //end loop over second leg
0726   }      //end loop over first leg
0727 
0728   if (m_verbosity && condResult) {
0729     LogDebug("L1TGlobal") << " pass(es) the correlation condition.\n" << std::endl;
0730   }
0731   return condResult;
0732 }
0733 
0734 /**
0735  * checkObjectParameter - Compare a single particle with a numbered condition
0736  *
0737  * @param iCondition: The number of the condition.
0738  * @param cand: The candidate to compare.
0739  *
0740  * @return: The result of the comparison (false if a condition does not exist)
0741  */
0742 
0743 const bool l1t::CorrThreeBodyCondition::checkObjectParameter(const int iCondition, const l1t::L1Candidate& cand) const {
0744   return true;
0745 }
0746 
0747 void l1t::CorrThreeBodyCondition::print(std::ostream& myCout) const {
0748   myCout << "Dummy Print for CorrThreeBodyCondition" << std::endl;
0749   m_gtCorrelationThreeBodyTemplate->print(myCout);
0750 
0751   ConditionEvaluation::print(myCout);
0752 }