Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508
/*
 * \file L1TRate.cc
 *
 * \author J. Pela, P. Musella
 *
 */

// L1TMonitor includes
#include "DQM/L1TMonitor/interface/L1TRate.h"
#include "DQM/L1TMonitor/interface/L1TOMDSHelper.h"

#include "DQMServices/Core/interface/DQMStore.h"

#include "DataFormats/Common/interface/ConditionsInEdm.h"  // Parameters associated to Run, LS and Event
#include "DataFormats/Luminosity/interface/LumiDetails.h"  // Luminosity Information
#include "DataFormats/Luminosity/interface/LumiSummary.h"  // Luminosity Information

#include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
#include "CondFormats/L1TObjects/interface/L1GtTriggerMenuFwd.h"
#include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
#include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"             // L1Gt - Masks
#include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h"  // L1Gt - Masks
#include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
#include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"

#include "TList.h"

using namespace edm;
using namespace std;

//_____________________________________________________________________
L1TRate::L1TRate(const ParameterSet& ps)
    : m_menuToken(esConsumes<edm::Transition::BeginRun>()),
      m_l1GtPfAlgoToken(esConsumes<edm::Transition::BeginRun>()),
      m_helperTokens(L1TMenuHelper::consumes<edm::Transition::BeginRun>(consumesCollector())),
      m_l1GtUtils(ps, consumesCollector(), false, *this) {
  m_maxNbins = 2500;  // Maximum LS for each run (for binning purposes)
  m_parameters = ps;

  // Mapping parameter input variables
  m_scalersSource_colLScal =
      consumes<LumiScalersCollection>(m_parameters.getParameter<InputTag>("inputTagScalersResults"));
  m_scalersSource_triggerScalers =
      consumes<Level1TriggerScalersCollection>(m_parameters.getParameter<InputTag>("inputTagScalersResults"));
  m_l1GtDataDaqInputTag =
      consumes<L1GlobalTriggerReadoutRecord>(m_parameters.getParameter<InputTag>("inputTagL1GtDataDaq"));
  m_verbose = m_parameters.getUntrackedParameter<bool>("verbose", false);
  m_refPrescaleSet = m_parameters.getParameter<int>("refPrescaleSet");
  m_lsShiftGTRates = m_parameters.getUntrackedParameter<int>("lsShiftGTRates", 0);

  // Getting which categories to monitor
  ParameterSet Categories = ps.getParameter<ParameterSet>("categories");
  m_inputCategories["Mu"] = Categories.getUntrackedParameter<bool>("Mu");
  m_inputCategories["EG"] = Categories.getUntrackedParameter<bool>("EG");
  m_inputCategories["IsoEG"] = Categories.getUntrackedParameter<bool>("IsoEG");
  m_inputCategories["Jet"] = Categories.getUntrackedParameter<bool>("Jet");
  m_inputCategories["CenJet"] = Categories.getUntrackedParameter<bool>("CenJet");
  m_inputCategories["ForJet"] = Categories.getUntrackedParameter<bool>("ForJet");
  m_inputCategories["TauJet"] = Categories.getUntrackedParameter<bool>("TauJet");
  m_inputCategories["ETM"] = Categories.getUntrackedParameter<bool>("ETM");
  m_inputCategories["ETT"] = Categories.getUntrackedParameter<bool>("ETT");
  m_inputCategories["HTT"] = Categories.getUntrackedParameter<bool>("HTT");
  m_inputCategories["HTM"] = Categories.getUntrackedParameter<bool>("HTM");

  // What to do if we want our output to be saved to a external file
  m_outputFile = ps.getUntrackedParameter<string>("outputFile", "");

  if (!m_outputFile.empty()) {
    cout << "L1T Monitoring histograms will be saved to " << m_outputFile.c_str() << endl;
  }

  bool disable = ps.getUntrackedParameter<bool>("disableROOToutput", false);
  if (disable) {
    m_outputFile = "";
  }
}

//_____________________________________________________________________
L1TRate::~L1TRate() {}

