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
< PFSCRegressionCalc > & 
gedRegression ()
 
void init (const edm::EventSetup &theEventSetup)
 
 PhotonEnergyCorrector (const edm::ParameterSet &config)
 
 ~PhotonEnergyCorrector ()
 

Private Attributes

edm::InputTag barrelEcalHits_
 
std::string candidateP4type_
 
edm::InputTag endcapEcalHits_
 
std::unique_ptr
< PFSCRegressionCalc
gedRegression_
 
double minR9Barrel_
 
double minR9Endcap_
 
EcalClusterFunctionBaseClassphotonEcalEnergyCorrFunction_
 
EnergyUncertaintyPhotonSpecificphotonUncertaintyCalculator_
 
EGEnergyCorrectorregressionCorrector_
 
EcalClusterFunctionBaseClassscCrackEnergyFunction_
 
EcalClusterFunctionBaseClassscEnergyErrorFunction_
 
EcalClusterFunctionBaseClassscEnergyFunction_
 
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 26 of file PhotonEnergyCorrector.h.

Constructor & Destructor Documentation

PhotonEnergyCorrector::PhotonEnergyCorrector ( const edm::ParameterSet config)

Definition at line 13 of file PhotonEnergyCorrector.cc.

References barrelEcalHits_, endcapEcalHits_, 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 
16  minR9Barrel_ = config.getParameter<double>("minR9Barrel");
17  minR9Endcap_ = config.getParameter<double>("minR9Endcap");
18  // get the geometry from the event setup:
19 
20  barrelEcalHits_ = config.getParameter<edm::InputTag>("barrelEcalHits");
21  endcapEcalHits_ = config.getParameter<edm::InputTag>("endcapEcalHits");
22  // candidateP4type_ = config.getParameter<std::string>("candidateP4type") ;
23 
24 
25  // function to extract f(eta) correction
26  scEnergyFunction_ = 0 ;
27  std::string superClusterFunctionName = config.getParameter<std::string>("superClusterEnergyCorrFunction") ;
28  scEnergyFunction_ = EcalClusterFunctionFactory::get()->create(superClusterFunctionName,config) ;
29 
30 
31  // function to extract corrections to cracks
33  std::string superClusterCrackFunctionName = config.getParameter<std::string>("superClusterCrackEnergyCorrFunction") ;
34  scCrackEnergyFunction_ = EcalClusterFunctionFactory::get()->create(superClusterCrackFunctionName,config) ;
35 
36 
37  // function to extract the error on the sc ecal correction
39  std::string superClusterErrorFunctionName = config.getParameter<std::string>("superClusterEnergyErrorFunction") ;
40  scEnergyErrorFunction_ = EcalClusterFunctionFactory::get()->create(superClusterErrorFunctionName,config) ;
41 
42 
43  // function to extract the error on the photon ecal correction
45  std::string photonEnergyFunctionName = config.getParameter<std::string>("photonEcalEnergyCorrFunction") ;
46  photonEcalEnergyCorrFunction_ = EcalClusterFunctionFactory::get()->create(photonEnergyFunctionName, config);
47  //ingredient for photon uncertainty
49 
50  if( config.existsAs<edm::ParameterSet>("regressionConfig") ) {
51  const edm::ParameterSet regr_conf =
52  config.getParameterSet("regressionConfig");
53  gedRegression_.reset(new PFSCRegressionCalc(regr_conf));
54  }
55 
56  // ingredient for energy regression
57  weightsfromDB_= config.getParameter<bool>("regressionWeightsFromDB");
58  w_file_ = config.getParameter<std::string>("energyRegressionWeightsFileLocation");
59  if (weightsfromDB_) w_db_ = config.getParameter<std::string>("energyRegressionWeightsDBLocation");
60  else w_db_ == "none" ;
62 
63 
64 }
std::unique_ptr< PFSCRegressionCalc > gedRegression_
T getParameter(std::string const &) const
EcalClusterFunctionBaseClass * scEnergyFunction_
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:184
EcalClusterFunctionBaseClass * photonEcalEnergyCorrFunction_
EGEnergyCorrector * regressionCorrector_
EcalClusterFunctionBaseClass * scEnergyErrorFunction_
EnergyUncertaintyPhotonSpecific * photonUncertaintyCalculator_
ParameterSet const & getParameterSet(std::string const &) const
EcalClusterFunctionBaseClass * scCrackEnergyFunction_
T get(const Candidate &c)
Definition: component.h:55
SCRegressionCalculator< BaselinePFSCRegression > PFSCRegressionCalc
PhotonEnergyCorrector::~PhotonEnergyCorrector ( )

