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 
9 
11 
13 
15 
17 #include "TVector3.h"
18 
19 #include <string>
20 #include <vector>
21 #include <typeinfo>
22 
23 template<typename T>
25  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
26  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
27  //minJetPt_(iConfig.getParameter<double> ("MinJetPt")),
28  //maxAbsJetEta_(iConfig.getParameter<double> ("MaxAbsJetEta")),
29  //minNJets_(iConfig.getParameter<unsigned int> ("MinNJets")),
30  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
31  //Only for VBF
32  //minSoftJetPt_(iConfig.getParameter< double > ("MinSoftJetPt")),
33  //minDeltaEta_(iConfig.getParameter< double > ("MinDeltaEta"))
34 {
35  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
36  produces<TCollectionVector> ();
37 }
38 
39 
40 template<typename T>
42 {
43  // do anything here that needs to be done at desctruction time
44  // (e.g. close files, deallocate resources etc.)
45 
46 }
47 
48 template<typename T>
51  desc.add<edm::InputTag> ("HltElectronTag", edm::InputTag("triggerFilterObjectWithRefs"));
52  desc.add<edm::InputTag> ("SourceJetTag", edm::InputTag("jetCollection"));
53  // desc.add<double> ("MinJetPt", 30.);
54  // desc.add<double> ("MaxAbsJetEta", 2.6);
55  // desc.add<unsigned int> ("MinNJets", 1);
56  desc.add<double> ("minDeltaR", 0.5);
57  //Only for VBF
58  // desc.add<double> ("MinSoftJetPt", 25.);
59  //desc.add<double> ("MinDeltaEta", -1.);
60  descriptions.add(std::string("hlt")+std::string(typeid(HLTJetCollectionsForElePlusJets<T>).name()),desc);
61 }
62 
63 //
64 // member functions
65 //
66 
67 
68 // ------------ method called to produce the data ------------
69 template <typename T>
70 void
72 {
73  using namespace edm;
74  using namespace std;
75 
76  typedef vector<T> TCollection;
77  typedef Ref<TCollection> TRef;
78 
79  typedef edm::RefVector<TCollection> TRefVector;
80 
81  typedef std::vector<edm::RefVector<std::vector<T>,T,edm::refhelper::FindUsingAdvance<std::vector<T>,T> > > TCollectionVector;
82 
84  iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
85 
86  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
87  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
88  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
89 
90  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
91  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
92 
93  //prepare the collection of 3-D vector for electron momenta
94  std::vector<TVector3> ElePs;
95 
96  if(!clusCands.empty()){ //try trigger cluster
97  for(size_t candNr=0;candNr<clusCands.size();candNr++){
98  TVector3 positionVector(
99  clusCands[candNr]->superCluster()->position().x(),
100  clusCands[candNr]->superCluster()->position().y(),
101  clusCands[candNr]->superCluster()->position().z());
102  ElePs.push_back(positionVector);
103  }
104  }else if(!eleCands.empty()){ // try trigger electrons
105  for(size_t candNr=0;candNr<eleCands.size();candNr++){
106  TVector3 positionVector(
107  eleCands[candNr]->superCluster()->position().x(),
108  eleCands[candNr]->superCluster()->position().y(),
109  eleCands[candNr]->superCluster()->position().z());
110  ElePs.push_back(positionVector);
111  }
112  }
113 
114  edm::Handle<TCollection> theJetCollectionHandle;
115  iEvent.getByLabel(sourceJetTag, theJetCollectionHandle);
116  //const TCollection* theJetCollection = theJetCollectionHandle.product();
117 
118  const TCollection & theJetCollection = *theJetCollectionHandle;
119 
120  //std::auto_ptr< TCollection > theFilteredJetCollection(new TCollection);
121 
122  std::auto_ptr < TCollectionVector > allSelections(new TCollectionVector());
123 
124  //bool foundSolution(false);
125 
126  for (unsigned int i = 0; i < ElePs.size(); i++) {
127 
128  // bool VBFJetPair = false;
129  //std::vector<int> store_jet;
130  TRefVector refVector;
131 
132  for (unsigned int j = 0; j < theJetCollection.size(); j++) {
133  TVector3 JetP(theJetCollection[j].px(), theJetCollection[j].py(), theJetCollection[j].pz());
134  double DR = ElePs[i].DeltaR(JetP);
135 
136  if (DR > minDeltaR_)
137  refVector.push_back(TRef(theJetCollectionHandle, j));
138  }
139  allSelections->push_back(refVector);
140  }
141 
142  //iEvent.put(theFilteredJetCollection);
143  iEvent.put(allSelections);
144 
145  return;
146 
147 }
int i
Definition: DBlmapReader.cc:9
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
HLTJetCollectionsForElePlusJets(const edm::ParameterSet &)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
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:361
void add(std::string const &label, ParameterSetDescription const &psetDescription)
Definition: DDAxes.h:10
long double T
virtual void produce(edm::Event &, const edm::EventSetup &)