CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/DQM/L1TMonitor/src/L1TRate.cc

Go to the documentation of this file.
00001 /*
00002  * \file L1TRate.cc
00003  *
00004  * $Date: 2011/11/15 10:41:00 $
00005  * $Revision: 1.11 $
00006  * \author J. Pela, P. Musella
00007  *
00008  */
00009 
00010 // L1TMonitor includes
00011 #include "DQM/L1TMonitor/interface/L1TRate.h"
00012 #include "DQM/L1TMonitor/interface/L1TMenuHelper.h"
00013 #include "DQM/L1TMonitor/interface/L1TOMDSHelper.h"
00014 
00015 #include "DQMServices/Core/interface/DQMStore.h"
00016 
00017 #include "DataFormats/Scalers/interface/LumiScalers.h"
00018 #include "DataFormats/Scalers/interface/Level1TriggerRates.h"
00019 #include "DataFormats/Scalers/interface/Level1TriggerScalers.h"
00020 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00021 #include "DataFormats/Common/interface/ConditionsInEdm.h" // Parameters associated to Run, LS and Event
00022 #include "DataFormats/Luminosity/interface/LumiDetails.h" // Luminosity Information
00023 #include "DataFormats/Luminosity/interface/LumiSummary.h" // Luminosity Information
00024 
00025 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenu.h"
00026 #include "CondFormats/L1TObjects/interface/L1GtTriggerMenuFwd.h"
00027 #include "CondFormats/L1TObjects/interface/L1GtPrescaleFactors.h"
00028 #include "CondFormats/L1TObjects/interface/L1GtTriggerMask.h"            // L1Gt - Masks
00029 #include "CondFormats/DataRecord/interface/L1GtTriggerMaskAlgoTrigRcd.h" // L1Gt - Masks
00030 #include "CondFormats/DataRecord/interface/L1GtTriggerMenuRcd.h"
00031 #include "CondFormats/DataRecord/interface/L1GtPrescaleFactorsAlgoTrigRcd.h"
00032 
00033 #include "DataFormats/Histograms/interface/MEtoEDMFormat.h"
00034 
00035 #include "TList.h"
00036 
00037 using namespace edm;
00038 using namespace std;
00039 
00040 //_____________________________________________________________________
00041 L1TRate::L1TRate(const ParameterSet & ps){
00042 
00043   m_maxNbins   = 2500; // Maximum LS for each run (for binning purposes)
00044   m_parameters = ps;
00045 
00046   // Mapping parameter input variables
00047   m_scalersSource       = m_parameters.getParameter         <InputTag>("inputTagScalersResults");
00048   m_l1GtDataDaqInputTag = m_parameters.getParameter         <InputTag>("inputTagL1GtDataDaq");
00049   m_verbose             = m_parameters.getUntrackedParameter<bool>    ("verbose",false);
00050   m_refPrescaleSet      = m_parameters.getParameter         <int>     ("refPrescaleSet");  
00051   m_lsShiftGTRates      = m_parameters.getUntrackedParameter<int>     ("lsShiftGTRates",0);
00052   
00053   // Getting which categories to monitor
00054   ParameterSet Categories     = ps.getParameter<ParameterSet>("categories");  
00055   m_inputCategories["Mu"]     = Categories.getUntrackedParameter<bool>("Mu"); 
00056   m_inputCategories["EG"]     = Categories.getUntrackedParameter<bool>("EG"); 
00057   m_inputCategories["IsoEG"]  = Categories.getUntrackedParameter<bool>("IsoEG"); 
00058   m_inputCategories["Jet"]    = Categories.getUntrackedParameter<bool>("Jet"); 
00059   m_inputCategories["CenJet"] = Categories.getUntrackedParameter<bool>("CenJet"); 
00060   m_inputCategories["ForJet"] = Categories.getUntrackedParameter<bool>("ForJet"); 
00061   m_inputCategories["TauJet"] = Categories.getUntrackedParameter<bool>("TauJet"); 
00062   m_inputCategories["ETM"]    = Categories.getUntrackedParameter<bool>("ETM"); 
00063   m_inputCategories["ETT"]    = Categories.getUntrackedParameter<bool>("ETT"); 
00064   m_inputCategories["HTT"]    = Categories.getUntrackedParameter<bool>("HTT"); 
00065   m_inputCategories["HTM"]    = Categories.getUntrackedParameter<bool>("HTM"); 
00066 
00067   // Inicializing Variables
00068   dbe = NULL;
00069 
00070   if (ps.getUntrackedParameter < bool > ("dqmStore", false)) {
00071     dbe = Service < DQMStore > ().operator->();
00072     dbe->setVerbose(0);
00073   }
00074   
00075   // What to do if we want our output to be saved to a external file
00076   m_outputFile = ps.getUntrackedParameter < string > ("outputFile", "");
00077   
00078   if (m_outputFile.size() != 0) {
00079     cout << "L1T Monitoring histograms will be saved to " << m_outputFile.c_str() << endl;
00080   }
00081   
00082   bool disable = ps.getUntrackedParameter < bool > ("disableROOToutput", false);
00083   if (disable) {m_outputFile = "";}
00084   
00085   if (dbe != NULL) {dbe->setCurrentFolder("L1T/L1TRate");}
00086   
00087 }
00088 
00089 //_____________________________________________________________________
00090 L1TRate::~L1TRate(){}
00091 
00092 //_____________________________________________________________________
00093 void L1TRate::beginJob(void){
00094 
00095   if (m_verbose) {cout << "[L1TRate:] Called beginJob." << endl;}
00096 
00097   // get hold of back-end interface
00098   DQMStore *dbe = 0;
00099   dbe = Service < DQMStore > ().operator->();
00100 
00101   if (dbe) {
00102     dbe->setCurrentFolder("L1T/L1TRate");
00103     dbe->rmdir("L1T/L1TRate");
00104   }
00105  
00106 }
00107 
00108 //_____________________________________________________________________
00109 void L1TRate::endJob(void){
00110 
00111   if (m_verbose) {cout << "[L1TRate:] Called endJob." << endl;}
00112 
00113   if (m_outputFile.size() != 0 && dbe)
00114     dbe->save(m_outputFile);
00115 
00116   return;
00117 }
00118 
00119 //_____________________________________________________________________
00120 // BeginRun
00121 //_____________________________________________________________________
00122 void L1TRate::beginRun(const edm::Run& run, const edm::EventSetup& iSetup){
00123 
00124   if (m_verbose) {cout << "[L1TRate:] Called beginRun." << endl;}
00125 
00126   ESHandle<L1GtTriggerMenu>     menuRcd;
00127   ESHandle<L1GtPrescaleFactors> l1GtPfAlgo;
00128 
00129   iSetup.get<L1GtTriggerMenuRcd>()            .get(menuRcd);
00130   iSetup.get<L1GtPrescaleFactorsAlgoTrigRcd>().get(l1GtPfAlgo);
00131 
00132   const L1GtTriggerMenu*     menu         = menuRcd   .product();
00133   const L1GtPrescaleFactors* m_l1GtPfAlgo = l1GtPfAlgo.product();
00134 
00135   // Initializing DQM Monitor Elements
00136   dbe->setCurrentFolder("L1T/L1TRate");
00137   m_ErrorMonitor = dbe->book1D("ErrorMonitor", "ErrorMonitor",5,0,5);
00138   m_ErrorMonitor->setBinLabel(1,"WARNING_DB_CONN_FAILED");        // Errors from L1TOMDSHelper
00139   m_ErrorMonitor->setBinLabel(2,"WARNING_DB_QUERY_FAILED");       // Errors from L1TOMDSHelper
00140   m_ErrorMonitor->setBinLabel(3,"WARNING_DB_INCORRECT_NBUNCHES"); // Errors from L1TOMDSHelper
00141   m_ErrorMonitor->setBinLabel(4,"WARNING_PY_MISSING_FIT");
00142   m_ErrorMonitor->setBinLabel(5,"UNKNOWN");
00143 
00144   // Retriving the list of prescale sets
00145   m_listsPrescaleFactors = &(m_l1GtPfAlgo->gtPrescaleFactors());
00146  
00147   // Getting Lowest Prescale Single Object Triggers from the menu
00148   L1TMenuHelper myMenuHelper = L1TMenuHelper(iSetup);
00149   m_selectedTriggers = myMenuHelper.getLUSOTrigger(m_inputCategories,m_refPrescaleSet);
00150 
00151   //-> Getting template fits for the algLo cross sections
00152   int srcAlgoXSecFit = m_parameters.getParameter<int>("srcAlgoXSecFit");
00153   if     (srcAlgoXSecFit == 0){getXSexFitsOMDS  (m_parameters);}
00154   else if(srcAlgoXSecFit == 1){getXSexFitsPython(m_parameters);}
00155 
00156   for (CItAlgo algo = menu->gtAlgorithmMap().begin(); algo!=menu->gtAlgorithmMap().end(); ++algo){
00157     m_algoBit[(algo->second).algoAlias()] = (algo->second).algoBitNumber();    
00158   }
00159 
00160   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
00161   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
00162 
00163   // Initializing DQM Monitor Elements
00164   for(map<string,string>::const_iterator i=m_selectedTriggers.begin() ; i!=m_selectedTriggers.end() ; i++){
00165 
00166     TString tCategory     = (*i).first;
00167     TString tTrigger      = (*i).second;
00168 
00169     TString tErrorMessage = "";  
00170     TF1*    tTestFunction;
00171 
00172     if(tTrigger != "Undefined" && m_templateFunctions.find(tTrigger) != m_templateFunctions.end()){
00173       tTestFunction = m_templateFunctions[tTrigger];
00174     }
00175     else if(tTrigger == "Undefined"){
00176       TString tFunc = "-1";
00177       tTestFunction = new TF1("FitParametrization_"+tTrigger,tFunc,0,double(m_maxNbins)-0.5);
00178     }
00179     else if(m_templateFunctions.find(tTrigger) == m_templateFunctions.end()){
00180       TString tFunc = "-1";
00181       tTestFunction = new TF1("FitParametrization_"+tTrigger,tFunc,0,double(m_maxNbins)-0.5);
00182       tErrorMessage = " (Undefined Test Function)";
00183     }
00184     else{
00185       TString tFunc = "-1";
00186       tTestFunction = new TF1("FitParametrization_"+tTrigger,tFunc,0,double(m_maxNbins)-0.5);
00187     }
00188 
00189     if(tTrigger != "Undefined"){
00190 
00191     if(myMenuHelper.getPrescaleByAlias(tCategory,tTrigger) != 1){
00192       tErrorMessage += " WARNING: Default Prescale = ";
00193       tErrorMessage += myMenuHelper.getPrescaleByAlias(tCategory,tTrigger);
00194     }
00195 
00196       if     (tCategory == "Mu"    && myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger) != 4294967295){
00197         tErrorMessage += " WARNING: Eta Range = ";
00198         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger);
00199       }
00200       else if(tCategory == "EG"    && myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger) != 32639){
00201         tErrorMessage += " WARNING: Eta Range = ";
00202         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger);
00203       }
00204       else if(tCategory == "IsoEG" && myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger) != 32639){
00205         tErrorMessage += " WARNING: Eta Range = ";
00206         tErrorMessage += myMenuHelper.getEtaRangeByAlias(tCategory,tTrigger);
00207       }
00208 
00209       if(tCategory == "Mu" && myMenuHelper.getQualityAlias(tCategory,tTrigger) != 240){
00210         tErrorMessage += " WARNING: Quality = ";
00211         tErrorMessage += myMenuHelper.getQualityAlias(tCategory,tTrigger);      
00212       }
00213 
00214     }
00215 
00216 
00217 
00218     dbe->setCurrentFolder("L1T/L1TRate/TriggerCrossSections");
00219     m_xSecVsInstLumi[tTrigger] = dbe->bookProfile(tCategory,
00220                                                   "Cross Sec. vs Inst. Lumi Algo: "+tTrigger+tErrorMessage,
00221                                                   m_maxNbins,
00222                                                   minInstantLuminosity,
00223                                                   maxInstantLuminosity,0,500); 
00224     m_xSecVsInstLumi[tTrigger] ->setAxisTitle("Instantaneous Luminosity [10^{30}cm^{-2}s^{-1}]" ,1);
00225     m_xSecVsInstLumi[tTrigger] ->setAxisTitle("Algorithm #sigma [#mu b]" ,2);
00226     m_xSecVsInstLumi[tTrigger] ->getTProfile()->GetListOfFunctions()->Add(tTestFunction);
00227     m_xSecVsInstLumi[tTrigger] ->getTProfile()->SetMarkerStyle(23);
00228 
00229     dbe->setCurrentFolder("L1T/L1TRate/Certification");
00230     m_xSecObservedToExpected[tTrigger] = dbe->book1D(tCategory, "Algo: "+tTrigger+tErrorMessage,m_maxNbins,-0.5,double(m_maxNbins)-0.5);
00231     m_xSecObservedToExpected[tTrigger] ->setAxisTitle("Lumi Section" ,1);
00232     m_xSecObservedToExpected[tTrigger] ->setAxisTitle("#sigma_{obs} / #sigma_{exp}" ,2);
00233 
00234   }  
00235 
00236 }
00237 
00238 //_____________________________________________________________________
00239 void L1TRate::endRun(const edm::Run& run, const edm::EventSetup& iSetup){
00240   if (m_verbose) {cout << "[L1TRate:] Called endRun." << endl;}
00241 }
00242 
00243 //_____________________________________________________________________
00244 void L1TRate::beginLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
00245 
00246   if (m_verbose) {cout << "[L1TRate:] Called beginLuminosityBlock at LS=" << lumiBlock.id().luminosityBlock() << endl;}
00247 
00248 }
00249 
00250 //_____________________________________________________________________
00251 void L1TRate::endLuminosityBlock(LuminosityBlock const& lumiBlock, EventSetup const& c) {
00252 
00253   int eventLS = lumiBlock.id().luminosityBlock();  
00254   if (m_verbose) {cout << "[L1TRate:] Called endLuminosityBlock at LS=" << eventLS << endl;}
00255 
00256   // We can certify LS -1 since we should have available:
00257   // gt rates: (current LS)-1
00258   // prescale: current LS
00259   // lumi    : current LS
00260   //eventLS--;
00261   
00262   // Checking if all necessary quantities are defined for our calculations
00263   bool isDefRate,isDefLumi,isDefPrescaleIndex;
00264   map<TString,double>* rates=0;
00265   double               lumi=0;
00266   int                  prescalesIndex=0;
00267 
00268   // Reseting MonitorElements so we can refill them
00269   for(map<string,string>::const_iterator i=m_selectedTriggers.begin() ; i!=m_selectedTriggers.end() ; i++){
00270     string tTrigger      = (*i).second;
00271     m_xSecObservedToExpected[tTrigger]->getTH1()->Reset("ICE");
00272     m_xSecVsInstLumi        [tTrigger]->getTH1()->Reset("ICE");
00273   }
00274     
00275   for(map<int,map<TString,double> >::iterator i=m_lsRates.begin() ; i!=m_lsRates.end() ; i++){
00276 
00277     unsigned int ls =  (*i).first;
00278     rates   = &(*i).second;
00279     isDefRate=true;
00280 
00281     if(m_lsLuminosity.find(ls)==m_lsLuminosity.end()){isDefLumi=false;}
00282     else{
00283       isDefLumi=true;
00284       lumi=m_lsLuminosity[ls];
00285     }
00286   
00287     if(m_lsPrescaleIndex.find(ls)==m_lsPrescaleIndex.end()){isDefPrescaleIndex=false;}
00288     else{
00289       isDefPrescaleIndex=true;
00290       prescalesIndex=m_lsPrescaleIndex[ls];
00291     }
00292     
00293     if(isDefRate && isDefLumi && isDefPrescaleIndex){
00294     
00295       const vector<int>& currentPrescaleFactors = (*m_listsPrescaleFactors).at(prescalesIndex);
00296      
00297       for(map<string,string>::const_iterator i=m_selectedTriggers.begin() ; i!=m_selectedTriggers.end() ; i++){
00298 
00299         string tTrigger      = (*i).second;
00300         TF1*   tTestFunction = (TF1*) m_xSecVsInstLumi[tTrigger]->getTProfile()->GetListOfFunctions()->First();
00301 
00302         // If trigger name is defined we get the rate fit parameters 
00303         if(tTrigger != "Undefined"){
00304 
00305           unsigned int   trigBit      = m_algoBit[tTrigger];
00306           double trigPrescale = currentPrescaleFactors[trigBit];
00307           double trigRate     = (*rates)[tTrigger];
00308 
00309           if(lumi!=0 && trigPrescale!=0 && trigRate!=0){
00310 
00311             double AlgoXSec              = (trigPrescale*trigRate)/lumi;
00312             double TemplateFunctionValue = tTestFunction->Eval(lumi);
00313 
00314             // Checking against Template function
00315             int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
00316             m_xSecObservedToExpected[tTrigger]->setBinContent(ibin,AlgoXSec/TemplateFunctionValue);
00317             m_xSecVsInstLumi        [tTrigger]->Fill(lumi,AlgoXSec);
00318   
00319             if(m_verbose){cout<<"[L1TRate:] ls="<<ls<<" Algo="<<tTrigger<<" XSec="<<AlgoXSec<<" Test="<<AlgoXSec/TemplateFunctionValue<<endl;}
00320 
00321           }
00322           else{
00323             int ibin = m_xSecObservedToExpected[tTrigger]->getTH1()->FindBin(ls);
00324             m_xSecObservedToExpected[tTrigger]->setBinContent(ibin,0.000001);
00325             if(m_verbose){cout << "[L1TRate:] Algo="<< tTrigger<< " XSec=Failed" << endl;}
00326           }
00327         }
00328       }  
00329     }    
00330   }
00331 }
00332 
00333 //_____________________________________________________________________
00334 void L1TRate::analyze(const Event & iEvent, const EventSetup & eventSetup){
00335 
00336   edm::Handle<L1GlobalTriggerReadoutRecord>   gtReadoutRecordData;
00337   edm::Handle<Level1TriggerScalersCollection> triggerScalers;
00338   edm::Handle<LumiScalersCollection>          colLScal;
00339  
00340   iEvent.getByLabel(m_l1GtDataDaqInputTag, gtReadoutRecordData);
00341   iEvent.getByLabel(m_scalersSource      , colLScal);
00342   iEvent.getByLabel(m_scalersSource      , triggerScalers);
00343 
00344   // Integers
00345   int  EventRun = iEvent.id().run();
00346   unsigned int eventLS  = iEvent.id().luminosityBlock();
00347 
00348   // Getting the trigger trigger rates from GT and buffering it
00349   if(triggerScalers.isValid()){
00350     
00351     Level1TriggerScalersCollection::const_iterator itL1TScalers = triggerScalers->begin();
00352     Level1TriggerRates trigRates(*itL1TScalers,EventRun);
00353     
00354     int gtLS = (*itL1TScalers).lumiSegmentNr()+m_lsShiftGTRates;
00355     
00356     // If we haven't got the data from this LS yet get it
00357     if(m_lsRates.find(gtLS)==m_lsRates.end()){
00358     
00359       if (m_verbose) {cout << "[L1TRate:] Buffering GT Rates for LS=" << gtLS << endl;}
00360       map<TString,double> bufferRate;
00361       
00362       // Buffer the rate informations for all selected bits
00363       for(map<string,string>::const_iterator i=m_selectedTriggers.begin() ; i!=m_selectedTriggers.end() ; i++){
00364 
00365         string tTrigger = (*i).second;
00366 
00367         // If trigger name is defined we store the rate
00368         if(tTrigger != "Undefined"){
00369 
00370           unsigned int   trigBit  = m_algoBit[tTrigger];
00371           double trigRate = trigRates.gtAlgoCountsRate()[trigBit]; 
00372   
00373           bufferRate[tTrigger] = trigRate;
00374         }
00375       }
00376       m_lsRates[gtLS] = bufferRate;
00377     }
00378   }
00379   
00380   // Getting from the SCAL the luminosity information and buffering it
00381   if(colLScal.isValid() && colLScal->size()){
00382     
00383     LumiScalersCollection::const_iterator itLScal = colLScal->begin();
00384     unsigned int scalLS  = itLScal->sectionNumber();
00385     
00386     // If we haven't got the data from this SCAL LS yet get it
00387     if(m_lsLuminosity.find(scalLS)==m_lsLuminosity.end()){
00388     
00389       if (m_verbose) {cout << "[L1TRate:] Buffering SCAL-HF Lumi for LS=" << scalLS << endl;}
00390       double instLumi       = itLScal->instantLumi();           // Getting Instant Lumi from HF (via SCAL)   
00391       double deadTimeNormHF = itLScal->deadTimeNormalization(); // Getting Dead Time Normalization from HF (via SCAL)
00392        
00393       // If HF Dead Time Corrections is requested we apply it
00394       // NOTE: By default this is assumed false since for now WbM fits do NOT assume this correction
00395       if(m_parameters.getUntrackedParameter<bool>("useHFDeadTimeNormalization",false)){
00396 
00397         // Protecting for deadtime = 0
00398         if(deadTimeNormHF==0){instLumi = 0;}
00399         else                 {instLumi = instLumi/deadTimeNormHF;}
00400       }
00401       // Buffering the luminosity information
00402       m_lsLuminosity[scalLS]=instLumi;
00403     }
00404   }
00405 
00406   // Getting the prescale index used when this event was triggered
00407   if(gtReadoutRecordData.isValid()){
00408     
00409     // If we haven't got the data from this LS yet get it
00410     if(m_lsPrescaleIndex.find(eventLS)==m_lsPrescaleIndex.end()){
00411       
00412       if (m_verbose) {cout << "[L1TRate:] Buffering Prescale Index for LS=" << eventLS << endl;}
00413 
00414       // Getting Final Decision Logic (FDL) Data from GT
00415       const vector<L1GtFdlWord>& gtFdlVectorData = gtReadoutRecordData->gtFdlVector();
00416 
00417       // Getting the index for the fdl data for this event
00418       int indexFDL=0;
00419       for(unsigned int i=0; i<gtFdlVectorData.size(); i++){
00420         if(gtFdlVectorData[i].bxInEvent()==0){indexFDL=i; break;}
00421       }
00422       
00423       int CurrentPrescalesIndex  = gtFdlVectorData[indexFDL].gtPrescaleFactorIndexAlgo();
00424       m_lsPrescaleIndex[eventLS] = CurrentPrescalesIndex;   
00425     }    
00426   }
00427 }
00428 
00429 //_____________________________________________________________________
00430 // function: getXSexFitsOMDS
00431 // Imputs: 
00432 //   * const edm::ParameterSet& ps = ParameterSet contaning necessary
00433 //     information for the OMDS data extraction
00434 // Outputs:
00435 //   * int error = Number of algos where you did not find a 
00436 //     corresponding fit 
00437 //_____________________________________________________________________
00438 bool L1TRate::getXSexFitsOMDS(const edm::ParameterSet& ps){
00439 
00440   bool noError = true;
00441 
00442   string oracleDB   = ps.getParameter<string>("oracleDB");
00443   string pathCondDB = ps.getParameter<string>("pathCondDB");
00444 
00445   L1TOMDSHelper myOMDSHelper;
00446   int conError;
00447   myOMDSHelper.connect(oracleDB,pathCondDB,conError);
00448  
00449   map<string,WbMTriggerXSecFit> wbmFits;
00450 
00451   if(conError == L1TOMDSHelper::NO_ERROR){
00452     int errorRetrive;
00453     wbmFits = myOMDSHelper.getWbMAlgoXsecFits(errorRetrive);
00454     
00455     // Filling errors if they exist
00456     if(errorRetrive != L1TOMDSHelper::NO_ERROR){
00457       noError = false;
00458       string eName = myOMDSHelper.enumToStringError(errorRetrive);
00459       m_ErrorMonitor->Fill(eName);
00460     }
00461   }else{
00462     noError = false;
00463     string eName = myOMDSHelper.enumToStringError(conError);
00464     m_ErrorMonitor->Fill(eName);
00465   }  
00466 
00467   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
00468   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
00469 
00470   // Getting rate fit parameters for all input triggers
00471   for(map<string,string>::const_iterator a=m_selectedTriggers.begin() ; a!=m_selectedTriggers.end() ; a++){
00472 
00473     string tTrigger = (*a).second;
00474 
00475     // If trigger name is defined we get the rate fit parameters 
00476     if(tTrigger != "Undefined"){
00477       
00478       if(wbmFits.find(tTrigger) != wbmFits.end()){
00479 
00480         WbMTriggerXSecFit tWbMParameters = wbmFits[tTrigger];
00481 
00482         vector<double> tParameters;
00483         tParameters.push_back(tWbMParameters.pm1);
00484         tParameters.push_back(tWbMParameters.p0);
00485         tParameters.push_back(tWbMParameters.p1);
00486         tParameters.push_back(tWbMParameters.p2);
00487         
00488         // Retriving and populating the m_templateFunctions array
00489         m_templateFunctions[tTrigger] = new TF1("FitParametrization_"+tWbMParameters.bitName,
00490                                                 tWbMParameters.fitFunction,
00491                                                 minInstantLuminosity,maxInstantLuminosity);
00492         m_templateFunctions[tTrigger] ->SetParameters(&tParameters[0]);
00493         m_templateFunctions[tTrigger] ->SetLineWidth(1);
00494         m_templateFunctions[tTrigger] ->SetLineColor(kRed);
00495 
00496       }else{noError = false;}
00497     }
00498   }
00499 
00500   return noError;
00501 
00502 }
00503 
00504 //_____________________________________________________________________
00505 // function: getXSexFitsPython
00506 // Imputs: 
00507 //   * const edm::ParameterSet& ps = ParameterSet contaning the fit
00508 //     functions and parameters for the selected triggers
00509 // Outputs:
00510 //   * int error = Number of algos where you did not find a 
00511 //     corresponding fit 
00512 //_____________________________________________________________________
00513 bool L1TRate::getXSexFitsPython(const edm::ParameterSet& ps){
00514 
00515   // error meaning
00516   bool noError = true;
00517 
00518   // Getting fit parameters
00519   std::vector<edm::ParameterSet>  m_fitParameters = ps.getParameter< vector<ParameterSet> >("fitParameters");
00520 
00521   double minInstantLuminosity = m_parameters.getParameter<double>("minInstantLuminosity");
00522   double maxInstantLuminosity = m_parameters.getParameter<double>("maxInstantLuminosity");
00523   
00524   // Getting rate fit parameters for all input triggers
00525   for(map<string,string>::const_iterator a=m_selectedTriggers.begin() ; a!=m_selectedTriggers.end() ; a++){
00526 
00527     string tTrigger = (*a).second;
00528 
00529     // If trigger name is defined we get the rate fit parameters 
00530     if(tTrigger != "Undefined"){
00531 
00532       bool foundFit = false;
00533 
00534       for(unsigned int b=0 ; b<m_fitParameters.size() ; b++){
00535         
00536         if(tTrigger == m_fitParameters[b].getParameter<string>("AlgoName")){
00537           
00538           TString        tAlgoName          = m_fitParameters[b].getParameter< string >        ("AlgoName");
00539           TString        tTemplateFunction  = m_fitParameters[b].getParameter< string >        ("TemplateFunction");
00540           vector<double> tParameters        = m_fitParameters[b].getParameter< vector<double> >("Parameters");
00541           
00542           // Retriving and populating the m_templateFunctions array
00543           m_templateFunctions[tTrigger] = new TF1("FitParametrization_"+tAlgoName,tTemplateFunction,
00544                                                   minInstantLuminosity,maxInstantLuminosity);
00545           m_templateFunctions[tTrigger] ->SetParameters(&tParameters[0]);
00546           m_templateFunctions[tTrigger] ->SetLineWidth(1);
00547           m_templateFunctions[tTrigger] ->SetLineColor(kRed);
00548 
00549           foundFit = true;
00550           break;
00551         }
00552 
00553         if(!foundFit){
00554           noError = false;
00555           string eName = "WARNING_PY_MISSING_FIT";
00556           m_ErrorMonitor->Fill(eName);
00557         }
00558       }
00559     }
00560   }
00561   
00562   return noError;
00563 
00564 }
00565 
00566 //define this as a plug-in
00567 DEFINE_FWK_MODULE(L1TRate);