//_____________________________________________________________________
// BeginRun
//_____________________________________________________________________
void L1TRate::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup& iSetup) {
  const L1GtTriggerMenu& menu = iSetup.getData(m_menuToken);
  const L1GtPrescaleFactors& l1GtPfAlgo = iSetup.getData(m_l1GtPfAlgoToken);

  // Initializing DQM Monitor Elements
  ibooker.setCurrentFolder("L1T/L1TRate");
  m_ErrorMonitor = ibooker.book1D("ErrorMonitor", "ErrorMonitor", 5, 0, 5);
  m_ErrorMonitor->setBinLabel(1, "WARNING_DB_CONN_FAILED");         // Errors from L1TOMDSHelper
  m_ErrorMonitor->setBinLabel(2, "WARNING_DB_QUERY_FAILED");        // Errors from L1TOMDSHelper
  m_ErrorMonitor->setBinLabel(3, "WARNING_DB_INCORRECT_NBUNCHES");  // Errors from L1TOMDSHelper
  m_ErrorMonitor->setBinLabel(4, "WARNING_PY_MISSING_FIT");
  m_ErrorMonitor->setBinLabel(5, "UNKNOWN");

  // Retriving the list of prescale sets
  m_listsPrescaleFactors = &(l1GtPfAlgo.gtPrescaleFactors());

  // Getting Lowest Prescale Single Object Triggers from the menu
  L1TMenuHelper myMenuHelper = L1TMenuHelper(iSetup, m_helperTokens);
  m_l1GtUtils.retrieveL1EventSetup(iSetup);
  m_selectedTriggers = myMenuHelper.getLUSOTrigger(m_inputCategories, m_refPrescaleSet, m_l1GtUtils);

  //-> Getting template fits for the algLo cross sections
  int srcAlgoXSecFit = m_parameters.getParameter<int>("srcAlgoXSecFit");
  if (srcAlgoXSecFit == 0) {
    getXSexFitsOMDS(m_parameters);
  } else if (srcAlgoXSecFit == 1) {
    getXSexFitsPython(m_parameters);
  }

  for (const auto& algo : menu.gtAlgorithmMap()) {
    m_algoBit[algo.second.algoAlias()] = algo.second.algoBitNumber();
  }

  double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
  double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");

  // Initializing DQM Monitor Elements
  for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
    TString tCategory = (*i).first;
    TString tTrigger = (*i).second;

    TString tErrorMessage = "";
    TF1* tTestFunction;

    if (tTrigger != "Undefined" && m_templateFunctions.find(tTrigger) != m_templateFunctions.end()) {
      tTestFunction = m_templateFunctions[tTrigger];
    } else if (tTrigger == "Undefined") {
      TString tFunc = "-1";
      tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
    } else if (m_templateFunctions.find(tTrigger) == m_templateFunctions.end()) {
      TString tFunc = "-1";
      tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
      tErrorMessage = " (Undefined Test Function)";
    } else {
      TString tFunc = "-1";
      tTestFunction = new TF1("FitParametrization_" + tTrigger, tFunc, 0, double(m_maxNbins) - 0.5);
    }

    if (tTrigger != "Undefined") {
      if (myMenuHelper.getPrescaleByAlias(tCategory, tTrigger) != 1) {
        tErrorMessage += " WARNING: Default Prescale = ";
        tErrorMessage += myMenuHelper.getPrescaleByAlias(tCategory, tTrigger);
      }

      if (tCategory == "Mu" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 4294967295) {
        tErrorMessage += " WARNING: Eta Range = ";
        tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
      } else if (tCategory == "EG" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 32639) {
        tErrorMessage += " WARNING: Eta Range = ";
        tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
      } else if (tCategory == "IsoEG" && myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger) != 32639) {
        tErrorMessage += " WARNING: Eta Range = ";
        tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory, tTrigger);
      }

      if (tCategory == "Mu" && myMenuHelper.getQualityAlias(tCategory, tTrigger) != 240) {
        tErrorMessage += " WARNING: Quality = ";
        tErrorMessage += myMenuHelper.getQualityAlias(tCategory, tTrigger);
      }
    }

    ibooker.setCurrentFolder("L1T/L1TRate/TriggerCrossSections");
    m_xSecVsInstLumi[tTrigger] = ibooker.bookProfile(tCategory,
                                                     "Cross Sec. vs Inst. Lumi Algo: " + tTrigger + tErrorMessage,
                                                     m_maxNbins,
                                                     minInstantLuminosity,
                                                     maxInstantLuminosity,
                                                     0,
                                                     500);
    m_xSecVsInstLumi[tTrigger]->setAxisTitle("Instantaneous Luminosity [10^{30}cm^{-2}s^{-1}]", 1);
    m_xSecVsInstLumi[tTrigger]->setAxisTitle("Algorithm #sigma [#mu b]", 2);
    m_xSecVsInstLumi[tTrigger]->getTProfile()->GetListOfFunctions()->Add(tTestFunction);
    m_xSecVsInstLumi[tTrigger]->getTProfile()->SetMarkerStyle(23);

    ibooker.setCurrentFolder("L1T/L1TRate/Certification");
    m_xSecObservedToExpected[tTrigger] =
        ibooker.book1D(tCategory, "Algo: " + tTrigger + tErrorMessage, m_maxNbins, -0.5, double(m_maxNbins) - 0.5);
    m_xSecObservedToExpected[tTrigger]->setAxisTitle("Lumi Section", 1);
    m_xSecObservedToExpected[tTrigger]->setAxisTitle("#sigma_{obs} / #sigma_{exp}", 2);
  }
}

