CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
PhotonEnergyCorrector Class Reference

#include <PhotonEnergyCorrector.h>

Public Member Functions

double applyCrackCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
 
void calculate (edm::Event &evt, reco::Photon &, int subdet, const reco::VertexCollection &vtxcol, const edm::EventSetup &iSetup)
 
std::unique_ptr
< ModifyObjectValueBase > & 
gedRegression ()
 
void init (const edm::EventSetup &theEventSetup)
 
 PhotonEnergyCorrector (const edm::ParameterSet &config, edm::ConsumesCollector &&iC)
 

Private Attributes

edm::InputTag barrelEcalHits_
 
edm::EDGetTokenT
< EcalRecHitCollection
barrelEcalHitsToken_
 
std::string candidateP4type_
 
edm::InputTag endcapEcalHits_
 
edm::EDGetTokenT
< EcalRecHitCollection
endcapEcalHitsToken_
 
std::unique_ptr
< ModifyObjectValueBase
gedRegression_
 
double minR9Barrel_
 
double minR9Endcap_
 
std::unique_ptr
< EcalClusterFunctionBaseClass
photonEcalEnergyCorrFunction_
 
std::unique_ptr
< EnergyUncertaintyPhotonSpecific
photonUncertaintyCalculator_
 
std::unique_ptr
< EGEnergyCorrector
regressionCorrector_
 
std::unique_ptr
< EcalClusterFunctionBaseClass
scCrackEnergyFunction_
 
std::unique_ptr
< EcalClusterFunctionBaseClass
scEnergyErrorFunction_
 
std::unique_ptr
< EcalClusterFunctionBaseClass
scEnergyFunction_
 
edm::ESHandle< CaloGeometrytheCaloGeom_
 
std::string w_db_
 
std::string w_file_
 
bool weightsfromDB_
 

Detailed Description

Author
Nancy Marinelli, U. of Notre Dame, US

Definition at line 27 of file PhotonEnergyCorrector.h.

Constructor & Destructor Documentation

PhotonEnergyCorrector::PhotonEnergyCorrector ( const edm::ParameterSet config,
edm::ConsumesCollector &&  iC 
)

Definition at line 13 of file PhotonEnergyCorrector.cc.

References barrelEcalHits_, barrelEcalHitsToken_, SurfaceDeformationFactory::create(), endcapEcalHits_, endcapEcalHitsToken_, edm::ParameterSet::existsAs(), gedRegression_, reco::get(), edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterSet(), minR9Barrel_, minR9Endcap_, photonEcalEnergyCorrFunction_, photonUncertaintyCalculator_, regressionCorrector_, scCrackEnergyFunction_, scEnergyErrorFunction_, scEnergyFunction_, AlCaHLTBitMon_QueryRunRegistry::string, w_db_, w_file_, and weightsfromDB_.

