CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoLuminosity/LumiProducer/plugins/LumiCalculator.cc

Go to the documentation of this file.
00001 #ifndef RecoLuminosity_LumiProducer_LumiCalculator_h
00002 #define RecoLuminosity_LumiProducer_LumiCalculator_h
00003 #include "FWCore/Framework/interface/EDAnalyzer.h"
00004 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00005 #include "FWCore/Framework/interface/Event.h"
00006 #include "FWCore/Framework/interface/Run.h"
00007 #include "FWCore/Framework/interface/LuminosityBlock.h"
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "DataFormats/Luminosity/interface/LumiSummary.h"
00010 #include "FWCore/Framework/interface/MakerMacros.h"
00011 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include <boost/regex.hpp>
00014 #include <iostream>
00015 #include <map>
00016 struct hltPerPathInfo{
00017   hltPerPathInfo():inputcount(0),outputcount(0),prescale(0){}
00018   unsigned int inputcount;
00019   unsigned int outputcount;
00020   unsigned int prescale;
00021 };
00022 struct l1PerBitInfo{
00023   l1PerBitInfo():count(0),prescale(0){}
00024   unsigned int count;
00025   unsigned int prescale;
00026 };
00027 struct MyPerLumiInfo{
00028   unsigned int lsnum;
00029   float livefraction;
00030   float intglumi;
00031   unsigned long long deadcount;
00032 };
00033 
00034 class LumiCalculator : public edm::EDAnalyzer{
00035 public:
00036   
00037   explicit LumiCalculator(edm::ParameterSet const& pset);
00038   virtual ~LumiCalculator();
00039 
00040 private:  
00041   virtual void beginJob();
00042   virtual void beginRun(const edm::Run& run, const edm::EventSetup& c);
00043   virtual void analyze(edm::Event const& e, edm::EventSetup const& c);
00044   virtual void endLuminosityBlock(edm::LuminosityBlock const& lumiBlock, 
00045                                   edm::EventSetup const& c);
00046   virtual void endRun(edm::Run const&, edm::EventSetup const&);
00047   virtual void endJob();
00048   std::vector<std::string> splitpathstr(const std::string& strValue,const std::string separator);
00049   HLTConfigProvider hltConfig_;
00050   std::multimap<std::string,std::string> trgpathMmap_;//key:hltpath,value:l1bit
00051   std::map<std::string,hltPerPathInfo> hltmap_;
00052   std::map<std::string,l1PerBitInfo> l1map_;//
00053   std::vector<MyPerLumiInfo> perrunlumiinfo_;
00054 private:
00055   edm::LogInfo* log_;
00056   bool showTrgInfo_;
00057   unsigned int currentlumi_;
00058 };//end class
00059 
00060 // -----------------------------------------------------------------
00061 
00062 LumiCalculator::LumiCalculator(edm::ParameterSet const& pset):log_( new edm::LogInfo("LumiReport")),currentlumi_(0){
00063   showTrgInfo_=pset.getUntrackedParameter<bool>("showTriggerInfo",false);
00064 }
00065 
00066 // -----------------------------------------------------------------
00067 
00068 LumiCalculator::~LumiCalculator(){
00069   delete log_; log_=0; 
00070 }
00071 
00072 // -----------------------------------------------------------------
00073 
00074 void LumiCalculator::analyze(edm::Event const& e,edm::EventSetup const&){
00075   
00076 }
00077 
00078 // -----------------------------------------------------------------
00079 
00080 void LumiCalculator::beginJob(){
00081   
00082 }
00083 
00084 // -----------------------------------------------------------------
00085 
00086 void LumiCalculator::beginRun(const edm::Run& run, const edm::EventSetup& c){
00087   //std::cout<<"I'm in run number "<<run.run()<<std::endl;
00088   //if(!hltConfig_.init("HLT")){
00089   // throw cms::Exception("HLT process cannot be initialized");
00090   //}
00091   bool changed(true);
00092   const std::string processname("HLT");
00093   if(!hltConfig_.init(run,c,processname,changed)){
00094     throw cms::Exception("HLT process cannot be initialized");
00095   }
00096   perrunlumiinfo_.clear();
00097   trgpathMmap_.clear();
00098   hltmap_.clear();
00099   l1map_.clear();
00100   //hltConfig_.dump("processName");
00101   //hltConfig_.dump("TableName");
00102   //hltConfig_.dump("Triggers");
00103   //hltConfig_.dump("Modules");  
00104   if(showTrgInfo_){
00105     *log_<<"======Trigger Configuration Overview======\n";
00106     *log_<<"Run "<<run.run()<<" Trigger Table : "<<hltConfig_.tableName()<<"\n";
00107   }
00108   unsigned int totaltrg=hltConfig_.size();
00109   for (unsigned int t=0;t<totaltrg;++t){
00110     std::string hltname(hltConfig_.triggerName(t));
00111     std::vector<std::string> numpathmodules=hltConfig_.moduleLabels(hltname);
00112     if(showTrgInfo_){
00113       *log_<<t<<" HLT path\t"<<hltname<<"\n";
00114     }
00115     hltPerPathInfo hlt;
00116     hlt.prescale=1;
00117     hltmap_.insert(std::make_pair(hltname,hlt));
00118     std::vector<std::string>::iterator hltpathBeg=numpathmodules.begin();
00119     std::vector<std::string>::iterator hltpathEnd=numpathmodules.end();
00120     unsigned int mycounter=0;
00121     for(std::vector<std::string>::iterator numpathmodule = hltpathBeg;
00122         numpathmodule!=hltpathEnd; ++numpathmodule ) {
00123       if (hltConfig_.moduleType(*numpathmodule) != "HLTLevel1GTSeed"){
00124         continue;
00125       }
00126       ++mycounter;
00127      
00128       edm::ParameterSet l1GTPSet=hltConfig_.modulePSet(*numpathmodule);
00129       std::string l1pathname=l1GTPSet.getParameter<std::string>("L1SeedsLogicalExpression");
00130       //*log_<<"numpathmodule "<< *numpathmodule <<" : l1pathname "<<l1pathname<<"\n";
00131       if(mycounter>1){
00132         if(showTrgInfo_){
00133           *log_<<"\tskip and erase previous seeds : multiple L1SeedsLogicalExpressions per hlt path\n";
00134         }
00135         //erase all previously calculated seeds for this path
00136         trgpathMmap_.erase(hltname);
00137         continue;
00138       }
00139       if(l1pathname.find("(")!=std::string::npos){
00140         if(showTrgInfo_){
00141           *log_<<"  L1SeedsLogicalExpression(Complex)\t"<<l1pathname<<"\n";
00142           *log_<<"\tskip:contain complex logic\n";
00143         }
00144         continue;
00145       }else if(l1pathname.find("OR")!=std::string::npos){
00146         if(showTrgInfo_){
00147           *log_<<"  L1SeedsLogicalExpression(ORed)\t"<<l1pathname<<"\n";
00148         }
00149         std::vector<std::string> seeds=splitpathstr(l1pathname," OR ");
00150         if(seeds.size()>2){
00151           if(showTrgInfo_){
00152             *log_<<"\tskip:contain >1 OR\n";
00153           }
00154           continue;
00155         }else{
00156           for(std::vector<std::string>::iterator i=seeds.begin();i!=seeds.end();++i){
00157             if(i->size()!=0 && showTrgInfo_)  *log_<<"\t\tseed: "<<*i<<"\n";
00158             if(i==seeds.begin()){//for now we take the first one from OR
00159                 trgpathMmap_.insert(std::make_pair(hltname,*i));
00160             }
00161           }
00162         }
00163       }else if (l1pathname.find("AND")!=std::string::npos){
00164         if(showTrgInfo_){
00165           *log_<<"  L1SeedsLogicalExpression(ANDed)\t"<< l1pathname<< "\n";
00166         }
00167         std::vector<std::string> seeds=splitpathstr(l1pathname," AND ");
00168         if(seeds.size()>2){
00169           if(showTrgInfo_){
00170             *log_<<"\tskip:contain >1 AND\n";
00171           }
00172           continue;
00173         }else{
00174           for(std::vector<std::string>::iterator i=seeds.begin();
00175               i!=seeds.end();++i){
00176             if(i->size()!=0 && showTrgInfo_) *log_<<"\t\tseed: "<<*i<<"\n";
00177             if(i==seeds.begin()){//for now we take the first one 
00178               trgpathMmap_.insert(std::make_pair(hltname,*i));
00179             }
00180           }
00181         }
00182       }else{
00183         if(showTrgInfo_){
00184           *log_<<"  L1SeedsLogicalExpression(ONE)\t"<< l1pathname<<"\n";
00185         }
00186         if(splitpathstr(l1pathname," NOT ").size()>1){
00187           if(showTrgInfo_){
00188             *log_<<"\tskip:contain NOT\n";
00189           }
00190           continue;
00191         }
00192         trgpathMmap_.insert(std::make_pair(hltname,l1pathname));
00193       } 
00194     }
00195   }
00196   if(showTrgInfo_){
00197     *log_<<"================\n";
00198   }
00199 }
00200 
00201 // -----------------------------------------------------------------
00202 void LumiCalculator::endLuminosityBlock(edm::LuminosityBlock const& lumiBlock, 
00203                                         edm::EventSetup const& c){
00207   //std::cout<<"I'm in lumi block "<<lumiBlock.id()<<std::endl;
00208   
00209   edm::Handle<LumiSummary> lumiSummary;
00210   lumiBlock.getByLabel("lumiProducer", lumiSummary);
00211   
00212   MyPerLumiInfo l;
00213   l.lsnum=lumiBlock.id().luminosityBlock();
00214 
00215   //
00216   //collect lumi. 
00217   //
00218   l.deadcount=lumiSummary->deadcount();
00219   l.intglumi=lumiSummary->avgInsDelLumi()*93.244;
00220   l.livefraction=lumiSummary->liveFrac();
00221   
00222   *log_<<"====== Lumi Section "<<lumiBlock.id().luminosityBlock()<<" ======\n";
00223   *log_<<"\t Luminosity "<<l.intglumi<<"\n";
00224   *log_<<"\t Dead count "<<l.deadcount<<"\n";
00225   *log_<<"\t Deadtime corrected Luminosity "<<l.intglumi*l.livefraction<<"\n";
00226 
00227   //
00228   //print correlated hlt-l1 info, only if you ask
00229   //
00230   if(showTrgInfo_){
00231     std::map<std::string,hltPerPathInfo>::iterator hltit;
00232     std::map<std::string,hltPerPathInfo>::iterator hltitBeg=hltmap_.begin();
00233     std::map<std::string,hltPerPathInfo>::iterator hltitEnd=hltmap_.end();
00234     
00235     typedef std::pair< std::multimap<std::string,std::string>::iterator,std::multimap<std::string,std::string>::iterator > TRGMAPIT;
00236     unsigned int c=0;
00237     for(hltit=hltitBeg;hltit!=hltitEnd;++hltit){
00238       std::string hltname=hltit->first;
00239       *log_<<c<<" HLT path  "<<hltname<<" , prescale : "<<hltit->second.prescale<<" , in : "<<hltit->second.inputcount<<" , out : "<<hltit->second.outputcount<<"\n";
00240       TRGMAPIT ppp;
00241       ppp=trgpathMmap_.equal_range(hltname);
00242       if(ppp.first==ppp.second){
00243         *log_<<"    no L1\n";
00244       }
00245       for(std::multimap<std::string,std::string>::iterator mit=ppp.first; mit!=ppp.second; ++mit){
00246         std::string l1name=mit->second;
00247         *log_<<"    L1 name : "<<l1name;
00248         LumiSummary::L1 l1result=lumiSummary->l1info(l1name);
00249         *log_<<" , count : "<<l1result.ratecount<<" , prescale : "<<l1result.prescale<<"\n";
00250         *log_<<"\n";
00251       }
00252       ++c;
00253     }
00254   }
00255   //
00256   //accumulate hlt counts. Absent for now
00257   //
00263   //
00264   //accumulate l1 counts
00265   //
00266   size_t n=lumiSummary->nTriggerLine();
00267   for(size_t i=0;i<n;++i){
00268     std::string l1bitname=lumiSummary->l1info(i).triggername;
00269     l1PerBitInfo t;
00270     if(currentlumi_==0){
00271       t.count=lumiSummary->l1info(i).ratecount;
00272       t.prescale=lumiSummary->l1info(i).prescale;
00273       l1map_.insert(std::make_pair(l1bitname,t));
00274     }else{
00275       std::map<std::string,l1PerBitInfo>::iterator it=l1map_.find(l1bitname);
00276       if(it!=l1map_.end()){
00277         it->second.count += lumiSummary->l1info(i).ratecount;
00278       }
00279     }
00280   }
00281   
00282   perrunlumiinfo_.push_back(l);
00283 
00284   ++currentlumi_;
00285 }
00286  
00287 // -----------------------------------------------------------------
00288 void LumiCalculator::endRun(edm::Run const& run, edm::EventSetup const& c){
00305   //std::cout<<"valid trigger lines "<<trgpathMmap_.size()<<std::endl;
00306   //std::cout<<"total lumi lines "<<perrunlumiinfo_.size()<<std::endl;
00307   std::vector<MyPerLumiInfo>::const_iterator lumiIt;
00308   std::vector<MyPerLumiInfo>::const_iterator lumiItBeg=perrunlumiinfo_.begin();
00309   std::vector<MyPerLumiInfo>::const_iterator lumiItEnd=perrunlumiinfo_.end(); 
00310   float recorded=0.0;
00311   
00312   *log_<<"================ Run Summary "<<run.run()<<"================\n";
00313   for(lumiIt=lumiItBeg;lumiIt!=lumiItEnd;++lumiIt){//loop over LS
00314     recorded += lumiIt->intglumi*lumiIt->livefraction;  
00315   }
00316   *log_<<"  CMS Recorded Lumi (e+27cm^-2) : "<<recorded<<"\n";
00317   *log_<<"  Effective Lumi (e+27cm^-2) per trigger path: "<<"\n\n";
00318   std::multimap<std::string,std::string>::iterator it;
00319   std::multimap<std::string,std::string>::iterator itBeg=trgpathMmap_.begin();
00320   std::multimap<std::string,std::string>::iterator itEnd=trgpathMmap_.end();
00321   unsigned int cc=0;
00322   for(it=itBeg;it!=itEnd;++it){
00323     *log_<<"  "<<cc<<"  "<<it->first<<" - "<<it->second<<" : ";
00324     ++cc;
00325     std::map<std::string,hltPerPathInfo>::const_iterator hltIt=hltmap_.find(it->first);
00326     if( hltIt==hltmap_.end() ){
00327       std::cout<<"HLT path "<<it->first<<" not found"<<std::endl;
00328       *log_<<"\n";
00329       continue;
00330     }
00331     std::map<std::string,l1PerBitInfo>::const_iterator l1It=l1map_.find(it->second);
00332     if( l1It==l1map_.end() ){
00333       std::cout<<"L1 bit "<<it->second<<" not found"<<std::endl;
00334       *log_<<"\n";
00335       continue;
00336     }
00337     unsigned int hltprescale=hltIt->second.prescale;
00338     unsigned int l1prescale=l1It->second.prescale;
00339     if( hltprescale!=0 && l1prescale!=0 ){
00340       float effectiveLumi=recorded/(hltprescale*l1prescale);
00341       *log_<<effectiveLumi<<"\n";
00342     }else{
00343       *log_<<"0 prescale exception\n";
00344       continue;
00345     }
00346     *log_<<"\n"; 
00347   }
00348 }
00349 
00350 
00351 // -----------------------------------------------------------------
00352 void LumiCalculator::endJob(){
00353 }
00354 
00355 std::vector<std::string>
00356 LumiCalculator::splitpathstr(const std::string& strValue,const std::string separator){
00357   std::vector<std::string> vecstrResult;
00358   boost::regex re(separator);
00359   boost::sregex_token_iterator p(strValue.begin(),strValue.end(),re,-1);
00360   boost::sregex_token_iterator end;
00361   while(p!=end){
00362     vecstrResult.push_back(*p++);
00363   }
00364   return vecstrResult;
00365 }
00366 
00367 DEFINE_FWK_MODULE(LumiCalculator);
00368 #endif