CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
JetCollectionForEleHT.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: JetCollectionForEleHT
4 // Class: JetCollectionForEleHT
5 //
13 //
14 // Original Author: Massimiliano Chiorboli,40 4-A01,+41227671535,
15 // Created: Mon Oct 4 11:57:35 CEST 2010
16 // $Id: JetCollectionForEleHT.cc,v 1.6 2011/04/28 06:22:41 gruen Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
27 
30 
33 
37 
41 
43 
45 
47 
49  hltElectronTag(iConfig.getParameter< edm::InputTag > ("HltElectronTag")),
50  sourceJetTag(iConfig.getParameter< edm::InputTag > ("SourceJetTag")),
51  minDeltaR_(iConfig.getParameter< double > ("minDeltaR"))
52 {
53  produces<reco::CaloJetCollection>();
54 }
55 
56 
57 
59 {
60 
61  // do anything here that needs to be done at desctruction time
62  // (e.g. close files, deallocate resources etc.)
63 
64 }
65 
66 void
69  desc.add<edm::InputTag>("HltElectronTag",edm::InputTag("triggerFilterObjectWithRefs"));
70  desc.add<edm::InputTag>("SourceJetTag",edm::InputTag("caloJetCollection"));
71  desc.add<double>("minDeltaR",0.5);
72  descriptions.add("hltJetCollectionForEleHT",desc);
73 }
74 
75 //
76 // member functions
77 //
78 
79 
80 // ------------ method called to produce the data ------------
81 // template <typename T>
82 void
84 {
85  using namespace edm;
86 
88  iEvent.getByLabel(hltElectronTag,PrevFilterOutput);
89 
90  //its easier on the if statement flow if I try everything at once, shouldnt add to timing
91  std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > clusCands;
92  PrevFilterOutput->getObjects(trigger::TriggerCluster,clusCands);
93  std::vector<edm::Ref<reco::ElectronCollection> > eleCands;
94  PrevFilterOutput->getObjects(trigger::TriggerElectron,eleCands);
95 
96  //prepare the collection of 3-D vector for electron momenta
97  std::vector<TVector3> ElePs;
98 
99  if(!clusCands.empty()){ //try trigger cluster
100  for(size_t candNr=0;candNr<clusCands.size();candNr++){
101  TVector3 positionVector(
102  clusCands[candNr]->superCluster()->position().x(),
103  clusCands[candNr]->superCluster()->position().y(),
104  clusCands[candNr]->superCluster()->position().z());
105  ElePs.push_back(positionVector);
106  }
107  }else if(!eleCands.empty()){ // try trigger electrons
108  for(size_t candNr=0;candNr<eleCands.size();candNr++){
109  TVector3 positionVector(
110  eleCands[candNr]->superCluster()->position().x(),
111  eleCands[candNr]->superCluster()->position().y(),
112  eleCands[candNr]->superCluster()->position().z());
113  ElePs.push_back(positionVector);
114  }
115  }
116 
117  edm::Handle<reco::CaloJetCollection> theCaloJetCollectionHandle;
118  iEvent.getByLabel(sourceJetTag, theCaloJetCollectionHandle);
119  const reco::CaloJetCollection* theCaloJetCollection = theCaloJetCollectionHandle.product();
120 
121  std::auto_ptr< reco::CaloJetCollection > theFilteredCaloJetCollection(new reco::CaloJetCollection);
122 
123  bool isOverlapping;
124 
125  for(unsigned int j=0; j<theCaloJetCollection->size(); j++) {
126 
127  isOverlapping = false;
128  for(unsigned int i=0; i<ElePs.size(); i++) {
129 
130  TVector3 JetP((*theCaloJetCollection)[j].px(), (*theCaloJetCollection)[j].py(), (*theCaloJetCollection)[j].pz());
131  double DR = ElePs[i].DeltaR(JetP);
132 
133  if(DR<minDeltaR_) {
134  isOverlapping = true;
135  break;
136  }
137  }
138 
139  if(!isOverlapping) theFilteredCaloJetCollection->push_back((*theCaloJetCollection)[j]);
140  }
141 
142  //do the filtering
143 
144  iEvent.put(theFilteredCaloJetCollection);
145 
146  return;
147 
148 }
149 
150 // ------------ method called once each job just before starting event loop ------------
151 void
153 {
154 }
155 
156 // ------------ method called once each job just after ending the event loop ------------
157 void
159 }
160 
161 //define this as a plug-in
163 
int i
Definition: DBlmapReader.cc:9
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
JetCollectionForEleHT(const edm::ParameterSet &)
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
virtual void produce(edm::Event &, const edm::EventSetup &)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
static int position[264][3]
Definition: ReadPGInfo.cc:509
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
x
Definition: VDTMath.h:216
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects