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  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  //prepare the collection of 3-D vector for electron momenta
86  std::vector<TVector3> ElePs;
87 
88  if(!clusCands.empty()){ //try trigger cluster
89  for(size_t candNr=0;candNr<clusCands.size();candNr++){
90  TVector3 positionVector(
91  clusCands[candNr]->superCluster()->position().x(),
92  clusCands[candNr]->superCluster()->position().y(),
93  clusCands[candNr]->superCluster()->position().z());
94  ElePs.push_back(positionVector);
95  }
96  }else if(!eleCands.empty()){ // try trigger electrons
97  for(size_t candNr=0;candNr<eleCands.size();candNr++){
98  TVector3 positionVector(
99  eleCands[candNr]->superCluster()->position().x(),
100  eleCands[candNr]->superCluster()->position().y(),
101  eleCands[candNr]->superCluster()->position().z());
102  ElePs.push_back(positionVector);
103  }
104  }
105 
106  edm::Handle<TCollection> theJetCollectionHandle;
107  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
108 
109  const TCollection & theJetCollection = *theJetCollectionHandle;
110 
111  std::auto_ptr< TCollection > theFilteredJetCollection(new TCollection);
112 
113  std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
114 
115  bool foundSolution(false);
116 
117  for (unsigned int i = 0; i < ElePs.size(); i++) {
118 
119  bool VBFJetPair = false;
120  std::vector<int> store_jet;
121  TRefVector refVector;
122 
123  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
124  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(),
125  theJetCollection[j].pz());
126  double DR = ElePs[i].DeltaR(JetP);
127 
128  if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) {
129  store_jet.push_back(j);
130  // The VBF part of the filter
131  if ( minDeltaEta_ > 0 ) {
132  for ( unsigned int k = j+1; k < theJetCollection.size(); k++ ) {
133  TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(),
134  theJetCollection[k].pz());
135  double softDR = ElePs[i].DeltaR(SoftJetP);
136 
137  if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_)
138  if ( std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_ ) {
139  store_jet.push_back(k);
140  VBFJetPair = true;
141  }
142  }
143  }
144  }
145 
146  }
147 
148  // Now remove duplicates from the jet collection to store
149  std::sort( store_jet.begin(), store_jet.end() );
150  store_jet.erase( unique( store_jet.begin(), store_jet.end() ), store_jet.end() );
151 
152  // Now save the cleaned jets
153  for ( unsigned int ijet = 0; ijet < store_jet.size(); ijet++ )
154  {
155  //store all selections
156  refVector.push_back(TRef(theJetCollectionHandle, store_jet.at(ijet)));
157  //store first selection which matches the criteria
158  if(!foundSolution)
159  theFilteredJetCollection->push_back(theJetCollection[store_jet.at(ijet)]);
160  }
161  //store all selections
162  allSelections->push_back(refVector);
163 
164  if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
165  foundSolution = true;
166  else if (VBFJetPair && minDeltaEta_ > 0)
167  foundSolution = true;
168  else if (!foundSolution)
169  theFilteredJetCollection->clear();
170 
171 
172  }
173 
174  iEvent.put(theFilteredJetCollection);
175 
176  return;
177 
178 }
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
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
int iEvent
Definition: GenABIO.cc:243
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)
virtual void produce(edm::Event &, const edm::EventSetup &)
Definition: DDAxes.h:10
long double T