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 
9 
11 
13 
15 
17 #include "TVector3.h"
18 
19 #include<typeinfo>
20 
21 template <typename T>
23  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
24  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
25  minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
26  maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
27  minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
28  minDeltaR_(iConfig.getParameter< double > ("minDeltaR")),
29  //Only for VBF
30  minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
31  minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
32 {
33  typedef std::vector<T> TCollection;
34  produces<TCollection>();
35 }
36 
37 
38 template <typename T>
40 {
41  // do anything here that needs to be done at desctruction time
42  // (e.g. close files, deallocate resources etc.)
43 
44 }
45 
46 template <typename T>
49  desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
50  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
51  desc.add<double> ("MinJetPt", 30.);
52  desc.add<double> ("MaxAbsJetEta", 2.6);
53  desc.add<unsigned int> ("MinNJets", 1);
54  desc.add<double> ("minDeltaR", 0.5);
55  //Only for VBF
56  desc.add<double> ("MinSoftJetPt", 25.);
57  desc.add<double> ("MinDeltaEta", -1.);
58  descriptions.add(std::string("hlt")+std::string(typeid(HLTJetCollForElePlusJets<T>).name()), desc);
59 }
60 
61 //
62 // member functions
63 //
64 
65 
66 // ------------ method called to produce the data ------------
67 template <typename T>
68 void
70 {
71  using namespace edm;
72  using namespace std;
73 
74  typedef vector<T> TCollection;
75  typedef Ref<TCollection> TRef;
76  typedef edm::RefVector<TCollection> TRefVector;
77  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
78 
80  iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
81 
82  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
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  //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 
110  edm::Handle<TCollection> theJetCollectionHandle;
111  iEvent.getByLabel(sourceJetTag, theJetCollectionHandle);
112 
113  const TCollection & theJetCollection = *theJetCollectionHandle;
114 
115  std::auto_ptr< TCollection > theFilteredJetCollection(new TCollection);
116 
117  std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
118 
119  bool foundSolution(false);
120 
121  for (unsigned int i = 0; i < ElePs.size(); i++) {
122 
123  bool VBFJetPair = false;
124  std::vector<int> store_jet;
125  TRefVector refVector;
126 
127  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
128  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(),
129  theJetCollection[j].pz());
130  double DR = ElePs[i].DeltaR(JetP);
131 
132  if (JetP.Pt() > minJetPt_ && std::abs(JetP.Eta()) < maxAbsJetEta_ && DR > minDeltaR_) {
133  store_jet.push_back(j);
134  // The VBF part of the filter
135  if ( minDeltaEta_ > 0 ) {
136  for ( unsigned int k = j+1; k < theJetCollection.size(); k++ ) {
137  TVector3 SoftJetP(theJetCollection[k].px(), theJetCollection[k].py(),
138  theJetCollection[k].pz());
139  double softDR = ElePs[i].DeltaR(SoftJetP);
140 
141  if (SoftJetP.Pt() > minSoftJetPt_ && std::abs(SoftJetP.Eta()) < maxAbsJetEta_ && softDR > minDeltaR_)
142  if ( std::abs(SoftJetP.Eta() - JetP.Eta()) > minDeltaEta_ ) {
143  store_jet.push_back(k);
144  VBFJetPair = true;
145  }
146  }
147  }
148  }
149 
150  }
151 
152  // Now remove duplicates from the jet collection to store
153  std::sort( store_jet.begin(), store_jet.end() );
154  store_jet.erase( unique( store_jet.begin(), store_jet.end() ), store_jet.end() );
155 
156  // Now save the cleaned jets
157  for ( unsigned int ijet = 0; ijet < store_jet.size(); ijet++ )
158  {
159  //store all selections
160  refVector.push_back(TRef(theJetCollectionHandle, store_jet.at(ijet)));
161  //store first selection which matches the criteria
162  if(!foundSolution)
163  theFilteredJetCollection->push_back(theJetCollection[store_jet.at(ijet)]);
164  }
165  //store all selections
166  allSelections->push_back(refVector);
167 
168  if (theFilteredJetCollection->size() >= minNJets_ && minDeltaEta_ < 0)
169  foundSolution = true;
170  else if (VBFJetPair && minDeltaEta_ > 0)
171  foundSolution = true;
172  else if (!foundSolution)
173  theFilteredJetCollection->clear();
174 
175 
176  }
177 
178  iEvent.put(theFilteredJetCollection);
179 
180  return;
181 
182 }
int i
Definition: DBlmapReader.cc:9
#define abs(x)
Definition: mlp_lapack.h:159
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
double double double z
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
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 &)
x
Definition: VDTMath.h:216
long double T