CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DQMOffline/Trigger/src/HLTTauDQMPathPlotter.cc

Go to the documentation of this file.
00001 #include "DQMOffline/Trigger/interface/HLTTauDQMPathPlotter.h"
00002 
00003 HLTTauDQMPathPlotter::HLTTauDQMPathPlotter( const edm::ParameterSet& ps, bool ref, std::string dqmBaseFolder ) {
00004     //Initialize Plotter
00005     name_ = "HLTTauDQMPathPlotter";
00006     
00007     //Process PSet
00008     try {
00009         triggerEventObject_   = ps.getUntrackedParameter<edm::InputTag>("TriggerEventObject");
00010         triggerTag_           = ps.getUntrackedParameter<std::string>("DQMFolder");
00011         triggerTagAlias_      = ps.getUntrackedParameter<std::string>("Alias","");
00012         filters_              = ps.getUntrackedParameter<std::vector<edm::ParameterSet> >("Filters");
00013         reference_            = ps.getUntrackedParameter<edm::ParameterSet>("Reference");
00014         refNTriggeredTaus_    = reference_.getUntrackedParameter<unsigned int>("NTriggeredTaus");
00015         refNTriggeredLeptons_ = reference_.getUntrackedParameter<unsigned int>("NTriggeredLeptons");
00016         refTauPt_             = reference_.getUntrackedParameter<double>("refTauPt",20);
00017         refLeptonPt_          = reference_.getUntrackedParameter<double>("refLeptonPt",15);
00018         dqmBaseFolder_        = dqmBaseFolder;
00019         doRefAnalysis_        = ref;
00020         validity_             = true;
00021     } catch ( cms::Exception &e ) {
00022         edm::LogInfo("HLTTauDQMPathPlotter::HLTTauDQMPathPlotter") << e.what() << std::endl;
00023         validity_ = false;
00024         return;
00025     }
00026     
00027     for ( std::vector<edm::ParameterSet>::const_iterator iter = filters_.begin(); iter != filters_.end(); ++iter ) {
00028         HLTTauDQMPlotter::FilterObject tmp(*iter);
00029         if (tmp.isValid()) filterObjs_.push_back(tmp);
00030     }
00031     
00032     if (store_) {
00033         //Create the histograms
00034         store_->setCurrentFolder(triggerTag());
00035         store_->removeContents();
00036         
00037         accepted_events = store_->book1D("TriggerBits","Accepted Events per Path;;entries",filterObjs_.size(),0,filterObjs_.size());
00038         for ( size_t k = 0; k < filterObjs_.size(); ++k ) {
00039             accepted_events->setBinLabel(k+1,filterObjs_[k].getAlias(),1);
00040         }
00041         if (doRefAnalysis_) {
00042             accepted_events_matched = store_->book1D("MatchedTriggerBits","Accepted+Matched Events per Path;;entries",filterObjs_.size()+1,0,filterObjs_.size()+1);
00043             accepted_events_matched->setBinLabel(1,"RefEvents",1);
00044             for ( size_t k = 0; k < filterObjs_.size(); ++k ) {
00045                 accepted_events_matched->setBinLabel(k+2,filterObjs_[k].getAlias(),1);
00046             }
00047         }
00048     }
00049 }
00050 
00051 HLTTauDQMPathPlotter::~HLTTauDQMPathPlotter() {
00052 }
00053 
00054 //
00055 // member functions
00056 //
00057 
00058 void HLTTauDQMPathPlotter::analyze( const edm::Event& iEvent, const edm::EventSetup& iSetup, const std::map<int,LVColl>& refC ) {
00059     using namespace std;
00060     using namespace edm;
00061     using namespace reco;
00062     using namespace l1extra;
00063     using namespace trigger;
00064     
00065     bool isGoodReferenceEvent = false;
00066     LVColl refTaus, refLeptons;
00067     
00068     if (doRefAnalysis_) {
00069         unsigned int highPtTaus = 0;
00070         unsigned int highPtElectrons = 0;
00071         unsigned int highPtMuons = 0;
00072         
00073         bool tau_ok = true;
00074         bool leptons_ok = true;
00075         
00076         std::map<int,LVColl>::const_iterator iref;
00077         
00078         //Tau reference
00079         iref = refC.find(15);
00080         if ( iref != refC.end() ) {
00081             for ( LVColl::const_iterator lvi = iref->second.begin(); lvi != iref->second.end(); ++lvi ) {
00082                 if ( lvi->Et() > refTauPt_ ) {
00083                     highPtTaus++;
00084                 }
00085                 refTaus.push_back(*lvi);
00086             }
00087         }
00088         if ( highPtTaus < refNTriggeredTaus_ ) tau_ok = false;
00089         
00090         //Electron reference
00091         iref = refC.find(11);
00092         if ( iref != refC.end() ) {
00093             for ( LVColl::const_iterator lvi = iref->second.begin(); lvi != iref->second.end(); ++lvi ) {
00094                 if ( lvi->Et() > refLeptonPt_ ) {
00095                     highPtElectrons++;
00096                 }
00097                 refLeptons.push_back(*lvi);
00098             }
00099         }
00100         if ( filterObjs_.size() > 0 && filterObjs_.back().leptonId() == 11 && highPtElectrons < refNTriggeredLeptons_ ) leptons_ok = false;
00101         
00102         //Muon reference
00103         iref = refC.find(13);
00104         if ( iref != refC.end() ) {
00105             for ( LVColl::const_iterator lvi = iref->second.begin(); lvi != iref->second.end(); ++lvi ) {
00106                 if ( lvi->Et() > refLeptonPt_ ) {
00107                     highPtMuons++;
00108                 }
00109                 refLeptons.push_back(*lvi);
00110             }
00111         }
00112         if ( filterObjs_.size() > 0 && filterObjs_.back().leptonId() == 13 && highPtMuons < refNTriggeredLeptons_ ) leptons_ok = false;
00113         
00114         if ( tau_ok && leptons_ok ) {
00115             accepted_events_matched->Fill(0.5);
00116             isGoodReferenceEvent = true;
00117         }
00118     }
00119     
00120     Handle<TriggerEventWithRefs> trigEv;
00121     bool gotTEV = iEvent.getByLabel(triggerEventObject_,trigEv) && trigEv.isValid();
00122     
00123     if (gotTEV) {
00124         //Loop through the filters
00125         for ( size_t i = 0; i < filterObjs_.size(); ++i ) {             
00126             size_t ID = trigEv->filterIndex(filterObjs_[i].getFilterName());
00127             
00128             if ( ID != trigEv->size() ) {
00129                 LVColl leptons = getFilterCollection(ID,filterObjs_[i].getLeptonType(),*trigEv);
00130                 LVColl taus = getFilterCollection(ID,filterObjs_[i].getTauType(),*trigEv);
00131                 //Fired
00132                 if ( leptons.size() >= filterObjs_[i].getNTriggeredLeptons() && taus.size() >= filterObjs_[i].getNTriggeredTaus() ) {
00133                     accepted_events->Fill(i+0.5);
00134                     //Now do the matching only though if we have a good reference event
00135                     if ( doRefAnalysis_ && isGoodReferenceEvent ) {
00136                         size_t nT = 0;
00137                         for ( size_t j = 0; j < taus.size(); ++j ) {
00138                             if( match(taus[j],refTaus,filterObjs_[i].getMatchDeltaR()).first ) nT++;
00139                         }
00140                         
00141                         size_t nL = 0;
00142                         for ( size_t j = 0; j < leptons.size(); ++j ) {
00143                             if ( match(leptons[j],refLeptons,filterObjs_[i].getMatchDeltaR()).first ) nL++;
00144                         }
00145                         
00146                         if ( nT >= filterObjs_[i].getNTriggeredTaus() && nL >= filterObjs_[i].getNTriggeredLeptons() ) {
00147                             accepted_events_matched->Fill(i+1.5);
00148                         }
00149                     }
00150                 }
00151             }
00152         }
00153     }
00154 }
00155 
00156 LVColl HLTTauDQMPathPlotter::getFilterCollection( size_t filterID, int id, const trigger::TriggerEventWithRefs& trigEv ) {
00157     using namespace trigger;
00158     LVColl out;
00159     
00160     if ( id == trigger::TriggerL1IsoEG || id == trigger::TriggerL1NoIsoEG ) {
00161         VRl1em obj;
00162         trigEv.getObjects(filterID,id,obj);
00163         for (size_t i=0;i<obj.size();++i)
00164             if (obj.at(i).isAvailable())
00165                 out.push_back(obj[i]->p4());
00166     }
00167     
00168     if ( id == trigger::TriggerL1Mu ) {
00169         VRl1muon obj;
00170         trigEv.getObjects(filterID,id,obj);
00171         for (size_t i=0;i<obj.size();++i)
00172             if (obj.at(i).isAvailable())
00173                 out.push_back(obj[i]->p4());
00174     }
00175     
00176     if ( id == trigger::TriggerMuon ) {
00177         VRmuon obj;
00178         trigEv.getObjects(filterID,id,obj);
00179         for (size_t i=0;i<obj.size();++i)
00180             if (obj.at(i).isAvailable())
00181                 out.push_back(obj[i]->p4());
00182     }
00183     
00184     if ( id == trigger::TriggerElectron ) {
00185         VRelectron obj;
00186         trigEv.getObjects(filterID,id,obj);
00187         for (size_t i=0;i<obj.size();++i)
00188             if (obj.at(i).isAvailable())
00189                 out.push_back(obj[i]->p4());
00190     }
00191     
00192     if ( id == trigger::TriggerL1TauJet ) {
00193         VRl1jet obj;
00194         trigEv.getObjects(filterID,id,obj);
00195         for (size_t i=0;i<obj.size();++i)
00196             if (obj.at(i).isAvailable())
00197                 out.push_back(obj[i]->p4());
00198         trigEv.getObjects(filterID,trigger::TriggerL1CenJet,obj);
00199         for (size_t i=0;i<obj.size();++i)
00200             if (obj.at(i).isAvailable())
00201                 out.push_back(obj[i]->p4());
00202         
00203     }
00204     
00205     if ( id == trigger::TriggerTau ) {
00206         VRjet obj;
00207         trigEv.getObjects(filterID,id,obj);
00208         for (size_t i = 0; i < obj.size(); ++i) {
00209             if (obj.at(i).isAvailable()) {
00210                 out.push_back(obj[i]->p4());
00211             }
00212         }
00213         VRpfjet pfjetobj;
00214         trigEv.getObjects(filterID,id,pfjetobj);
00215         for (size_t i = 0; i < pfjetobj.size(); ++i) {
00216             if (pfjetobj.at(i).isAvailable()) {
00217                 out.push_back(pfjetobj[i]->p4());
00218             }
00219         }
00220         VRpftau pftauobj;
00221         trigEv.getObjects(filterID,id,pftauobj);
00222         for (size_t i = 0; i < pftauobj.size(); ++i) {
00223             if (pftauobj.at(i).isAvailable()) {
00224                 out.push_back(pftauobj[i]->p4());
00225             }
00226         }
00227     }
00228 
00229     return out;
00230 }