void L1TRate::dqmBeginRun(edm::Run const&, edm::EventSetup const&) {
  //
  if (m_verbose) {
    cout << "[L1TRate:] Called beginRun." << endl;
  }
}
//_____________________________________________________________________
void L1TRate::beginLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
  if (m_verbose) {
    cout << "[L1TRate:] Called beginLuminosityBlock at LS=" << lumiBlock.id().luminosityBlock() << endl;
  }
}

//_____________________________________________________________________
void L1TRate::endLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
  int eventLS = lumiBlock.id().luminosityBlock();
  if (m_verbose) {
    cout << "[L1TRate:] Called endLuminosityBlock at LS=" << eventLS << endl;
  }

  // We can certify LS -1 since we should have available:
  // gt rates: (current LS)-1
  // prescale: current LS
  // lumi    : current LS
  //eventLS--;

  // Checking if all necessary quantities are defined for our calculations
  bool isDefRate, isDefLumi, isDefPrescaleIndex;
  map<TString, double>* rates = nullptr;
  double lumi = 0;
  int prescalesIndex = 0;

  // Reseting MonitorElements so we can refill them
  for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
    string tTrigger = (*i).second;
    m_xSecObservedToExpected[tTrigger]->getTH1()->Reset("ICE");
    m_xSecVsInstLumi[tTrigger]->getTH1()->Reset("ICE");
  }

  for (map<int, map<TString, double> >::iterator i = m_lsRates.begin(); i != m_lsRates.end(); i++) {
    unsigned int ls = (*i).first;
    rates = &(*i).second;
    isDefRate = true;

    if (m_lsLuminosity.find(ls) == m_lsLuminosity.end()) {
      isDefLumi = false;
    } else {
      isDefLumi = true;
      lumi = m_lsLuminosity[ls];
    }

    if (m_lsPrescaleIndex.find(ls) == m_lsPrescaleIndex.end()) {
      isDefPrescaleIndex = false;
    } else {
      isDefPrescaleIndex = true;
      prescalesIndex = m_lsPrescaleIndex[ls];
    }

    if (isDefRate && isDefLumi && isDefPrescaleIndex) {
      const vector<int>& currentPrescaleFactors = (*m_listsPrescaleFactors).at(prescalesIndex);

      for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
        string tTrigger = (*i).second;
        TF1* tTestFunction = (TF1*)m_xSecVsInstLumi[tTrigger]->getTProfile()->GetListOfFunctions()->First();

        // If trigger name is defined we get the rate fit parameters
        if (tTrigger != "Undefined") {
          unsigned int trigBit = m_algoBit[tTrigger];
          double trigPrescale = currentPrescaleFactors[trigBit];
          double trigRate = (*rates)[tTrigger];

          if (lumi != 0 && trigPrescale != 0 && trigRate != 0) {
            double AlgoXSec = (trigPrescale * trigRate) / lumi;
            double TemplateFunctionValue = tTestFunction->Eval(lumi);

            // Checking against Template function
            int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
            m_xSecObservedToExpected[tTrigger]->setBinContent(ibin, AlgoXSec / TemplateFunctionValue);
            m_xSecVsInstLumi[tTrigger]->Fill(lumi, AlgoXSec);

            if (m_verbose) {
              cout << "[L1TRate:] ls=" << ls << " Algo=" << tTrigger << " XSec=" << AlgoXSec
                   << " Test=" << AlgoXSec / TemplateFunctionValue << endl;
            }

          } else {
            int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
            m_xSecObservedToExpected[tTrigger]->setBinContent(ibin, 0.000001);
            if (m_verbose) {
              cout << "[L1TRate:] Algo=" << tTrigger << " XSec=Failed" << endl;
            }
          }
        }
      }
    }
  }
}

