CMS 3D CMS Logo

HLTJetCollForElePlusJets.cc
Go to the documentation of this file.
2 
4 
6 
11 #include "TVector3.h"
13 
14 template <typename T>
16  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
17  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
18  minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
19  maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
20  minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
21  minDeltaR_(iConfig.getParameter< double > ("minDeltaR")),
22  //Only for VBF
23  minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
24  minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
25 {
26  typedef std::vector<T> TCollection;
27  m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
28  m_theJetToken = consumes<TCollection>(sourceJetTag);
29  produces<TCollection>();
30 }
31 
32 
33 template <typename T>
35 {
36  // do anything here that needs to be done at desctruction time
37  // (e.g. close files, deallocate resources etc.)
38 
39 }
40 
41 template <typename T>
44  desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
45  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
46  desc.add<double> ("MinJetPt", 30.);
47  desc.add<double> ("MaxAbsJetEta", 2.6);
48  desc.add<unsigned int> ("MinNJets", 1);
49  desc.add<double> ("minDeltaR", 0.5);
50  //Only for VBF
51  desc.add<double> ("MinSoftJetPt", 25.);
52  desc.add<double> ("MinDeltaEta", -1.);
53  descriptions.add(defaultModuleLabel<HLTJetCollForElePlusJets<T>>(), desc);
54 }
55 
56 //
57 // member functions
58 //
59 
60 
61 // ------------ method called to produce the data ------------
62 template <typename T>
63 void
65 {
66  using namespace edm;
67  using namespace std;
68 
69  typedef vector<T> TCollection;
70  typedef Ref<TCollection> TRef;
71  typedef edm::RefVector<TCollection> TRefVector;
72  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
73 
75  iEvent.getByToken(m_theElectronToken,PrevFilterOutput);
76 
77  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
78  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
79  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
80  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
81 
82  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
83  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
84 
85  trigger::VRphoton photonCands;
86  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
87 
88  //prepare the collection of 3-D vector for electron momenta
89  std::vector<TVector3> ElePs;
90 
91  if(!clusCands.empty()){ //try trigger cluster
92  for(auto & clusCand : clusCands){
93  TVector3 positionVector(
94  clusCand->superCluster()->position().x(),
95  clusCand->superCluster()->position().y(),
96  clusCand->superCluster()->position().z());
97  ElePs.push_back(positionVector);
98  }
99  }else if(!eleCands.empty()){ // try trigger electrons
100  for(auto & eleCand : eleCands){
101  TVector3 positionVector(
102  eleCand->superCluster()->position().x(),
103  eleCand->superCluster()->position().y(),
104  eleCand->superCluster()->position().z());
105  ElePs.push_back(positionVector);
106  }
107  }
108  else if(!photonCands.empty()){ // try trigger photons
109  for(auto & photonCand : photonCands){
110  TVector3 positionVector(
111  photonCand->superCluster()->position().x(),
112  photonCand->superCluster()->position().y(),
113  photonCand->superCluster()->position().z());
114  ElePs.push_back(positionVector);
115  }
116  }
117 
118  edm::Handle<TCollection> theJetCollectionHandle;
119  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
120 
121  const TCollection & theJetCollection = *theJetCollectionHandle;
122 
123  std::unique_ptr< TCollection > theFilteredJetCollection(new TCollection);
124 
125  std::unique_ptr < TCollectionVector > allSelections(new TCollectionVector());
126 
127  bool foundSolution(false);
128 
129  for (auto & EleP : ElePs) {
130 
131  bool VBFJetPair = false;
132  std::vector<int> store_jet;
133  TRefVector refVector;
134 
135  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
136  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(),
137  theJetCollection[j].pz());
138  double DR = EleP.DeltaR(JetP);
139 
140  if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) {
141  store_jet.push_back(j);
142  // The VBF part of the filter
143  if ( minDeltaEta_ > 0 ) {
144  for ( unsigned int k = j+1; k < theJetCollection.size(); k++ ) {
145  TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(),
146  theJetCollection[k].pz());
147  double softDR = EleP.DeltaR(SoftJetP);
148 
149  if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_)
150  if ( std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_ ) {
151  store_jet.push_back(k);
152  VBFJetPair = true;
153  }
154  }
155  }
156  }
157 
158  }
159 
160  // Now remove duplicates from the jet collection to store
161  std::sort( store_jet.begin(), store_jet.end() );
162  store_jet.erase( unique( store_jet.begin(), store_jet.end() ), store_jet.end() );
163 
164  // Now save the cleaned jets
165  for (int & ijet : store_jet)
166  {
167  //store all selections
168  refVector.push_back(TRef(theJetCollectionHandle, ijet));
169  //store first selection which matches the criteria
170  if(!foundSolution)
171  theFilteredJetCollection->push_back(theJetCollection[ijet]);
172  }
173  //store all selections
174  allSelections->push_back(refVector);
175 
176  if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
177  foundSolution = true;
178  else if (VBFJetPair && minDeltaEta_ > 0)
179  foundSolution = true;
180  else if (!foundSolution)
181  theFilteredJetCollection->clear();
182 
183 
184  }
185 
186  iEvent.put(std::move(theFilteredJetCollection));
187 
188  return;
189 
190 }
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
void getObjects(Vids &ids, VRphoton &refs) const
various physics-level getters:
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
edm::EDGetTokenT< std::vector< T > > m_theJetToken
std::string defaultModuleLabel()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
int iEvent
Definition: GenABIO.cc:230
def unique(seq, keepstr=True)
Definition: tier0.py:24
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
HLTJetCollForElePlusJets(const edm::ParameterSet &)
int k[5][pyjets_maxn]
void add(std::string const &label, ParameterSetDescription const &psetDescription)
void produce(edm::Event &, const edm::EventSetup &) override
HLT enums.
long double T
std::vector< reco::RecoEcalCandidateRef > VRphoton
def move(src, dest)
Definition: eostools.py:510