CMS 3D CMS Logo

ElectronRegressionEnergyProducer.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // user include files
6 
12 
15 
25 
26 //
27 // class declaration
28 //
29 
31 public:
34 private:
35  bool filter(edm::Event&, const edm::EventSetup&) override;
36 
37  // ----------member data ---------------------------
40 
43 
46 
48 
51 
54 
56 
59 
60 };
61 
62 
64  printDebug_ = iConfig.getUntrackedParameter<bool>("printDebug", false);
65  electronToken_ = consumes<reco::GsfElectronCollection>(iConfig.getParameter<edm::InputTag>("electronTag"));
66 
67  regressionInputFile_ = iConfig.getParameter<std::string>("regressionInputFile");
68  energyRegressionType_ = iConfig.getParameter<uint32_t>("energyRegressionType");
69 
70  nameEnergyReg_ = iConfig.getParameter<std::string>("nameEnergyReg");
71  nameEnergyErrorReg_ = iConfig.getParameter<std::string>("nameEnergyErrorReg");
72 
73  recHitCollectionEBToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitCollectionEB"));
74  recHitCollectionEEToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitCollectionEE"));
75 
76  hVertexToken_ = consumes<reco::VertexCollection>(edm::InputTag("offlinePrimaryVertices"));
77  hRhoKt6PFJetsToken_ = consumes<double>(edm::InputTag("kt6PFJets","rho"));
78 
79  produces<edm::ValueMap<double> >(nameEnergyReg_);
80  produces<edm::ValueMap<double> >(nameEnergyErrorReg_);
81 
83 
84  //set regression type
86  if (energyRegressionType_ == 1) type = ElectronEnergyRegressionEvaluate::kNoTrkVar;
87  else if (energyRegressionType_ == 2) type = ElectronEnergyRegressionEvaluate::kWithSubCluVar;
88  else if (energyRegressionType_ == 3) type = ElectronEnergyRegressionEvaluate::kWithTrkVarV1;
89  else if (energyRegressionType_ == 4) type = ElectronEnergyRegressionEvaluate::kWithTrkVarV2;
90 
91  //load weights and initialize
92  regressionEvaluator->initialize(regressionInputFile_,type);
93 
94  geomInitialized_ = false;
95 
96 }
97 
98 
100 {
101  delete regressionEvaluator;
102 }
103 
104 // ------------ method called on each new Event ------------
106 
108 
109  if (!geomInitialized_) {
110  edm::ESHandle<CaloTopology> theCaloTopology;
111  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
112  ecalTopology_ = & (*theCaloTopology);
113 
114  edm::ESHandle<CaloGeometry> theCaloGeometry;
115  iSetup.get<CaloGeometryRecord>().get(theCaloGeometry);
116  caloGeometry_ = & (*theCaloGeometry);
117  geomInitialized_ = true;
118  }
119 
120  std::unique_ptr<edm::ValueMap<double> > regrEnergyMap(new edm::ValueMap<double>() );
121  edm::ValueMap<double>::Filler energyFiller(*regrEnergyMap);
122 
123  std::unique_ptr<edm::ValueMap<double> > regrEnergyErrorMap(new edm::ValueMap<double>() );
124  edm::ValueMap<double>::Filler energyErrorFiller(*regrEnergyErrorMap);
125 
127  iEvent.getByToken(electronToken_,egCollection);
128  const reco::GsfElectronCollection egCandidates = (*egCollection.product());
129 
130  std::vector<double> energyValues;
131  std::vector<double> energyErrorValues;
132  energyValues.reserve(egCollection->size());
133  energyErrorValues.reserve(egCollection->size());
134  //
135  //**************************************************************************
136  // Rechits
137  //**************************************************************************
140  iEvent.getByToken( recHitCollectionEBToken_, pEBRecHits );
141  iEvent.getByToken( recHitCollectionEEToken_, pEERecHits );
142 
143  //**************************************************************************
144  //Get Number of Vertices
145  //**************************************************************************
147  iEvent.getByToken(hVertexToken_,hVertexProduct);
148  const reco::VertexCollection inVertices = *(hVertexProduct.product());
149 
150  // loop through all vertices
151  Int_t nvertices = 0;
152  for (reco::VertexCollection::const_iterator inV = inVertices.begin();
153  inV != inVertices.end(); ++inV) {
154 
155  // pass these vertex cuts
156  if (inV->ndof() >= 4
157  && inV->position().Rho() <= 2.0
158  && fabs(inV->z()) <= 24.0
159  ) {
160  nvertices++;
161  }
162  }
163 
164  //**************************************************************************
165  //Get Rho
166  //**************************************************************************
167  double rho = 0;
168  edm::Handle<double> hRhoKt6PFJets;
169  iEvent.getByToken(hRhoKt6PFJetsToken_, hRhoKt6PFJets);
170  rho = (*hRhoKt6PFJets);
171 
172 
173  for ( reco::GsfElectronCollection::const_iterator egIter = egCandidates.begin();
174  egIter != egCandidates.end(); ++egIter) {
175 
176  const EcalRecHitCollection * recHits=nullptr;
177  if(egIter->isEB())
178  recHits = pEBRecHits.product();
179  else
180  recHits = pEERecHits.product();
181 
182  SuperClusterHelper mySCHelper(&(*egIter),recHits,ecalTopology_,caloGeometry_);
183 
184  double energy=regressionEvaluator->calculateRegressionEnergy(&(*egIter),
185  mySCHelper,
186  rho,nvertices,
187  printDebug_);
188 
190  mySCHelper,
191  rho,nvertices,
192  printDebug_);
193 
194  energyValues.push_back(energy);
195  energyErrorValues.push_back(error);
196 
197  }
198 
199  energyFiller.insert( egCollection, energyValues.begin(), energyValues.end() );
200  energyFiller.fill();
201 
202  energyErrorFiller.insert( egCollection, energyErrorValues.begin(), energyErrorValues.end() );
203  energyErrorFiller.fill();
204 
205  iEvent.put(std::move(regrEnergyMap),nameEnergyReg_);
206  iEvent.put(std::move(regrEnergyErrorMap),nameEnergyErrorReg_);
207 
208  return true;
209 
210 }
211 
212 
213 //define this as a plug-in
215 
216 
217 
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
ElectronEnergyRegressionEvaluate * regressionEvaluator
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:53
double calculateRegressionEnergy(const reco::GsfElectron *ele, SuperClusterHelper &mySCHelper, double rho, double nvertices, bool printDebug=false)
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEEToken_
std::vector< GsfElectron > GsfElectronCollection
collection of GsfElectron objects
int iEvent
Definition: GenABIO.cc:230
ElectronRegressionEnergyProducer(const edm::ParameterSet &)
bool filter(edm::Event &, const edm::EventSetup &) override
edm::EDGetTokenT< reco::VertexCollection > hVertexToken_
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEBToken_
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:59
double calculateRegressionEnergyUncertainty(const reco::GsfElectron *ele, SuperClusterHelper &mySCHelper, double rho, double nvertices, bool printDebug=false)
edm::EDGetTokenT< reco::GsfElectronCollection > electronToken_
def move(src, dest)
Definition: eostools.py:510