//_____________________________________________________________________
void L1TRate::analyze(const Event& iEvent, const EventSetup& eventSetup) {
  edm::Handle<L1GlobalTriggerReadoutRecord> gtReadoutRecordData;
  edm::Handle<Level1TriggerScalersCollection> triggerScalers;
  edm::Handle<LumiScalersCollection> colLScal;

  iEvent.getByToken(m_l1GtDataDaqInputTag, gtReadoutRecordData);
  iEvent.getByToken(m_scalersSource_colLScal, colLScal);
  iEvent.getByToken(m_scalersSource_triggerScalers, triggerScalers);

  // Integers
  int EventRun = iEvent.id().run();
  unsigned int eventLS = iEvent.id().luminosityBlock();

  // Getting the trigger trigger rates from GT and buffering it
  if (triggerScalers.isValid()) {
    Level1TriggerScalersCollection::const_iterator itL1TScalers = triggerScalers->begin();
    Level1TriggerRates trigRates(*itL1TScalers, EventRun);

    int gtLS = (*itL1TScalers).lumiSegmentNr() + m_lsShiftGTRates;

    // If we haven't got the data from this LS yet get it
    if (m_lsRates.find(gtLS) == m_lsRates.end()) {
      if (m_verbose) {
        cout << "[L1TRate:] Buffering GT Rates for LS=" << gtLS << endl;
      }
      map<TString, double> bufferRate;

      // Buffer the rate informations for all selected bits
      for (map<string, string>::const_iterator i = m_selectedTriggers.begin(); i != m_selectedTriggers.end(); i++) {
        string tTrigger = (*i).second;

        // If trigger name is defined we store the rate
        if (tTrigger != "Undefined") {
          unsigned int trigBit = m_algoBit[tTrigger];
          double trigRate = trigRates.gtAlgoCountsRate()[trigBit];

          bufferRate[tTrigger] = trigRate;
        }
      }
      m_lsRates[gtLS] = bufferRate;
    }
  }

  // Getting from the SCAL the luminosity information and buffering it
  if (colLScal.isValid() && !colLScal->empty()) {
    LumiScalersCollection::const_iterator itLScal = colLScal->begin();
    unsigned int scalLS = itLScal->sectionNumber();

    // If we haven't got the data from this SCAL LS yet get it
    if (m_lsLuminosity.find(scalLS) == m_lsLuminosity.end()) {
      if (m_verbose) {
        cout << "[L1TRate:] Buffering SCAL-HF Lumi for LS=" << scalLS << endl;
      }
      double instLumi = itLScal->instantLumi();                  // Getting Instant Lumi from HF (via SCAL)
      double deadTimeNormHF = itLScal->deadTimeNormalization();  // Getting Dead Time Normalization from HF (via SCAL)

      // If HF Dead Time Corrections is requested we apply it
      // NOTE: By default this is assumed false since for now WbM fits do NOT assume this correction
      if (m_parameters.getUntrackedParameter<bool>("useHFDeadTimeNormalization", false)) {
        // Protecting for deadtime = 0
        if (deadTimeNormHF == 0) {
          instLumi = 0;
        } else {
          instLumi = instLumi / deadTimeNormHF;
        }
      }
      // Buffering the luminosity information
      m_lsLuminosity[scalLS] = instLumi;
    }
  }

  // Getting the prescale index used when this event was triggered
  if (gtReadoutRecordData.isValid()) {
    // If we haven't got the data from this LS yet get it
    if (m_lsPrescaleIndex.find(eventLS) == m_lsPrescaleIndex.end()) {
      if (m_verbose) {
        cout << "[L1TRate:] Buffering Prescale Index for LS=" << eventLS << endl;
      }

      // Getting Final Decision Logic (FDL) Data from GT
      const vector<L1GtFdlWord>& gtFdlVectorData = gtReadoutRecordData->gtFdlVector();

      // Getting the index for the fdl data for this event
      int indexFDL = 0;
      for (unsigned int i = 0; i < gtFdlVectorData.size(); i++) {
        if (gtFdlVectorData[i].bxInEvent() == 0) {
          indexFDL = i;
          break;
        }
      }

      int CurrentPrescalesIndex = gtFdlVectorData[indexFDL].gtPrescaleFactorIndexAlgo();
      m_lsPrescaleIndex[eventLS] = CurrentPrescalesIndex;
    }
  }
}

