Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-02-14 14:20:32

0001 #include "HLTrigger/JetMET/interface/HLTJetCollForElePlusJets.h"
0002 
0003 #include "FWCore/Framework/interface/MakerMacros.h"
0004 
0005 #include "DataFormats/Common/interface/Handle.h"
0006 
0007 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h"
0008 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
0009 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0010 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0011 #include "TVector3.h"
0012 #include "HLTrigger/HLTcore/interface/defaultModuleLabel.h"
0013 
0014 template <typename T>
0015 HLTJetCollForElePlusJets<T>::HLTJetCollForElePlusJets(const edm::ParameterSet& iConfig)
0016     : hltElectronTag(iConfig.getParameter<edm::InputTag>("HltElectronTag")),
0017       sourceJetTag(iConfig.getParameter<edm::InputTag>("SourceJetTag")),
0018       minJetPt_(iConfig.getParameter<double>("MinJetPt")),
0019       maxAbsJetEta_(iConfig.getParameter<double>("MaxAbsJetEta")),
0020       minNJets_(iConfig.getParameter<unsigned int>("MinNJets")),
0021       minDeltaR_(iConfig.getParameter<double>("minDeltaR")),
0022       //Only for VBF
0023       minSoftJetPt_(iConfig.getParameter<double>("MinSoftJetPt")),
0024       minDeltaEta_(iConfig.getParameter<double>("MinDeltaEta")) {
0025   typedef std::vector<T> TCollection;
0026   m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
0027   m_theJetToken = consumes<TCollection>(sourceJetTag);
0028   produces<TCollection>();
0029 }
0030 
0031 template <typename T>
0032 HLTJetCollForElePlusJets<T>::~HLTJetCollForElePlusJets() {
0033   // do anything here that needs to be done at desctruction time
0034   // (e.g. close files, deallocate resources etc.)
0035 }
0036 
0037 template <typename T>
0038 void HLTJetCollForElePlusJets<T>::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0039   edm::ParameterSetDescription desc;
0040   desc.add<edm::InputTag>("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
0041   desc.add<edm::InputTag>("SourceJetTag", edm::InputTag("jetCollection"));
0042   desc.add<double>("MinJetPt", 30.);
0043   desc.add<double>("MaxAbsJetEta", 2.6);
0044   desc.add<unsigned int>("MinNJets", 1);
0045   desc.add<double>("minDeltaR", 0.5);
0046   //Only for VBF
0047   desc.add<double>("MinSoftJetPt", 25.);
0048   desc.add<double>("MinDeltaEta", -1.);
0049   descriptions.add(defaultModuleLabel<HLTJetCollForElePlusJets<T>>(), desc);
0050 }
0051 
0052 //
0053 // member functions
0054 //
0055 
0056 // ------------ method called to produce the data  ------------
0057 template <typename T>
0058 void HLTJetCollForElePlusJets<T>::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0059   using namespace edm;
0060   using namespace std;
0061 
0062   typedef vector<T> TCollection;
0063   typedef Ref<TCollection> TRef;
0064   typedef edm::RefVector<TCollection> TRefVector;
0065   typedef std::vector<edm::RefVector<std::vector<T>, T, edm::refhelper::FindUsingAdvance<std::vector<T>, T>>>
0066       TCollectionVector;
0067 
0068   edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput;
0069   iEvent.getByToken(m_theElectronToken, PrevFilterOutput);
0070 
0071   //its easier on the if statement flow if I try everything at once, shouldnt add to timing
0072   // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
0073   std::vector<edm::Ref<reco::RecoEcalCandidateCollection>> clusCands;
0074   PrevFilterOutput->getObjects(trigger::TriggerCluster, clusCands);
0075 
0076   std::vector<edm::Ref<reco::ElectronCollection>> eleCands;
0077   PrevFilterOutput->getObjects(trigger::TriggerElectron, eleCands);
0078 
0079   trigger::VRphoton photonCands;
0080   PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
0081 
0082   //prepare the collection of 3-D vector for electron momenta
0083   std::vector<TVector3> ElePs;
0084 
0085   if (!clusCands.empty()) {  //try trigger cluster
0086     for (auto& clusCand : clusCands) {
0087       TVector3 positionVector(clusCand->superCluster()->position().x(),
0088                               clusCand->superCluster()->position().y(),
0089                               clusCand->superCluster()->position().z());
0090       ElePs.push_back(positionVector);
0091     }
0092   } else if (!eleCands.empty()) {  // try trigger electrons
0093     for (auto& eleCand : eleCands) {
0094       TVector3 positionVector(eleCand->superCluster()->position().x(),
0095                               eleCand->superCluster()->position().y(),
0096                               eleCand->superCluster()->position().z());
0097       ElePs.push_back(positionVector);
0098     }
0099   } else if (!photonCands.empty()) {  // try trigger photons
0100     for (auto& photonCand : photonCands) {
0101       TVector3 positionVector(photonCand->superCluster()->position().x(),
0102                               photonCand->superCluster()->position().y(),
0103                               photonCand->superCluster()->position().z());
0104       ElePs.push_back(positionVector);
0105     }
0106   }
0107 
0108   edm::Handle<TCollection> theJetCollectionHandle;
0109   iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
0110 
0111   const TCollection& theJetCollection = *theJetCollectionHandle;
0112 
0113   std::unique_ptr<TCollection> theFilteredJetCollection(new TCollection);
0114 
0115   std::unique_ptr<TCollectionVector> allSelections(new TCollectionVector());
0116 
0117   bool foundSolution(false);
0118 
0119   for (auto& EleP : ElePs) {
0120     bool VBFJetPair = false;
0121     std::vector<int> store_jet;
0122     TRefVector refVector;
0123 
0124     for (unsigned int j = 0; j < theJetCollection.size(); j++) {
0125       TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
0126       double DR = EleP.DeltaR(JetP);
0127 
0128       if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) {
0129         store_jet.push_back(j);
0130         // The VBF part of the filter
0131         if (minDeltaEta_ > 0) {
0132           for (unsigned int k = j + 1; k < theJetCollection.size(); k++) {
0133             TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(), theJetCollection[k].pz());
0134             double softDR = EleP.DeltaR(SoftJetP);
0135 
0136             if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_)
0137               if (std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_) {
0138                 store_jet.push_back(k);
0139                 VBFJetPair = true;
0140               }
0141           }
0142         }
0143       }
0144     }
0145 
0146     // Now remove duplicates from the jet collection to store
0147     std::sort(store_jet.begin(), store_jet.end());
0148     store_jet.erase(unique(store_jet.begin(), store_jet.end()), store_jet.end());
0149 
0150     // Now save the cleaned jets
0151     for (int& ijet : store_jet) {
0152       //store all selections
0153       refVector.push_back(TRef(theJetCollectionHandle, ijet));
0154       //store first selection which matches the criteria
0155       if (!foundSolution)
0156         theFilteredJetCollection->push_back(theJetCollection[ijet]);
0157     }
0158     //store all selections
0159     allSelections->push_back(refVector);
0160 
0161     if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
0162       foundSolution = true;
0163     else if (VBFJetPair && minDeltaEta_ > 0)
0164       foundSolution = true;
0165     else if (!foundSolution)
0166       theFilteredJetCollection->clear();
0167   }
0168 
0169   iEvent.put(std::move(theFilteredJetCollection));
0170 
0171   return;
0172 }