13  {
14 
15  minR9Barrel_ = config.getParameter<double>("minR9Barrel");
16  minR9Endcap_ = config.getParameter<double>("minR9Endcap");
17  // get the geometry from the event setup:
18 
19  barrelEcalHits_ = config.getParameter<edm::InputTag>("barrelEcalHits");
20  endcapEcalHits_ = config.getParameter<edm::InputTag>("endcapEcalHits");
23 
24  // candidateP4type_ = config.getParameter<std::string>("candidateP4type") ;
25 
26 
27  // function to extract f(eta) correction
28  scEnergyFunction_ = 0 ;
29  std::string superClusterFunctionName = config.getParameter<std::string>("superClusterEnergyCorrFunction") ;
30  scEnergyFunction_.reset(EcalClusterFunctionFactory::get()->create(superClusterFunctionName,config));
31 
32 
33  // function to extract corrections to cracks
35  std::string superClusterCrackFunctionName = config.getParameter<std::string>("superClusterCrackEnergyCorrFunction") ;
36  scCrackEnergyFunction_.reset(EcalClusterFunctionFactory::get()->create(superClusterCrackFunctionName,config));
37 
38 
39  // function to extract the error on the sc ecal correction
41  std::string superClusterErrorFunctionName = config.getParameter<std::string>("superClusterEnergyErrorFunction") ;
42  scEnergyErrorFunction_.reset(EcalClusterFunctionFactory::get()->create(superClusterErrorFunctionName,config));
43 
44 
45  // function to extract the error on the photon ecal correction
47  std::string photonEnergyFunctionName = config.getParameter<std::string>("photonEcalEnergyCorrFunction") ;
48  photonEcalEnergyCorrFunction_.reset(EcalClusterFunctionFactory::get()->create(photonEnergyFunctionName, config));
49  //ingredient for photon uncertainty
51 
52  if( config.existsAs<edm::ParameterSet>("regressionConfig") ) {
53  const edm::ParameterSet regr_conf =
54  config.getParameterSet("regressionConfig");
55  const std::string& mname = regr_conf.getParameter<std::string>("modifierName");
56  ModifyObjectValueBase* regr = ModifyObjectValueFactory::get()->create(mname,regr_conf);
57  gedRegression_.reset(regr);
58  } else {
59  gedRegression_.reset(nullptr);
60  }
61 
62  // ingredient for energy regression
63  weightsfromDB_= config.getParameter<bool>("regressionWeightsFromDB");
64  w_file_ = config.getParameter<std::string>("energyRegressionWeightsFileLocation");
65  if (weightsfromDB_) w_db_ = config.getParameter<std::string>("energyRegressionWeightsDBLocation");
66  else w_db_ = "none" ;
68 
69 
70 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
std::unique_ptr< EcalClusterFunctionBaseClass > scCrackEnergyFunction_
edm::EDGetTokenT< EcalRecHitCollection > endcapEcalHitsToken_
std::unique_ptr< EnergyUncertaintyPhotonSpecific > photonUncertaintyCalculator_
std::unique_ptr< EcalClusterFunctionBaseClass > scEnergyFunction_
edm::EDGetTokenT< EcalRecHitCollection > barrelEcalHitsToken_
std::unique_ptr< EcalClusterFunctionBaseClass > photonEcalEnergyCorrFunction_
ParameterSet const & getParameterSet(std::string const &) const
std::unique_ptr< EcalClusterFunctionBaseClass > scEnergyErrorFunction_
std::unique_ptr< ModifyObjectValueBase > gedRegression_
std::unique_ptr< EGEnergyCorrector > regressionCorrector_
SurfaceDeformation * create(int type, const std::vector< double > &params)
T get(const Candidate &c)
Definition: component.h:55

Member Function Documentation

double PhotonEnergyCorrector::applyCrackCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass crackCorrectionFunction 
)

Definition at line 185 of file PhotonEnergyCorrector.cc.

References reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), EcalClusterFunctionBaseClass::getValue(), and reco::SuperCluster::rawEnergy().

Referenced by calculate().

186  {
187 
188 
189  double crackcor = 1.;
190 
191  for(reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
192 
193  const reco::CaloClusterPtr cc = *cIt;
194  crackcor *= ( (cl.rawEnergy() +
195  cc->energy()*(crackCorrectionFunction->getValue(*cc)-1.)) /
196  cl.rawEnergy() );
197  }// loop on BCs
198 
199 
200  return crackcor;
201 
202 }
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:75
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:78
void PhotonEnergyCorrector::calculate ( edm::Event evt,
reco::Photon thePhoton,
int  subdet,
const reco::VertexCollection vtxcol,
const edm::EventSetup iSetup 
)

Definition at line 94 of file PhotonEnergyCorrector.cc.

References applyCrackCorrection(), barrelEcalHitsToken_, reco::Photon::e5x5(), DetId::Ecal, reco::Photon::ecal_photons, reco::Photon::ecal_standard, EcalBarrel, EcalEndcap, endcapEcalHitsToken_, edm::false, gedRegression_, reco::Photon::getCorrectedEnergy(), reco::Photon::getCorrectedEnergyError(), minR9Barrel_, minR9Endcap_, photonEcalEnergyCorrFunction_, photonUncertaintyCalculator_, reco::Photon::r9(), reco::Photon::regression1, reco::Photon::regression2, regressionCorrector_, scCrackEnergyFunction_, scEnergyFunction_, reco::Photon::setCorrectedEnergy(), reco::Photon::superCluster(), theCaloGeom_, w_file_, and weightsfromDB_.

