CMS 3D CMS Logo

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