CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/HLTrigger/special/src/HLTPixelIsolTrackFilter.cc

Go to the documentation of this file.
00001 #include "HLTrigger/special/interface/HLTPixelIsolTrackFilter.h"
00002 
00003 #include "DataFormats/HcalIsolatedTrack/interface/IsolatedPixelTrackCandidate.h"
00004 
00005 #include "DataFormats/Common/interface/Handle.h"
00006 
00007 #include "DataFormats/Common/interface/RefToBase.h"
00008 
00009 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00010 
00011 #include "DataFormats/HLTReco/interface/TriggerTypeDefs.h"
00012 
00013 //#include "DataFormats/HLTReco/interface/HLTFilterObject.h"
00014 
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00017 
00018 HLTPixelIsolTrackFilter::HLTPixelIsolTrackFilter(const edm::ParameterSet& iConfig)
00019 {
00020   candTag_             = iConfig.getParameter<edm::InputTag> ("candTag");
00021   hltGTseedlabel_      = iConfig.getParameter<edm::InputTag> ("L1GTSeedLabel");
00022   minDeltaPtL1Jet_     = iConfig.getParameter<double> ("MinDeltaPtL1Jet");
00023   minpttrack_          = iConfig.getParameter<double> ("MinPtTrack");
00024   maxptnearby_         = iConfig.getParameter<double> ("MaxPtNearby");
00025   maxetatrack_         = iConfig.getParameter<double> ("MaxEtaTrack");
00026   minetatrack_         = iConfig.getParameter<double> ("MinEtaTrack");
00027   filterE_             = iConfig.getParameter<bool> ("filterTrackEnergy");
00028   minEnergy_           = iConfig.getParameter<double>("MinEnergyTrack");
00029   nMaxTrackCandidates_ = iConfig.getParameter<int>("NMaxTrackCandidates");
00030   dropMultiL2Event_    = iConfig.getParameter<bool> ("DropMultiL2Event");
00031   //register your products
00032   produces<trigger::TriggerFilterObjectWithRefs>();
00033 }
00034 
00035 HLTPixelIsolTrackFilter::~HLTPixelIsolTrackFilter(){}
00036 
00037 bool HLTPixelIsolTrackFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00038 {
00039 
00040   // The Filter object
00041   std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00042 
00043   // Ref to Candidate object to be recorded in filter object
00044   edm::Ref<reco::IsolatedPixelTrackCandidateCollection> candref;
00045 
00046   // get hold of filtered candidates
00047   edm::Handle<reco::IsolatedPixelTrackCandidateCollection> recotrackcands;
00048   iEvent.getByLabel(candTag_,recotrackcands);
00049 
00050   //Filtering
00051 
00052   //find leading L1 jet:
00053   double ptTriggered  = -10;
00054   double etaTriggered = -100;
00055   double phiTriggered = -100;
00056 
00057   edm::Handle<trigger::TriggerFilterObjectWithRefs> l1trigobj;
00058   iEvent.getByLabel(hltGTseedlabel_, l1trigobj);
00059   
00060   std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1tauobjref;
00061   std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1jetobjref;
00062   std::vector< edm::Ref<l1extra::L1JetParticleCollection> > l1forjetobjref;
00063   
00064   l1trigobj->getObjects(trigger::TriggerL1TauJet, l1tauobjref);
00065   l1trigobj->getObjects(trigger::TriggerL1CenJet, l1jetobjref);
00066   l1trigobj->getObjects(trigger::TriggerL1ForJet, l1forjetobjref);
00067   
00068   for (unsigned int p=0; p<l1tauobjref.size(); p++)
00069     {
00070       if (l1tauobjref[p]->pt()>ptTriggered)
00071         {
00072           ptTriggered  = l1tauobjref[p]->pt(); 
00073           phiTriggered = l1tauobjref[p]->phi();
00074           etaTriggered = l1tauobjref[p]->eta();
00075         }
00076     }
00077   for (unsigned int p=0; p<l1jetobjref.size(); p++)
00078     {
00079       if (l1jetobjref[p]->pt()>ptTriggered)
00080         {
00081           ptTriggered  = l1jetobjref[p]->pt();
00082           phiTriggered = l1jetobjref[p]->phi();
00083           etaTriggered = l1jetobjref[p]->eta();
00084         }
00085     }
00086   for (unsigned int p=0; p<l1forjetobjref.size(); p++)
00087     {
00088       if (l1forjetobjref[p]->pt()>ptTriggered)
00089         {
00090           ptTriggered=l1forjetobjref[p]->pt();
00091           phiTriggered=l1forjetobjref[p]->phi();
00092           etaTriggered=l1forjetobjref[p]->eta();
00093         }
00094     }
00095 
00096   int n=0;
00097   for (unsigned int i=0; i<recotrackcands->size(); i++)
00098     {
00099       candref = edm::Ref<reco::IsolatedPixelTrackCandidateCollection>(recotrackcands, i);
00100       
00101       // cut on deltaPT
00102       if (ptTriggered-candref->pt()<minDeltaPtL1Jet_) continue;
00103       
00104       // select on transverse momentum
00105       if (!filterE_&&(candref->maxPtPxl()<maxptnearby_)&&
00106           (candref->pt()>minpttrack_)&&fabs(candref->track()->eta())<maxetatrack_&&fabs(candref->track()->eta())>minetatrack_)
00107         {
00108           filterproduct->addObject(trigger::TriggerTrack, candref);
00109           n++;
00110         }
00111 
00112       // select on momentum
00113       if (filterE_){
00114         if ((candref->maxPtPxl()<maxptnearby_)&&((candref->pt())*cosh(candref->track()->eta())>minEnergy_)&&fabs(candref->track()->eta())<maxetatrack_&&fabs(candref->track()->eta())>minetatrack_)
00115           {
00116             filterproduct->addObject(trigger::TriggerTrack, candref);
00117             n++;
00118           }
00119       }
00120 
00121       // stop looping over tracks if max number is reached
00122       if(!dropMultiL2Event_ && n>=nMaxTrackCandidates_) break; 
00123 
00124     } // loop over tracks
00125   
00126   
00127   bool accept(n>0);
00128 
00129   if( dropMultiL2Event_ && n>nMaxTrackCandidates_ ) accept=false;
00130 
00131   filterproduct->addCollectionTag(candTag_);
00132 
00133   iEvent.put(filterproduct);
00134 
00135   return accept;
00136 
00137 }
00138