Referenced by PhotonProducer::fillPhotonCollection(), and GEDPhotonProducer::fillPhotonCollection().

94  {
95 
96  double phoEcalEnergy = -9999.;
97  double phoEcalEnergyError = -9999.;
98  double phoRegr1Energy = -9999.;
99  double phoRegr1EnergyError = -9999.;
100  theCaloGeom_->getSubdetectorGeometry(DetId::Ecal, subdet);
101 
102  double minR9=0;
103  if (subdet==EcalBarrel) {
104  minR9=minR9Barrel_;
105  } else if (subdet==EcalEndcap) {
106  minR9=minR9Endcap_;
107  }
108 
110 
112  if ( thePhoton.r9() > minR9 ) {
113  // f(eta) correction to e5x5
114  double deltaE = scEnergyFunction_->getValue(*(thePhoton.superCluster()), 1);
115  float e5x5=thePhoton.e5x5();
116  if (subdet==EcalBarrel) e5x5 = e5x5 * (1.0 + deltaE/thePhoton.superCluster()->rawEnergy() );
117  phoEcalEnergy = e5x5 + thePhoton.superCluster()->preshowerEnergy() ;
118  } else {
119  phoEcalEnergy = thePhoton.superCluster()->energy();
120  }
121  // store the value in the Photon.h
122  thePhoton.setCorrectedEnergy( reco::Photon::ecal_standard, phoEcalEnergy, phoEcalEnergyError, false);
123 
125 
126  if ( thePhoton.r9() > minR9 ) {
127 
128 
129 
130  // f(eta) correction to e5x5
131  double deltaE = scEnergyFunction_->getValue(*(thePhoton.superCluster()), 1);
132  float e5x5=thePhoton.e5x5();
133  if (subdet==EcalBarrel) e5x5 = e5x5 * (1.0 + deltaE/thePhoton.superCluster()->rawEnergy() );
134  phoEcalEnergy = e5x5 + thePhoton.superCluster()->preshowerEnergy() ;
135  // add correction for cracks
136  phoEcalEnergy *= scCrackEnergyFunction_->getValue(*(thePhoton.superCluster()));
137  phoEcalEnergyError = photonUncertaintyCalculator_->computePhotonEnergyUncertainty_highR9(thePhoton.superCluster()->eta(), thePhoton.superCluster()->phiWidth()/thePhoton.superCluster()->etaWidth(), phoEcalEnergy);
138  } else {
139 
140 
141  // correction for low r9
142  phoEcalEnergy = photonEcalEnergyCorrFunction_->getValue(*(thePhoton.superCluster()), 1);
143  phoEcalEnergy *= applyCrackCorrection(*(thePhoton.superCluster()), scCrackEnergyFunction_.get());
144  phoEcalEnergyError = photonUncertaintyCalculator_->computePhotonEnergyUncertainty_lowR9(thePhoton.superCluster()->eta(), thePhoton.superCluster()->phiWidth()/thePhoton.superCluster()->etaWidth(), phoEcalEnergy);
145  }
146 
147 
148  // store the value in the Photon.h
149  thePhoton.setCorrectedEnergy( reco::Photon::ecal_photons, phoEcalEnergy, phoEcalEnergyError, false);
150 
152  //
153  if ( ( weightsfromDB_ && !gedRegression_) || ( !weightsfromDB_ && !(w_file_ == "none") ) ) {
154  std::pair<double,double> cor = regressionCorrector_->CorrectedEnergyWithError(thePhoton, vtxcol, lazyTools, iSetup);
155  phoRegr1Energy = cor.first;
156  phoRegr1EnergyError = cor.second;
157  // store the value in the Photon.h
158  thePhoton.setCorrectedEnergy( reco::Photon::regression1, phoRegr1Energy, phoRegr1EnergyError, false);
159  }
160 
161  if( gedRegression_ ) {
162  gedRegression_->modifyObject(thePhoton); // uses regression2 slot
163  // force regresions1 and 2 to be the same (no reason to be different)
167  false );
168  }
169 
170  /*
171  std::cout << " ------------------------- " << std::endl;
172  std::cout << " Corrector " << std::endl;
173  std::cout << " P4 Type " << thePhoton.getCandidateP4type() << " candidate p4 " << thePhoton.p4() << std::endl;
174  std::cout << " photon ecalEnergy " << thePhoton.getCorrectedEnergy(reco::Photon::ecal_photons) << " error " << thePhoton.getCorrectedEnergyError(reco::Photon::ecal_photons) << std::endl;
175  std::cout << " ecal p4 from accessor " << thePhoton.p4(reco::Photon::ecal_photons) << std::endl;
176  std::cout << " ------------------------- " << std::endl;
177  std::cout << " reg1 energy " << thePhoton.getCorrectedEnergy(reco::Photon::regression1) << " error " << thePhoton.getCorrectedEnergyError(reco::Photon::regression1) << std::endl;
178  std::cout << " New p4 from regression " << thePhoton.p4(reco::Photon::regression1) << std::endl;
179  std::cout << " ------------------------- " << std::endl;
180  */
181 
182 
183 }
edm::ESHandle< CaloGeometry > theCaloGeom_
reco::SuperClusterRef superCluster() const
Ref to SuperCluster.
std::unique_ptr< EcalClusterFunctionBaseClass > scCrackEnergyFunction_
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
float e5x5() const
Definition: Photon.h:221
double applyCrackCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
edm::EDGetTokenT< EcalRecHitCollection > endcapEcalHitsToken_
std::unique_ptr< EnergyUncertaintyPhotonSpecific > photonUncertaintyCalculator_
std::unique_ptr< EcalClusterFunctionBaseClass > scEnergyFunction_
edm::EDGetTokenT< EcalRecHitCollection > barrelEcalHitsToken_
float getCorrectedEnergyError(P4type type) const
std::unique_ptr< EcalClusterFunctionBaseClass > photonEcalEnergyCorrFunction_
float getCorrectedEnergy(P4type type) const
std::unique_ptr< ModifyObjectValueBase > gedRegression_
float r9() const
Definition: Photon.h:227
std::unique_ptr< EGEnergyCorrector > regressionCorrector_
volatile std::atomic< bool > shutdown_flag false
std::unique_ptr<ModifyObjectValueBase>& PhotonEnergyCorrector::gedRegression ( )
inline

