Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:53

0001 /*
0002  * \file L1TRate.cc
0003  *
0004  * \author J. Pela, P. Musella
0005  *
0006  */
0007 
0008 // L1TMonitor includes
0009 #include "DQM/L1TMonitor/interface/L1TRate.h"
0010 #include "DQM/L1TMonitor/interface/L1TOMDSHelper.h"
0011 
0012 #include "DQMServices/Core/interface/DQMStore.h"
0013 
0014 #include "DataFormats/Common/interface/ConditionsInEdm.h"  // Parameters associated to Run, LS and Event
0015 #include "DataFormats/Luminosity/interface/LumiDetails.h"  // Luminosity Information
0016 #include "DataFormats/Luminosity/interface/LumiSummary.h"  // Luminosity Information
0017 
0018 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
0019 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenuFwd.h"
0020 #include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
0021 #include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"             // L1Gt - Masks
0022 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h"  // L1Gt - Masks
0023 #include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
0024 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"
0025 
0026 #include "TList.h"
0027 
0028 using namespace edm;
0029 using namespace std;
0030 
0031 //_____________________________________________________________________
0032 L1TRate::L1TRate(const ParameterSet& ps)
0033     : m_menuToken(esConsumes<edm::Transition::BeginRun>()),
0034       m_l1GtPfAlgoToken(esConsumes<edm::Transition::BeginRun>()),
0035       m_helperTokens(L1TMenuHelper::consumes<edm::Transition::BeginRun>(consumesCollector())),
0036       m_l1GtUtils(ps, consumesCollector(), false, *this) {
0037   m_maxNbins = 2500;  // Maximum LS for each run (for binning purposes)
0038   m_parameters = ps;
0039 
0040   // Mapping parameter input variables
0041   m_scalersSource_colLScal =
0042       consumes<LumiScalersCollection>(m_parameters.getParameter<InputTag>("inputTagScalersResults"));
0043   m_scalersSource_triggerScalers =
0044       consumes<Level1TriggerScalersCollection>(m_parameters.getParameter<InputTag>("inputTagScalersResults"));
0045   m_l1GtDataDaqInputTag =
0046       consumes<L1GlobalTriggerReadoutRecord>(m_parameters.getParameter<InputTag>("inputTagL1GtDataDaq"));
0047   m_verbose = m_parameters.getUntrackedParameter<bool>("verbose", false);
0048   m_refPrescaleSet = m_parameters.getParameter<int>("refPrescaleSet");
0049   m_lsShiftGTRates = m_parameters.getUntrackedParameter<int>("lsShiftGTRates", 0);
0050 
0051   // Getting which categories to monitor
0052   ParameterSet Categories = ps.getParameter<ParameterSet>("categories");
0053   m_inputCategories["Mu"] = Categories.getUntrackedParameter<bool>("Mu");
0054   m_inputCategories["EG"] = Categories.getUntrackedParameter<bool>("EG");
0055   m_inputCategories["IsoEG"] = Categories.getUntrackedParameter<bool>("IsoEG");
0056   m_inputCategories["Jet"] = Categories.getUntrackedParameter<bool>("Jet");
0057   m_inputCategories["CenJet"] = Categories.getUntrackedParameter<bool>("CenJet");
0058   m_inputCategories["ForJet"] = Categories.getUntrackedParameter<bool>("ForJet");
0059   m_inputCategories["TauJet"] = Categories.getUntrackedParameter<bool>("TauJet");
0060   m_inputCategories["ETM"] = Categories.getUntrackedParameter<bool>("ETM");
0061   m_inputCategories["ETT"] = Categories.getUntrackedParameter<bool>("ETT");
0062   m_inputCategories["HTT"] = Categories.getUntrackedParameter<bool>("HTT");
0063   m_inputCategories["HTM"] = Categories.getUntrackedParameter<bool>("HTM");
0064 
0065   // What to do if we want our output to be saved to a external file
0066   m_outputFile = ps.getUntrackedParameter<string>("outputFile", "");
0067 
0068   if (!m_outputFile.empty()) {
0069     cout << "L1T Monitoring histograms will be saved to " << m_outputFile.c_str() << endl;
0070   }
0071 
0072   bool disable = ps.getUntrackedParameter<bool>("disableROOToutput", false);
0073   if (disable) {
0074     m_outputFile = "";
0075   }
0076 }
0077 
0078 //_____________________________________________________________________
0079 L1TRate::~L1TRate() {}
0080 
0081 //_____________________________________________________________________
0082 // BeginRun
0083 //_____________________________________________________________________
0084 void L1TRate::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup& iSetup) {
0085   const L1GtTriggerMenu& menu = iSetup.getData(m_menuToken);
0086   const L1GtPrescaleFactors& l1GtPfAlgo = iSetup.getData(m_l1GtPfAlgoToken);
0087 
0088   // Initializing DQM Monitor Elements
0089   ibooker.setCurrentFolder("L1T/L1TRate");
0090   m_ErrorMonitor = ibooker.book1D("ErrorMonitor", "ErrorMonitor", 5, 0, 5);
0091   m_ErrorMonitor->setBinLabel(1, "WARNING_DB_CONN_FAILED");         // Errors from L1TOMDSHelper
0092   m_ErrorMonitor->setBinLabel(2, "WARNING_DB_QUERY_FAILED");        // Errors from L1TOMDSHelper
0093   m_ErrorMonitor->setBinLabel(3, "WARNING_DB_INCORRECT_NBUNCHES");  // Errors from L1TOMDSHelper
0094   m_ErrorMonitor->setBinLabel(4, "WARNING_PY_MISSING_FIT");
0095   m_ErrorMonitor->setBinLabel(5, "UNKNOWN");
0096 
0097   // Retriving the list of prescale sets
0098   m_listsPrescaleFactors = &(l1GtPfAlgo.gtPrescaleFactors());
0099 
0100   // Getting Lowest Prescale Single Object Triggers from the menu
0101   L1TMenuHelper myMenuHelper = L1TMenuHelper(iSetup, m_helperTokens);
0102   m_l1GtUtils.retrieveL1EventSetup(iSetup);
0103   m_selectedTriggers = myMenuHelper.getLUSOTrigger(m_inputCategories, m_refPrescaleSet, m_l1GtUtils);
0104 
0105   //-> Getting template fits for the algLo cross sections
0106   int srcAlgoXSecFit = m_parameters.getParameter<int>("srcAlgoXSecFit");
0107   if (srcAlgoXSecFit == 0) {
0108     getXSexFitsOMDS(m_parameters);
0109   } else if (srcAlgoXSecFit == 1) {
0110     getXSexFitsPython(m_parameters);
0111   }
0112 
0113   for (const auto& algo : menu.gtAlgorithmMap()) {
0114     m_algoBit[algo.second.algoAlias()] = algo.second.algoBitNumber();
0115   }
0116 
0117   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
0118   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
0119 
0120   // Initializing DQM Monitor Elements
0121   for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
0122     TString tCategory = (*i).first;
0123     TString tTrigger = (*i).second;
0124 
0125     TString tErrorMessage = "";
0126     TF1* tTestFunction;
0127 
0128     if (tTrigger != "Undefined" && m_templateFunctions.find(tTrigger) != m_templateFunctions.end()) {
0129       tTestFunction = m_templateFunctions[tTrigger];
0130     } else if (tTrigger == "Undefined") {
0131       TString tFunc = "-1";
0132       tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
0133     } else if (m_templateFunctions.find(tTrigger) == m_templateFunctions.end()) {
0134       TString tFunc = "-1";
0135       tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
0136       tErrorMessage = " (Undefined Test Function)";
0137     } else {
0138       TString tFunc = "-1";
0139       tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
0140     }
0141 
0142     if (tTrigger != "Undefined") {
0143       if (myMenuHelper.getPrescaleByAlias(tCategory, tTrigger) != 1) {
0144         tErrorMessage += " WARNING: Default Prescale = ";
0145         tErrorMessage += myMenuHelper.getPrescaleByAlias(tCategory, tTrigger);
0146       }
0147 
0148       if (tCategory == "Mu" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 4294967295) {
0149         tErrorMessage += " WARNING: Eta Range = ";
0150         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
0151       } else if (tCategory == "EG" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 32639) {
0152         tErrorMessage += " WARNING: Eta Range = ";
0153         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
0154       } else if (tCategory == "IsoEG" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 32639) {
0155         tErrorMessage += " WARNING: Eta Range = ";
0156         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
0157       }
0158 
0159       if (tCategory == "Mu" && myMenuHelper.getQualityAlias(tCategory, tTrigger) != 240) {
0160         tErrorMessage += " WARNING: Quality = ";
0161         tErrorMessage += myMenuHelper.getQualityAlias(tCategory, tTrigger);
0162       }
0163     }
0164 
0165     ibooker.setCurrentFolder("L1T/L1TRate/TriggerCrossSections");
0166     m_xSecVsInstLumi[tTrigger] = ibooker.bookProfile(tCategory,
0167                                                      "Cross Sec. vs Inst. Lumi Algo: " + tTrigger + tErrorMessage,
0168                                                      m_maxNbins,
0169                                                      minInstantLuminosity,
0170                                                      maxInstantLuminosity,
0171                                                      0,
0172                                                      500);
0173     m_xSecVsInstLumi[tTrigger]->setAxisTitle("Instantaneous Luminosity [10^{30}cm^{-2}s^{-1}]", 1);
0174     m_xSecVsInstLumi[tTrigger]->setAxisTitle("Algorithm #sigma [#mu b]", 2);
0175     m_xSecVsInstLumi[tTrigger]->getTProfile()->GetListOfFunctions()->Add(tTestFunction);
0176     m_xSecVsInstLumi[tTrigger]->getTProfile()->SetMarkerStyle(23);
0177 
0178     ibooker.setCurrentFolder("L1T/L1TRate/Certification");
0179     m_xSecObservedToExpected[tTrigger] =
0180         ibooker.book1D(tCategory, "Algo: " + tTrigger + tErrorMessage, m_maxNbins, -0.5, double(m_maxNbins) - 0.5);
0181     m_xSecObservedToExpected[tTrigger]->setAxisTitle("Lumi Section", 1);
0182     m_xSecObservedToExpected[tTrigger]->setAxisTitle("#sigma_{obs} / #sigma_{exp}", 2);
0183   }
0184 }
0185 
0186 void L1TRate::dqmBeginRun(edm::Run const&, edm::EventSetup const&) {
0187   //
0188   if (m_verbose) {
0189     cout << "[L1TRate:] Called beginRun." << endl;
0190   }
0191 }
0192 //_____________________________________________________________________
0193 void L1TRate::beginLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
0194   if (m_verbose) {
0195     cout << "[L1TRate:] Called beginLuminosityBlock at LS=" << lumiBlock.id().luminosityBlock() << endl;
0196   }
0197 }
0198 
0199 //_____________________________________________________________________
0200 void L1TRate::endLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
0201   int eventLS = lumiBlock.id().luminosityBlock();
0202   if (m_verbose) {
0203     cout << "[L1TRate:] Called endLuminosityBlock at LS=" << eventLS << endl;
0204   }
0205 
0206   // We can certify LS -1 since we should have available:
0207   // gt rates: (current LS)-1
0208   // prescale: current LS
0209   // lumi    : current LS
0210   //eventLS--;
0211 
0212   // Checking if all necessary quantities are defined for our calculations
0213   bool isDefRate, isDefLumi, isDefPrescaleIndex;
0214   map<TString, double>* rates = nullptr;
0215   double lumi = 0;
0216   int prescalesIndex = 0;
0217 
0218   // Reseting MonitorElements so we can refill them
0219   for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
0220     string tTrigger = (*i).second;
0221     m_xSecObservedToExpected[tTrigger]->getTH1()->Reset("ICE");
0222     m_xSecVsInstLumi[tTrigger]->getTH1()->Reset("ICE");
0223   }
0224 
0225   for (map<int, map<TString, double> >::iterator i = m_lsRates.begin(); i != m_lsRates.end(); i++) {
0226     unsigned int ls = (*i).first;
0227     rates = &(*i).second;
0228     isDefRate = true;
0229 
0230     if (m_lsLuminosity.find(ls) == m_lsLuminosity.end()) {
0231       isDefLumi = false;
0232     } else {
0233       isDefLumi = true;
0234       lumi = m_lsLuminosity[ls];
0235     }
0236 
0237     if (m_lsPrescaleIndex.find(ls) == m_lsPrescaleIndex.end()) {
0238       isDefPrescaleIndex = false;
0239     } else {
0240       isDefPrescaleIndex = true;
0241       prescalesIndex = m_lsPrescaleIndex[ls];
0242     }
0243 
0244     if (isDefRate && isDefLumi && isDefPrescaleIndex) {
0245       const vector<int>& currentPrescaleFactors = (*m_listsPrescaleFactors).at(prescalesIndex);
0246 
0247       for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
0248         string tTrigger = (*i).second;
0249         TF1* tTestFunction = (TF1*)m_xSecVsInstLumi[tTrigger]->getTProfile()->GetListOfFunctions()->First();
0250 
0251         // If trigger name is defined we get the rate fit parameters
0252         if (tTrigger != "Undefined") {
0253           unsigned int trigBit = m_algoBit[tTrigger];
0254           double trigPrescale = currentPrescaleFactors[trigBit];
0255           double trigRate = (*rates)[tTrigger];
0256 
0257           if (lumi != 0 && trigPrescale != 0 && trigRate != 0) {
0258             double AlgoXSec = (trigPrescale * trigRate) / lumi;
0259             double TemplateFunctionValue = tTestFunction->Eval(lumi);
0260 
0261             // Checking against Template function
0262             int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
0263             m_xSecObservedToExpected[tTrigger]->setBinContent(ibin, AlgoXSec / TemplateFunctionValue);
0264             m_xSecVsInstLumi[tTrigger]->Fill(lumi, AlgoXSec);
0265 
0266             if (m_verbose) {
0267               cout << "[L1TRate:] ls=" << ls << " Algo=" << tTrigger << " XSec=" << AlgoXSec
0268                    << " Test=" << AlgoXSec / TemplateFunctionValue << endl;
0269             }
0270 
0271           } else {
0272             int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
0273             m_xSecObservedToExpected[tTrigger]->setBinContent(ibin, 0.000001);
0274             if (m_verbose) {
0275               cout << "[L1TRate:] Algo=" << tTrigger << " XSec=Failed" << endl;
0276             }
0277           }
0278         }
0279       }
0280     }
0281   }
0282 }
0283 
0284 //_____________________________________________________________________
0285 void L1TRate::analyze(const Event& iEvent, const EventSetup& eventSetup) {
0286   edm::Handle<L1GlobalTriggerReadoutRecord> gtReadoutRecordData;
0287   edm::Handle<Level1TriggerScalersCollection> triggerScalers;
0288   edm::Handle<LumiScalersCollection> colLScal;
0289 
0290   iEvent.getByToken(m_l1GtDataDaqInputTag, gtReadoutRecordData);
0291   iEvent.getByToken(m_scalersSource_colLScal, colLScal);
0292   iEvent.getByToken(m_scalersSource_triggerScalers, triggerScalers);
0293 
0294   // Integers
0295   int EventRun = iEvent.id().run();
0296   unsigned int eventLS = iEvent.id().luminosityBlock();
0297 
0298   // Getting the trigger trigger rates from GT and buffering it
0299   if (triggerScalers.isValid()) {
0300     Level1TriggerScalersCollection::const_iterator itL1TScalers = triggerScalers->begin();
0301     Level1TriggerRates trigRates(*itL1TScalers, EventRun);
0302 
0303     int gtLS = (*itL1TScalers).lumiSegmentNr() + m_lsShiftGTRates;
0304 
0305     // If we haven't got the data from this LS yet get it
0306     if (m_lsRates.find(gtLS) == m_lsRates.end()) {
0307       if (m_verbose) {
0308         cout << "[L1TRate:] Buffering GT Rates for LS=" << gtLS << endl;
0309       }
0310       map<TString, double> bufferRate;
0311 
0312       // Buffer the rate informations for all selected bits
0313       for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
0314         string tTrigger = (*i).second;
0315 
0316         // If trigger name is defined we store the rate
0317         if (tTrigger != "Undefined") {
0318           unsigned int trigBit = m_algoBit[tTrigger];
0319           double trigRate = trigRates.gtAlgoCountsRate()[trigBit];
0320 
0321           bufferRate[tTrigger] = trigRate;
0322         }
0323       }
0324       m_lsRates[gtLS] = bufferRate;
0325     }
0326   }
0327 
0328   // Getting from the SCAL the luminosity information and buffering it
0329   if (colLScal.isValid() && !colLScal->empty()) {
0330     LumiScalersCollection::const_iterator itLScal = colLScal->begin();
0331     unsigned int scalLS = itLScal->sectionNumber();
0332 
0333     // If we haven't got the data from this SCAL LS yet get it
0334     if (m_lsLuminosity.find(scalLS) == m_lsLuminosity.end()) {
0335       if (m_verbose) {
0336         cout << "[L1TRate:] Buffering SCAL-HF Lumi for LS=" << scalLS << endl;
0337       }
0338       double instLumi = itLScal->instantLumi();                  // Getting Instant Lumi from HF (via SCAL)
0339       double deadTimeNormHF = itLScal->deadTimeNormalization();  // Getting Dead Time Normalization from HF (via SCAL)
0340 
0341       // If HF Dead Time Corrections is requested we apply it
0342       // NOTE: By default this is assumed false since for now WbM fits do NOT assume this correction
0343       if (m_parameters.getUntrackedParameter<bool>("useHFDeadTimeNormalization", false)) {
0344         // Protecting for deadtime = 0
0345         if (deadTimeNormHF == 0) {
0346           instLumi = 0;
0347         } else {
0348           instLumi = instLumi / deadTimeNormHF;
0349         }
0350       }
0351       // Buffering the luminosity information
0352       m_lsLuminosity[scalLS] = instLumi;
0353     }
0354   }
0355 
0356   // Getting the prescale index used when this event was triggered
0357   if (gtReadoutRecordData.isValid()) {
0358     // If we haven't got the data from this LS yet get it
0359     if (m_lsPrescaleIndex.find(eventLS) == m_lsPrescaleIndex.end()) {
0360       if (m_verbose) {
0361         cout << "[L1TRate:] Buffering Prescale Index for LS=" << eventLS << endl;
0362       }
0363 
0364       // Getting Final Decision Logic (FDL) Data from GT
0365       const vector<L1GtFdlWord>& gtFdlVectorData = gtReadoutRecordData->gtFdlVector();
0366 
0367       // Getting the index for the fdl data for this event
0368       int indexFDL = 0;
0369       for (unsigned int i = 0; i < gtFdlVectorData.size(); i++) {
0370         if (gtFdlVectorData[i].bxInEvent() == 0) {
0371           indexFDL = i;
0372           break;
0373         }
0374       }
0375 
0376       int CurrentPrescalesIndex = gtFdlVectorData[indexFDL].gtPrescaleFactorIndexAlgo();
0377       m_lsPrescaleIndex[eventLS] = CurrentPrescalesIndex;
0378     }
0379   }
0380 }
0381 
0382 //_____________________________________________________________________
0383 // function: getXSexFitsOMDS
0384 // Imputs:
0385 //   * const edm::ParameterSet& ps = ParameterSet contaning necessary
0386 //     information for the OMDS data extraction
0387 // Outputs:
0388 //   * int error = Number of algos where you did not find a
0389 //     corresponding fit
0390 //_____________________________________________________________________
0391 bool L1TRate::getXSexFitsOMDS(const edm::ParameterSet& ps) {
0392   bool noError = true;
0393 
0394   string oracleDB = ps.getParameter<string>("oracleDB");
0395   string pathCondDB = ps.getParameter<string>("pathCondDB");
0396 
0397   L1TOMDSHelper myOMDSHelper;
0398   int conError;
0399   myOMDSHelper.connect(oracleDB, pathCondDB, conError);
0400 
0401   map<string, WbMTriggerXSecFit> wbmFits;
0402 
0403   if (conError == L1TOMDSHelper::NO_ERROR) {
0404     int errorRetrive;
0405     wbmFits = myOMDSHelper.getWbMAlgoXsecFits(errorRetrive);
0406 
0407     // Filling errors if they exist
0408     if (errorRetrive != L1TOMDSHelper::NO_ERROR) {
0409       noError = false;
0410       string eName = myOMDSHelper.enumToStringError(errorRetrive);
0411       m_ErrorMonitor->Fill(eName);
0412     }
0413   } else {
0414     noError = false;
0415     string eName = myOMDSHelper.enumToStringError(conError);
0416     m_ErrorMonitor->Fill(eName);
0417   }
0418 
0419   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
0420   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
0421 
0422   // Getting rate fit parameters for all input triggers
0423   for (map<string, string>::const_iterator a = m_selectedTriggers.begin(); a != m_selectedTriggers.end(); a++) {
0424     string tTrigger = (*a).second;
0425 
0426     // If trigger name is defined we get the rate fit parameters
0427     if (tTrigger != "Undefined") {
0428       if (wbmFits.find(tTrigger) != wbmFits.end()) {
0429         WbMTriggerXSecFit tWbMParameters = wbmFits[tTrigger];
0430 
0431         vector<double> tParameters;
0432         tParameters.push_back(tWbMParameters.pm1);
0433         tParameters.push_back(tWbMParameters.p0);
0434         tParameters.push_back(tWbMParameters.p1);
0435         tParameters.push_back(tWbMParameters.p2);
0436 
0437         // Retriving and populating the m_templateFunctions array
0438         m_templateFunctions[tTrigger] = new TF1("FitParametrization_" + tWbMParameters.bitName,
0439                                                 tWbMParameters.fitFunction,
0440                                                 minInstantLuminosity,
0441                                                 maxInstantLuminosity);
0442         m_templateFunctions[tTrigger]->SetParameters(&tParameters[0]);
0443         m_templateFunctions[tTrigger]->SetLineWidth(1);
0444         m_templateFunctions[tTrigger]->SetLineColor(kRed);
0445 
0446       } else {
0447         noError = false;
0448       }
0449     }
0450   }
0451 
0452   return noError;
0453 }
0454 
0455 //_____________________________________________________________________
0456 // function: getXSexFitsPython
0457 // Imputs:
0458 //   * const edm::ParameterSet& ps = ParameterSet contaning the fit
0459 //     functions and parameters for the selected triggers
0460 // Outputs:
0461 //   * int error = Number of algos where you did not find a
0462 //     corresponding fit
0463 //_____________________________________________________________________
0464 bool L1TRate::getXSexFitsPython(const edm::ParameterSet& ps) {
0465   // error meaning
0466   bool noError = true;
0467 
0468   // Getting fit parameters
0469   std::vector<edm::ParameterSet> m_fitParameters = ps.getParameter<vector<ParameterSet> >("fitParameters");
0470 
0471   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
0472   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
0473 
0474   // Getting rate fit parameters for all input triggers
0475   for (map<string, string>::const_iterator a = m_selectedTriggers.begin(); a != m_selectedTriggers.end(); a++) {
0476     string tTrigger = (*a).second;
0477 
0478     // If trigger name is defined we get the rate fit parameters
0479     if (tTrigger != "Undefined") {
0480       bool foundFit = false;
0481 
0482       for (unsigned int b = 0; b < m_fitParameters.size(); b++) {
0483         if (tTrigger == m_fitParameters[b].getParameter<string>("AlgoName")) {
0484           TString tAlgoName = m_fitParameters[b].getParameter<string>("AlgoName");
0485           TString tTemplateFunction = m_fitParameters[b].getParameter<string>("TemplateFunction");
0486           vector<double> tParameters = m_fitParameters[b].getParameter<vector<double> >("Parameters");
0487 
0488           // Retriving and populating the m_templateFunctions array
0489           m_templateFunctions[tTrigger] =
0490               new TF1("FitParametrization_" + tAlgoName, tTemplateFunction, minInstantLuminosity, maxInstantLuminosity);
0491           m_templateFunctions[tTrigger]->SetParameters(&tParameters[0]);
0492           m_templateFunctions[tTrigger]->SetLineWidth(1);
0493           m_templateFunctions[tTrigger]->SetLineColor(kRed);
0494 
0495           foundFit = true;
0496           break;
0497         }
0498       }
0499       if (!foundFit) {
0500         noError = false;
0501         string eName = "WARNING_PY_MISSING_FIT";
0502         m_ErrorMonitor->Fill(eName);
0503       }
0504     }
0505   }
0506 
0507   return noError;
0508 }