CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HLTJetCollectionsForElePlusJets.cc
Go to the documentation of this file.
2 
4 
7 
9 
11 
13 #include "TVector3.h"
14 
15 #include <string>
16 #include <vector>
17 #include <typeinfo>
18 
19 template<typename T>
21  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
22  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
23  //minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
24  //maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
25  //minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
26  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
27  //Only for VBF
28  //minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
29  //minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
30 {
31  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
32  m_theElectronToken = consumes<trigger::TriggerFilterObjectWithRefs>(hltElectronTag);
33  m_theJetToken = consumes<std::vector<T>>(sourceJetTag);
34  produces<TCollectionVector> ();
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(HLTJetCollectionsForElePlusJets<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 
77  typedef edm::RefVector<TCollection> TRefVector;
78 
79  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
80 
82  iEvent.getByToken(m_theElectronToken,PrevFilterOutput);
83 
84  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
85  // Electrons can be stored as objects of types TriggerCluster, TriggerElectron, or TriggerPhoton
86  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
87  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
88 
89  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
90  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
91 
92  trigger::VRphoton photonCands;
93  PrevFilterOutput->getObjects(trigger::TriggerPhoton, photonCands);
94 
95  //prepare the collection of 3-D vector for electron momenta
96  std::vector<TVector3> ElePs;
97 
98  if(!clusCands.empty()){ //try trigger cluster
99  for(size_t candNr=0;candNr<clusCands.size();candNr++){
100  TVector3 positionVector(
101  clusCands[candNr]->superCluster()->position().x(),
102  clusCands[candNr]->superCluster()->position().y(),
103  clusCands[candNr]->superCluster()->position().z());
104  ElePs.push_back(positionVector);
105  }
106  }else if(!eleCands.empty()){ // try trigger electrons
107  for(size_t candNr=0;candNr<eleCands.size();candNr++){
108  TVector3 positionVector(
109  eleCands[candNr]->superCluster()->position().x(),
110  eleCands[candNr]->superCluster()->position().y(),
111  eleCands[candNr]->superCluster()->position().z());
112  ElePs.push_back(positionVector);
113  }
114  }
115  else if(!photonCands.empty()){ // try trigger photons
116  for(size_t candNr=0;candNr<photonCands.size();candNr++){
117  TVector3 positionVector(
118  photonCands[candNr]->superCluster()->position().x(),
119  photonCands[candNr]->superCluster()->position().y(),
120  photonCands[candNr]->superCluster()->position().z());
121  ElePs.push_back(positionVector);
122  }
123  }
124 
125  edm::Handle<TCollection> theJetCollectionHandle;
126  iEvent.getByToken(m_theJetToken, theJetCollectionHandle);
127 
128  const TCollection & theJetCollection = *theJetCollectionHandle;
129 
130  //std::auto_ptr< TCollection > theFilteredJetCollection(new TCollection);
131 
132  std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
133 
134  //bool foundSolution(false);
135 
136  for (unsigned int i = 0; i < ElePs.size(); i++) {
137 
138  // bool VBFJetPair = false;
139  //std::vector<int> store_jet;
140  TRefVector refVector;
141 
142  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
143  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
144  double DR = ElePs[i].DeltaR(JetP);
145 
146  if (DR > minDeltaR_)
147  refVector.push_back(TRef(theJetCollectionHandle, j));
148  }
149  allSelections->push_back(refVector);
150  }
151 
152  //iEvent.put(theFilteredJetCollection);
153  iEvent.put(allSelections);
154 
155  return;
156 
157 }
edm::EDGetTokenT< std::vector< T > > m_theJetToken
int i
Definition: DBlmapReader.cc:9
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
HLTJetCollectionsForElePlusJets(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:116
int j
Definition: DBlmapReader.cc:9
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< trigger::TriggerFilterObjectWithRefs > m_theElectronToken
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: DDAxes.h:10
long double T
std::vector< reco::RecoEcalCandidateRef > VRphoton
virtual void produce(edm::Event &, const edm::EventSetup &)