CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/DQMOffline/Trigger/src/HLTMuonOverlap.cc

Go to the documentation of this file.
00001 
00005 //
00006 #include "DQMOffline/Trigger/interface/HLTMuonOverlap.h"
00007 
00008 // Collaborating Class Header
00009 #include "FWCore/Framework/interface/MakerMacros.h"
00010 #include "FWCore/Framework/interface/Frameworkfwd.h"
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 #include "DataFormats/Common/interface/Handle.h"
00015 #include "DataFormats/Common/interface/RefToBase.h"
00016 
00017 #include "DataFormats/Common/interface/TriggerResults.h"
00018 #include "FWCore/Common/interface/TriggerNames.h"
00019 
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 
00022 
00023 HLTMuonOverlap::HLTMuonOverlap(const edm::ParameterSet& pset)
00024 {
00025   TrigResultsIn=true;
00026   Ntp = 0;
00027   Nall_trig = Nevents = 0;
00028   theCrossSection = pset.getParameter<double>("CrossSection");
00029   // Convert it already into /nb/s)
00030   theLuminosity = pset.getUntrackedParameter<double>("Luminosity",1.e+32)*1.e-33;
00031   TrigResLabel_=pset.getParameter<edm::InputTag>("TriggerResultLabel");
00032   init_=false;
00033   size=0;
00034 }
00035 
00036 
00037 
00038 void HLTMuonOverlap::begin() {
00039 }
00040 
00041 
00042 
00043 void HLTMuonOverlap::analyze(const edm::Event & event ) {
00044   using namespace edm;  
00045 
00046   if (!TrigResultsIn)return;
00047   ++Nevents;
00048   
00049   Handle<TriggerResults> trigRes;
00050   event.getByLabel(TrigResLabel_, trigRes);
00051   if (!trigRes.isValid()){
00052   edm::InputTag triggerResultsLabelFU(TrigResLabel_.label(),TrigResLabel_.instance(), "FU");
00053   event.getByLabel(triggerResultsLabelFU,trigRes);
00054   if(!trigRes.isValid()) {
00055     LogWarning("HLTMuonVal")<< "No trigger Results";
00056     TrigResultsIn=false;
00057     // Do nothing
00058     return;
00059    }
00060   }
00061   size=trigRes->size();
00062   LogTrace("HLTMuonVal")<< "Ntp="<<Ntp<<" Size of trigger results="<<size;
00063   const edm::TriggerNames & triggerNames = event.triggerNames(*trigRes);
00064 
00065   if(Ntp)
00066     assert(Ntp == size);
00067   else
00068     Ntp = size;
00069   
00070   // has any trigger fired this event?
00071   if(trigRes->accept())++Nall_trig;
00072   else return;
00073   LogTrace("HLTMuonVal")<<" This event has fired ";
00074   // loop over all paths, get trigger decision
00075   for(unsigned i = 0; i != size; ++i)
00076     {
00077       std::string name = triggerNames.triggerName(i);
00078       LogTrace("HLTMuonVal") << name << " has decision "<<trigRes->accept(i);
00079       fired[name] = trigRes->accept(i);
00080       if(fired[name]) ++(Ntrig[name]);
00081     }
00082   
00083   // NOTE: WE SHOULD MAKE THIS A SYMMETRIC MATRIX...
00084   // double-loop over all paths, get trigger overlaps
00085   for(unsigned i = 0; i != size; ++i)
00086     {
00087       std::string name = triggerNames.triggerName(i);
00088       if(!fired[name])continue;
00089       
00090       bool correlation = false;
00091       
00092       for(unsigned j = 0; j != size; ++j)
00093         {
00094           // skip the same name; 
00095           // this entry correponds to events triggered by single trigger
00096           if(i == j) continue;
00097           std::string name2 = triggerNames.triggerName(j);
00098           if(fired[name2])
00099             {
00100               correlation = true;
00101               ++(Ncross[name][name2]);
00102             }
00103         } // loop over j-trigger paths
00104       
00105       if(!correlation) // events triggered by single trigger
00106         ++(Ncross[name][name]);
00107       
00108     } //  // loop over i-trigger paths
00109   
00110 }
00111 
00112 
00113 void HLTMuonOverlap::finish()
00114 {
00115   using namespace edm;
00116   if (!TrigResultsIn || Nevents == 0 )return;
00117   LogVerbatim("HLTMuonVal") << " Total trigger fraction: " << Nall_trig << "/" << Nevents
00118                             << " events (" << 100*Nall_trig/Nevents<<"%), the Rate="<< Nall_trig*theLuminosity*theCrossSection/Nevents << "Hz ";
00119   
00120   LogVerbatim("HLTMuonVal") << " Individual path rates: " ;
00121   typedef trigPath::iterator It;
00122   int ix = 1;
00123   for(It it = Ntrig.begin(); 
00124       it != Ntrig.end(); ++it, ++ix)
00125     {
00126       LogVerbatim("HLTMuonVal") << " Trigger path \"" << it->first << "\": " 
00127                                 << it->second << "/"
00128                                 << Nevents << " events, Rate=" << (it->second)*theLuminosity*theCrossSection/Nevents 
00129                                 << "Hz " ;
00130     }
00131     
00132   LogVerbatim("HLTMuonVal") << " Trigger path correlations: " ;
00133   typedef std::map<std::string, trigPath>::iterator IIt;
00134     
00135   ix = 1;
00136   for(IIt it = Ncross.begin(); 
00137       it != Ncross.end(); ++it, ++ix)
00138     { // loop over first trigger of pair
00139         
00140       trigPath & cross = it->second;
00141       int iy = 1;
00142       for(It it2 = cross.begin(); 
00143           it2 != cross.end(); ++it2, ++iy)
00144         { // loop over second trigger of pair
00145           // skip symmetric pairs: 1st pass does "path1", "path2";
00146           // 2nd pass should skip "path2", "path1"
00147           if(it->first > it2->first)continue;
00148 
00149           // if second trigger = first trigger, 
00150           // this corresponds to unique rate (ie. no correlation)
00151           if(it->first == it2->first)
00152             LogVerbatim("HLTMuonVal") << " \"" << it->first << "\"" 
00153                                       << " (unique rate): "<< it2->second << "/"
00154                                       << Nevents << " events, Rate=" 
00155                                       << (it2->second)*theLuminosity*theCrossSection/Nevents 
00156                                       << "Hz ";
00157           else
00158             LogVerbatim("HLTMuonVal") << " \"" << it->first << "\""<< " x \"" 
00159                                       << it2->first << "\": "<< it2->second 
00160                                       << "/"<< Nevents << " events, Rate=" 
00161                                       << (it2->second)*theLuminosity*theCrossSection/Nevents 
00162                                       << "Hz ";
00163 
00164         }
00165     }
00166 
00167 }