Definition at line 67 of file PhotonEnergyCorrector.cc.

References photonUncertaintyCalculator_, and regressionCorrector_.

67  {
68  delete regressionCorrector_;
70 }
EGEnergyCorrector * regressionCorrector_
EnergyUncertaintyPhotonSpecific * photonUncertaintyCalculator_

Member Function Documentation

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

Definition at line 195 of file PhotonEnergyCorrector.cc.

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

Referenced by calculate().

196  {
197 
198 
199  double crackcor = 1.;
200 
201  for(reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
202 
203  const reco::CaloClusterPtr cc = *cIt;
204  crackcor *= ( (cl.rawEnergy() +
205  cc->energy()*(crackCorrectionFunction->getValue(*cc)-1.)) /
206  cl.rawEnergy() );
207  }// loop on BCs
208 
209 
210  return crackcor;
211 
212 }
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 100 of file PhotonEnergyCorrector.cc.

References applyCrackCorrection(), barrelEcalHits_, EnergyUncertaintyPhotonSpecific::computePhotonEnergyUncertainty_highR9(), EnergyUncertaintyPhotonSpecific::computePhotonEnergyUncertainty_lowR9(), EGEnergyCorrector::CorrectedEnergyWithError(), reco::Photon::e5x5(), DetId::Ecal, reco::Photon::ecal_photons, reco::Photon::ecal_standard, EcalBarrel, EcalEndcap, endcapEcalHits_, gedRegression_, EcalClusterFunctionBaseClass::getValue(), minR9Barrel_, minR9Endcap_, photonEcalEnergyCorrFunction_, photonUncertaintyCalculator_, reco::Photon::r9(), reco::Photon::regression1, regressionCorrector_, scCrackEnergyFunction_, scEnergyFunction_, reco::Photon::setCorrectedEnergy(), reco::Photon::superCluster(), theCaloGeom_, w_file_, and weightsfromDB_.

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