//_____________________________________________________________________
// function: getXSexFitsOMDS
// Imputs:
//   * const edm::ParameterSet& ps = ParameterSet contaning necessary
//     information for the OMDS data extraction
// Outputs:
//   * int error = Number of algos where you did not find a
//     corresponding fit
//_____________________________________________________________________
bool L1TRate::getXSexFitsOMDS(const edm::ParameterSet& ps) {
  bool noError = true;

  string oracleDB = ps.getParameter<string>("oracleDB");
  string pathCondDB = ps.getParameter<string>("pathCondDB");

  L1TOMDSHelper myOMDSHelper;
  int conError;
  myOMDSHelper.connect(oracleDB, pathCondDB, conError);

  map<string, WbMTriggerXSecFit> wbmFits;

  if (conError == L1TOMDSHelper::NO_ERROR) {
    int errorRetrive;
    wbmFits = myOMDSHelper.getWbMAlgoXsecFits(errorRetrive);

    // Filling errors if they exist
    if (errorRetrive != L1TOMDSHelper::NO_ERROR) {
      noError = false;
      string eName = myOMDSHelper.enumToStringError(errorRetrive);
      m_ErrorMonitor->Fill(eName);
    }
  } else {
    noError = false;
    string eName = myOMDSHelper.enumToStringError(conError);
    m_ErrorMonitor->Fill(eName);
  }

  double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
  double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");

  // Getting rate fit parameters for all input triggers
  for (map<string, string>::const_iterator a = m_selectedTriggers.begin(); a != m_selectedTriggers.end(); a++) {
    string tTrigger = (*a).second;

    // If trigger name is defined we get the rate fit parameters
    if (tTrigger != "Undefined") {
      if (wbmFits.find(tTrigger) != wbmFits.end()) {
        WbMTriggerXSecFit tWbMParameters = wbmFits[tTrigger];

        vector<double> tParameters;
        tParameters.push_back(tWbMParameters.pm1);
        tParameters.push_back(tWbMParameters.p0);
        tParameters.push_back(tWbMParameters.p1);
        tParameters.push_back(tWbMParameters.p2);

        // Retriving and populating the m_templateFunctions array
        m_templateFunctions[tTrigger] = new TF1("FitParametrization_" + tWbMParameters.bitName,
                                                tWbMParameters.fitFunction,
                                                minInstantLuminosity,
                                                maxInstantLuminosity);
        m_templateFunctions[tTrigger]->SetParameters(&tParameters[0]);
        m_templateFunctions[tTrigger]->SetLineWidth(1);
        m_templateFunctions[tTrigger]->SetLineColor(kRed);

      } else {
        noError = false;
      }
    }
  }

  return noError;
}

//_____________________________________________________________________
// function: getXSexFitsPython
// Imputs:
//   * const edm::ParameterSet& ps = ParameterSet contaning the fit
//     functions and parameters for the selected triggers
// Outputs:
//   * int error = Number of algos where you did not find a
//     corresponding fit
//_____________________________________________________________________
bool L1TRate::getXSexFitsPython(const edm::ParameterSet& ps) {
  // error meaning
  bool noError = true;

  // Getting fit parameters
  std::vector<edm::ParameterSet> m_fitParameters = ps.getParameter<vector<ParameterSet> >("fitParameters");

  double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
  double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");

  // Getting rate fit parameters for all input triggers
  for (map<string, string>::const_iterator a = m_selectedTriggers.begin(); a != m_selectedTriggers.end(); a++) {
    string tTrigger = (*a).second;

    // If trigger name is defined we get the rate fit parameters
    if (tTrigger != "Undefined") {
      bool foundFit = false;

      for (unsigned int b = 0; b < m_fitParameters.size(); b++) {
        if (tTrigger == m_fitParameters[b].getParameter<string>("AlgoName")) {
          TString tAlgoName = m_fitParameters[b].getParameter<string>("AlgoName");
          TString tTemplateFunction = m_fitParameters[b].getParameter<string>("TemplateFunction");
          vector<double> tParameters = m_fitParameters[b].getParameter<vector<double> >("Parameters");

          // Retriving and populating the m_templateFunctions array
          m_templateFunctions[tTrigger] =
              new TF1("FitParametrization_" + tAlgoName, tTemplateFunction, minInstantLuminosity, maxInstantLuminosity);
          m_templateFunctions[tTrigger]->SetParameters(&tParameters[0]);
          m_templateFunctions[tTrigger]->SetLineWidth(1);
          m_templateFunctions[tTrigger]->SetLineColor(kRed);

          foundFit = true;
          break;
        }
      }
      if (!foundFit) {
        noError = false;
        string eName = "WARNING_PY_MISSING_FIT";
        m_ErrorMonitor->Fill(eName);
      }
    }
  }

  return noError;
}