Definition at line 33 of file PhotonEnergyCorrector.h.

References gedRegression_.

Referenced by GEDPhotonProducer::produce().

34  { return gedRegression_; }
std::unique_ptr< ModifyObjectValueBase > gedRegression_
void PhotonEnergyCorrector::init ( const edm::EventSetup theEventSetup)

Definition at line 72 of file PhotonEnergyCorrector.cc.

References edm::EventSetup::get(), photonEcalEnergyCorrFunction_, photonUncertaintyCalculator_, regressionCorrector_, scCrackEnergyFunction_, scEnergyErrorFunction_, scEnergyFunction_, theCaloGeom_, w_db_, w_file_, and weightsfromDB_.

Referenced by GEDPhotonProducer::produce().

72  {
73  theEventSetup.get<CaloGeometryRecord>().get(theCaloGeom_);
74 
75 
76  scEnergyFunction_->init(theEventSetup);
77  scCrackEnergyFunction_->init(theEventSetup);
78  scEnergyErrorFunction_->init(theEventSetup);
79  photonEcalEnergyCorrFunction_->init(theEventSetup);
80 
81  if ( weightsfromDB_ ) {
82  if (!regressionCorrector_->IsInitialized()) regressionCorrector_->Initialize(theEventSetup,w_db_,weightsfromDB_);
83  }
84  if ( !weightsfromDB_ && !(w_file_ == "none") ) {
85  if (!regressionCorrector_->IsInitialized()) regressionCorrector_->Initialize(theEventSetup,w_file_,weightsfromDB_);
86  }
87 
88 
89  photonUncertaintyCalculator_->init(theEventSetup);
90 
91 }
edm::ESHandle< CaloGeometry > theCaloGeom_
std::unique_ptr< EcalClusterFunctionBaseClass > scCrackEnergyFunction_
std::unique_ptr< EnergyUncertaintyPhotonSpecific > photonUncertaintyCalculator_
std::unique_ptr< EcalClusterFunctionBaseClass > scEnergyFunction_
std::unique_ptr< EcalClusterFunctionBaseClass > photonEcalEnergyCorrFunction_
const T & get() const
Definition: EventSetup.h:56
std::unique_ptr< EcalClusterFunctionBaseClass > scEnergyErrorFunction_
std::unique_ptr< EGEnergyCorrector > regressionCorrector_