100  {
101 
102  double phoEcalEnergy = -9999.;
103  double phoEcalEnergyError = -9999.;
104  double phoRegr1Energy = -9999.;
105  double phoRegr1EnergyError = -9999.;
106  theCaloGeom_->getSubdetectorGeometry(DetId::Ecal, subdet);
107 
108  double minR9=0;
109  if (subdet==EcalBarrel) {
110  minR9=minR9Barrel_;
111  } else if (subdet==EcalEndcap) {
112  minR9=minR9Endcap_;
113  }
114 
115 
116 
117  EcalClusterLazyTools lazyTools(evt, iSetup, barrelEcalHits_,endcapEcalHits_);
118 
119 
120 
122  if ( thePhoton.r9() > minR9 ) {
123  // f(eta) correction to e5x5
124  double deltaE = scEnergyFunction_->getValue(*(thePhoton.superCluster()), 1);
125  float e5x5=thePhoton.e5x5();
126  if (subdet==EcalBarrel) e5x5 = e5x5 * (1.0 + deltaE/thePhoton.superCluster()->rawEnergy() );
127  phoEcalEnergy = e5x5 + thePhoton.superCluster()->preshowerEnergy() ;
128  } else {
129  phoEcalEnergy = thePhoton.superCluster()->energy();
130  }
131  // store the value in the Photon.h
132  thePhoton.setCorrectedEnergy( reco::Photon::ecal_standard, phoEcalEnergy, phoEcalEnergyError, false);
133 
135 
136  if ( thePhoton.r9() > minR9 ) {
137 
138 
139 
140  // f(eta) correction to e5x5
141  double deltaE = scEnergyFunction_->getValue(*(thePhoton.superCluster()), 1);
142  float e5x5=thePhoton.e5x5();
143  if (subdet==EcalBarrel) e5x5 = e5x5 * (1.0 + deltaE/thePhoton.superCluster()->rawEnergy() );
144  phoEcalEnergy = e5x5 + thePhoton.superCluster()->preshowerEnergy() ;
145  // add correction for cracks
146  phoEcalEnergy *= scCrackEnergyFunction_->getValue(*(thePhoton.superCluster()));
147  phoEcalEnergyError = photonUncertaintyCalculator_->computePhotonEnergyUncertainty_highR9(thePhoton.superCluster()->eta(), thePhoton.superCluster()->phiWidth()/thePhoton.superCluster()->etaWidth(), phoEcalEnergy);
148  } else {
149 
150 
151  // correction for low r9
152  phoEcalEnergy = photonEcalEnergyCorrFunction_->getValue(*(thePhoton.superCluster()), 1);
153  phoEcalEnergy *= applyCrackCorrection(*(thePhoton.superCluster()), scCrackEnergyFunction_);
154  phoEcalEnergyError = photonUncertaintyCalculator_->computePhotonEnergyUncertainty_lowR9(thePhoton.superCluster()->eta(), thePhoton.superCluster()->phiWidth()/thePhoton.superCluster()->etaWidth(), phoEcalEnergy);
155  }
156 
157 
158  // store the value in the Photon.h
159  thePhoton.setCorrectedEnergy( reco::Photon::ecal_photons, phoEcalEnergy, phoEcalEnergyError, false);
160 
162  //
163  if ( ( weightsfromDB_ && !gedRegression_) || ( !weightsfromDB_ && !(w_file_ == "none") ) ) {
164  std::pair<double,double> cor = regressionCorrector_->CorrectedEnergyWithError(thePhoton, vtxcol, lazyTools, iSetup);
165  phoRegr1Energy = cor.first;
166  phoRegr1EnergyError = cor.second;
167  // store the value in the Photon.h
168  thePhoton.setCorrectedEnergy( reco::Photon::regression1, phoRegr1Energy, phoRegr1EnergyError, false);
169  }
170 
171  if( gedRegression_ ) {
172  gedRegression_->varCalc()->setEvent(evt);
173  std::pair<float,float> cor = gedRegression_->getCorrectionWithErrors(*(thePhoton.superCluster()));
174  phoRegr1Energy = cor.first*thePhoton.superCluster()->correctedEnergy();
175  phoRegr1EnergyError = cor.second*thePhoton.superCluster()->correctedEnergy();
176  // store the value in the Photon.h
177  thePhoton.setCorrectedEnergy( reco::Photon::regression1, phoRegr1Energy, phoRegr1EnergyError, false);
178  }
179 
180  /*
181  std::cout << " ------------------------- " << std::endl;
182  std::cout << " Corrector " << std::endl;
183  std::cout << " P4 Type " << thePhoton.getCandidateP4type() << " candidate p4 " << thePhoton.p4() << std::endl;
184  std::cout << " photon ecalEnergy " << thePhoton.getCorrectedEnergy(reco::Photon::ecal_photons) << " error " << thePhoton.getCorrectedEnergyError(reco::Photon::ecal_photons) << std::endl;
185  std::cout << " ecal p4 from accessor " << thePhoton.p4(reco::Photon::ecal_photons) << std::endl;
186  std::cout << " ------------------------- " << std::endl;
187  std::cout << " reg1 energy " << thePhoton.getCorrectedEnergy(reco::Photon::regression1) << " error " << thePhoton.getCorrectedEnergyError(reco::Photon::regression1) << std::endl;
188  std::cout << " New p4 from regression " << thePhoton.p4(reco::Photon::regression1) << std::endl;
189  std::cout << " ------------------------- " << std::endl;
190  */
191 
192 
193 }
std::unique_ptr< PFSCRegressionCalc > gedRegression_
EcalClusterFunctionBaseClass * scEnergyFunction_
edm::ESHandle< CaloGeometry > theCaloGeom_
double computePhotonEnergyUncertainty_highR9(double eta, double brem, double energy)
float e5x5() const
Definition: Photon.h:185
double applyCrackCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
reco::SuperClusterRef superCluster() const
Ref to SuperCluster.
Definition: Photon.cc:58
EcalClusterFunctionBaseClass * photonEcalEnergyCorrFunction_
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
std::pair< double, double > CorrectedEnergyWithError(const reco::Photon &p, const reco::VertexCollection &vtxcol, EcalClusterLazyTools &clustertools, const edm::EventSetup &es)
EGEnergyCorrector * regressionCorrector_
void setCorrectedEnergy(P4type type, float E, float dE, bool toCand=true)
Definition: Photon.cc:92
double computePhotonEnergyUncertainty_lowR9(double eta, double brem, double energy)
EnergyUncertaintyPhotonSpecific * photonUncertaintyCalculator_
EcalClusterFunctionBaseClass * scCrackEnergyFunction_
float r9() const
Definition: Photon.h:191
std::unique_ptr<PFSCRegressionCalc>& PhotonEnergyCorrector::gedRegression ( )
inline

