CMS 3D CMS Logo

HLTJetCollectionsForElePlusJets.cc
Go to the documentation of this file.
1 #include <string>
2 #include <vector>
3 
4 #include "TVector3.h"
5 
14 
15 
16 template<typename T>
18  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
19  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
20  //minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
21  //maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
22  //minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
23  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
24  //Only for VBF
25  //minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
26  //minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
27 {
28  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
29  m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
30  m_theJetToken = consumes<std::vector<T>>(sourceJetTag);
31  produces<TCollectionVector> ();
32 }
33 
34 
35 template<typename T>
37 {
38  // do anything here that needs to be done at desctruction time
39  // (e.g. close files, deallocate resources etc.)
40 
41 }
42 
43 template<typename T>
46  desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
47  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
48  // desc.add<double> ("MinJetPt", 30.);
49  // desc.add<double> ("MaxAbsJetEta", 2.6);
50  // desc.add<unsigned int> ("MinNJets", 1);
51  desc.add<double> ("minDeltaR", 0.5);
52  //Only for VBF
53  // desc.add<double> ("MinSoftJetPt", 25.);
54  //desc.add<double> ("MinDeltaEta", -1.);
56 }
57 
58 //
59 // member functions
60 //
61 
62 
63 // ------------ method called to produce the data ------------
64 template <typename T>
65 void
67 {
68  using namespace edm;
69  using namespace std;
70 
71  typedef vector<T> TCollection;
72  typedef Ref<TCollection> TRef;
73 
74  typedef edm::RefVector<TCollection> TRefVector;
75 
76  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
77 
79  iEvent.getByToken(m_theElectronToken,PrevFilterOutput);
80 
81  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
82  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
83  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
84  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
85 
86  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
87  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
88 
89  trigger::VRphoton photonCands;
90  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
91 
92  //prepare the collection of 3-D vector for electron momenta
93  std::vector<TVector3> ElePs;
94 
95  if(!clusCands.empty()){ //try trigger cluster
96  for(auto & clusCand : clusCands){
97  TVector3 positionVector(
98  clusCand->superCluster()->position().x(),
99  clusCand->superCluster()->position().y(),
100  clusCand->superCluster()->position().z());
101  ElePs.push_back(positionVector);
102  }
103  }else if(!eleCands.empty()){ // try trigger electrons
104  for(auto & eleCand : eleCands){
105  TVector3 positionVector(
106  eleCand->superCluster()->position().x(),
107  eleCand->superCluster()->position().y(),
108  eleCand->superCluster()->position().z());
109  ElePs.push_back(positionVector);
110  }
111  }
112  else if(!photonCands.empty()){ // try trigger photons
113  for(auto & photonCand : photonCands){
114  TVector3 positionVector(
115  photonCand->superCluster()->position().x(),
116  photonCand->superCluster()->position().y(),
117  photonCand->superCluster()->position().z());
118  ElePs.push_back(positionVector);
119  }
120  }
121 
122  edm::Handle<TCollection> theJetCollectionHandle;
123  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
124 
125  const TCollection & theJetCollection = *theJetCollectionHandle;
126 
127  //std::unique_ptr< TCollection > theFilteredJetCollection(new TCollection);
128 
129  std::unique_ptr < TCollectionVector > allSelections(new TCollectionVector());
130 
131  //bool foundSolution(false);
132 
133  for (auto & EleP : ElePs) {
134 
135  // bool VBFJetPair = false;
136  //std::vector<int> store_jet;
137  TRefVector refVector;
138 
139  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
140  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
141  double DR = EleP.DeltaR(JetP);
142 
143  if (DR > minDeltaR_)
144  refVector.push_back(TRef(theJetCollectionHandle, j));
145  }
146  allSelections->push_back(refVector);
147  }
148 
149  //iEvent.put(std::move(theFilteredJetCollection));
150  iEvent.put(std::move(allSelections));
151 
152  return;
153 
154 }
std::string defaultModuleLabel()
edm::EDGetTokenT< std::vector< T > > m_theJetToken
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:460
HLTJetCollectionsForElePlusJets(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:230
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
long double T
std::vector< reco::RecoEcalCandidateRef > VRphoton
virtual void produce(edm::Event &, const edm::EventSetup &)
def move(src, dest)
Definition: eostools.py:510