CMS 3D CMS Logo

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