CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch1/src/HLTrigger/JetMET/src/HLTJetCollForElePlusJets.cc

Go to the documentation of this file.
00001 #include "HLTrigger/JetMET/interface/HLTJetCollForElePlusJets.h"
00002 
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 
00005 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00006 #include "DataFormats/EgammaCandidates/interface/Electron.h"
00007 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00008 #include "DataFormats/JetReco/interface/PFJetCollection.h"
00009 
00010 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
00011 
00012 #include "DataFormats/Common/interface/Handle.h"
00013 
00014 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00015 
00016 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00017 #include "TVector3.h"
00018 
00019 #include<typeinfo>
00020 
00021 template <typename T>
00022 HLTJetCollForElePlusJets<T>::HLTJetCollForElePlusJets(const edm::ParameterSet& iConfig):
00023   hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
00024   sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
00025   minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
00026   maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
00027   minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
00028   minDeltaR_(iConfig.getParameter< double > ("minDeltaR")),
00029   //Only for VBF
00030   minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
00031   minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
00032 {
00033   typedef std::vector<T> TCollection;
00034   produces<TCollection>();
00035 }
00036 
00037 
00038 template <typename T>
00039 HLTJetCollForElePlusJets<T>::~HLTJetCollForElePlusJets()
00040 {
00041    // do anything here that needs to be done at desctruction time
00042    // (e.g. close files, deallocate resources etc.)
00043 
00044 }
00045 
00046 template <typename T>
00047 void HLTJetCollForElePlusJets<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00048     edm::ParameterSetDescription desc;
00049     desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
00050     desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
00051     desc.add<double> ("MinJetPt", 30.);
00052     desc.add<double> ("MaxAbsJetEta", 2.6);
00053     desc.add<unsigned int> ("MinNJets", 1);
00054     desc.add<double> ("minDeltaR", 0.5);
00055     //Only for VBF
00056     desc.add<double> ("MinSoftJetPt", 25.);
00057     desc.add<double> ("MinDeltaEta", -1.);    
00058     descriptions.add(std::string("hlt")+std::string(typeid(HLTJetCollForElePlusJets<T>).name()), desc);
00059 }
00060 
00061 //
00062 // member functions
00063 //
00064 
00065 
00066 // ------------ method called to produce the data  ------------
00067 template <typename T>
00068 void
00069 HLTJetCollForElePlusJets<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00070 {
00071   using namespace edm;
00072   using namespace std;
00073 
00074   typedef vector<T> TCollection;
00075   typedef Ref<TCollection> TRef;
00076   typedef edm::RefVector<TCollection> TRefVector;
00077   typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
00078 
00079   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
00080   iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
00081  
00082   //its easier on the if statement flow if I try everything at once, shouldnt add to timing
00083   std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
00084   PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
00085 
00086   std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
00087   PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
00088   
00089   //prepare the collection of 3-D vector for electron momenta
00090   std::vector<TVector3> ElePs;
00091 
00092   if(!clusCands.empty()){ //try trigger cluster
00093     for(size_t candNr=0;candNr<clusCands.size();candNr++){
00094       TVector3 positionVector(
00095                               clusCands[candNr]->superCluster()->position().x(),
00096                               clusCands[candNr]->superCluster()->position().y(),
00097                               clusCands[candNr]->superCluster()->position().z());
00098       ElePs.push_back(positionVector);
00099     }
00100   }else if(!eleCands.empty()){ // try trigger electrons
00101     for(size_t candNr=0;candNr<eleCands.size();candNr++){
00102       TVector3 positionVector(
00103                               eleCands[candNr]->superCluster()->position().x(),
00104                               eleCands[candNr]->superCluster()->position().y(),
00105                               eleCands[candNr]->superCluster()->position().z());
00106       ElePs.push_back(positionVector);
00107     }
00108   }
00109   
00110   edm::Handle<TCollection> theJetCollectionHandle;
00111   iEvent.getByLabel(sourceJetTag, theJetCollectionHandle);
00112   
00113   const TCollection & theJetCollection = *theJetCollectionHandle;
00114   
00115   std::auto_ptr< TCollection >  theFilteredJetCollection(new TCollection);
00116   
00117   std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
00118   
00119   bool foundSolution(false);
00120 
00121   for (unsigned int i = 0; i < ElePs.size(); i++) {
00122     
00123     bool VBFJetPair = false;
00124     std::vector<int> store_jet;
00125     TRefVector refVector;
00126     
00127     for (unsigned int j = 0; j < theJetCollection.size(); j++) {
00128       TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(),
00129                     theJetCollection[j].pz());
00130       double DR = ElePs[i].DeltaR(JetP);
00131       
00132       if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) {
00133         store_jet.push_back(j);
00134         // The VBF part of the filter
00135         if ( minDeltaEta_ > 0 ) {
00136           for ( unsigned int k = j+1; k < theJetCollection.size(); k++ ) {
00137             TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(),
00138                               theJetCollection[k].pz());
00139             double softDR = ElePs[i].DeltaR(SoftJetP);
00140             
00141             if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_)
00142               if ( std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_ ) {
00143                 store_jet.push_back(k);
00144                 VBFJetPair = true;
00145               }
00146           }
00147         }
00148       }
00149       
00150     }
00151     
00152     // Now remove duplicates from the jet collection to store
00153     std::sort( store_jet.begin(), store_jet.end() );
00154     store_jet.erase( unique( store_jet.begin(), store_jet.end() ), store_jet.end() );
00155     
00156     // Now save the cleaned jets
00157     for ( unsigned int ijet = 0; ijet < store_jet.size(); ijet++ )
00158       {
00159         //store all selections
00160         refVector.push_back(TRef(theJetCollectionHandle, store_jet.at(ijet)));
00161         //store first selection which matches the criteria
00162         if(!foundSolution)
00163           theFilteredJetCollection->push_back(theJetCollection[store_jet.at(ijet)]);
00164       }
00165     //store all selections
00166     allSelections->push_back(refVector);
00167     
00168     if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
00169       foundSolution = true;
00170     else if (VBFJetPair && minDeltaEta_ > 0)
00171       foundSolution = true;
00172     else if (!foundSolution)
00173       theFilteredJetCollection->clear();
00174     
00175     
00176   }
00177   
00178   iEvent.put(theFilteredJetCollection);
00179   
00180   return;
00181   
00182 }