Member Data Documentation

edm::InputTag PhotonEnergyCorrector::barrelEcalHits_
private

Definition at line 55 of file PhotonEnergyCorrector.h.

Referenced by PhotonEnergyCorrector().

edm::EDGetTokenT<EcalRecHitCollection> PhotonEnergyCorrector::barrelEcalHitsToken_
private

Definition at line 57 of file PhotonEnergyCorrector.h.

Referenced by calculate(), and PhotonEnergyCorrector().

std::string PhotonEnergyCorrector::candidateP4type_
private

Definition at line 45 of file PhotonEnergyCorrector.h.

edm::InputTag PhotonEnergyCorrector::endcapEcalHits_
private

Definition at line 56 of file PhotonEnergyCorrector.h.

Referenced by PhotonEnergyCorrector().

edm::EDGetTokenT<EcalRecHitCollection> PhotonEnergyCorrector::endcapEcalHitsToken_
private

Definition at line 58 of file PhotonEnergyCorrector.h.

Referenced by calculate(), and PhotonEnergyCorrector().

std::unique_ptr<ModifyObjectValueBase> PhotonEnergyCorrector::gedRegression_
private

Definition at line 51 of file PhotonEnergyCorrector.h.

Referenced by calculate(), gedRegression(), and PhotonEnergyCorrector().

double PhotonEnergyCorrector::minR9Barrel_
private

Definition at line 52 of file PhotonEnergyCorrector.h.

Referenced by calculate(), and PhotonEnergyCorrector().

double PhotonEnergyCorrector::minR9Endcap_
private

Definition at line 53 of file PhotonEnergyCorrector.h.

Referenced by calculate(), and PhotonEnergyCorrector().

std::unique_ptr<EcalClusterFunctionBaseClass> PhotonEnergyCorrector::photonEcalEnergyCorrFunction_
private

Definition at line 50 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

std::unique_ptr<EnergyUncertaintyPhotonSpecific> PhotonEnergyCorrector::photonUncertaintyCalculator_
private

Definition at line 60 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

std::unique_ptr<EGEnergyCorrector> PhotonEnergyCorrector::regressionCorrector_
private

Definition at line 46 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

std::unique_ptr<EcalClusterFunctionBaseClass> PhotonEnergyCorrector::scCrackEnergyFunction_
private

Definition at line 48 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

std::unique_ptr<EcalClusterFunctionBaseClass> PhotonEnergyCorrector::scEnergyErrorFunction_
private

Definition at line 49 of file PhotonEnergyCorrector.h.

Referenced by init(), and PhotonEnergyCorrector().

std::unique_ptr<EcalClusterFunctionBaseClass> PhotonEnergyCorrector::scEnergyFunction_
private

Definition at line 47 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

edm::ESHandle<CaloGeometry> PhotonEnergyCorrector::theCaloGeom_
private

Definition at line 54 of file PhotonEnergyCorrector.h.

Referenced by calculate(), and init().

std::string PhotonEnergyCorrector::w_db_
private

Definition at line 44 of file PhotonEnergyCorrector.h.

Referenced by init(), and PhotonEnergyCorrector().

std::string PhotonEnergyCorrector::w_file_
private

Definition at line 43 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().

bool PhotonEnergyCorrector::weightsfromDB_
private

Definition at line 42 of file PhotonEnergyCorrector.h.

Referenced by calculate(), init(), and PhotonEnergyCorrector().