Definition at line 33 of file PhotonEnergyCorrector.h.

References gedRegression_.

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

Definition at line 74 of file PhotonEnergyCorrector.cc.

References gedRegression_, edm::EventSetup::get(), EnergyUncertaintyPhotonSpecific::init(), EcalClusterFunctionBaseClass::init(), EGEnergyCorrector::Initialize(), EGEnergyCorrector::IsInitialized(), photonEcalEnergyCorrFunction_, photonUncertaintyCalculator_, regressionCorrector_, scCrackEnergyFunction_, scEnergyErrorFunction_, scEnergyFunction_, theCaloGeom_, w_db_, w_file_, and weightsfromDB_.

74  {
75  theEventSetup.get<CaloGeometryRecord>().get(theCaloGeom_);
76 
77 
78  scEnergyFunction_->init(theEventSetup);
79  scCrackEnergyFunction_->init(theEventSetup);
80  scEnergyErrorFunction_->init(theEventSetup);
81  photonEcalEnergyCorrFunction_->init(theEventSetup);
82 
83  if ( weightsfromDB_ ) {
85  }
86  if ( !weightsfromDB_ && !(w_file_ == "none") ) {
88  }
89 
90 
91  photonUncertaintyCalculator_->init(theEventSetup);
92 
93  if( gedRegression_ )
94  gedRegression_->update(theEventSetup);
95 
96 
97 }
std::unique_ptr< PFSCRegressionCalc > gedRegression_
void Initialize(const edm::EventSetup &iSetup, std::string regweights, bool weightsFromDB=false)
EcalClusterFunctionBaseClass * scEnergyFunction_
edm::ESHandle< CaloGeometry > theCaloGeom_
EcalClusterFunctionBaseClass * photonEcalEnergyCorrFunction_
EGEnergyCorrector * regressionCorrector_
EcalClusterFunctionBaseClass * scEnergyErrorFunction_
EnergyUncertaintyPhotonSpecific * photonUncertaintyCalculator_
EcalClusterFunctionBaseClass * scCrackEnergyFunction_
const T & get() const
Definition: EventSetup.h:55
Bool_t IsInitialized() const
virtual void init(const edm::EventSetup &es)=0
void init(const edm::EventSetup &theEventSetup)

Member Data Documentation

edm::InputTag PhotonEnergyCorrector::barrelEcalHits_
private

Definition at line 55 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 calculate(), and PhotonEnergyCorrector().

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

Definition at line 51 of file PhotonEnergyCorrector.h.

Referenced by calculate(), gedRegression(), init(), 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().

EcalClusterFunctionBaseClass* PhotonEnergyCorrector::photonEcalEnergyCorrFunction_
private

Definition at line 50 of file PhotonEnergyCorrector.h.

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

EnergyUncertaintyPhotonSpecific* PhotonEnergyCorrector::photonUncertaintyCalculator_
private
EGEnergyCorrector* PhotonEnergyCorrector::regressionCorrector_
private
EcalClusterFunctionBaseClass* PhotonEnergyCorrector::scCrackEnergyFunction_
private

Definition at line 48 of file PhotonEnergyCorrector.h.

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

EcalClusterFunctionBaseClass* PhotonEnergyCorrector::scEnergyErrorFunction_
private

Definition at line 49 of file PhotonEnergyCorrector.h.

Referenced by init(), and PhotonEnergyCorrector().

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().