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, edm::ConsumesCollector &&iC)
 
 ~PhotonEnergyCorrector ()
 

Private Attributes

edm::InputTag barrelEcalHits_
 
edm::EDGetTokenT
< EcalRecHitCollection
barrelEcalHitsToken_
 
std::string candidateP4type_
 
edm::InputTag endcapEcalHits_
 
edm::EDGetTokenT
< EcalRecHitCollection
endcapEcalHitsToken_
 
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,
edm::ConsumesCollector &&  iC 
)

Definition at line 13 of file PhotonEnergyCorrector.cc.

References barrelEcalHits_, barrelEcalHitsToken_, 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_ = 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_ = 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_ = 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_ = 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  gedRegression_.reset(new PFSCRegressionCalc(regr_conf));
56  }
57 
58  // ingredient for energy regression
59  weightsfromDB_= config.getParameter<bool>("regressionWeightsFromDB");
60  w_file_ = config.getParameter<std::string>("energyRegressionWeightsFileLocation");
61  if (weightsfromDB_) w_db_ = config.getParameter<std::string>("energyRegressionWeightsDBLocation");
62  else w_db_ = "none" ;
64 
65 
66 }
std::unique_ptr< PFSCRegressionCalc > gedRegression_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
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:186
edm::EDGetTokenT< EcalRecHitCollection > endcapEcalHitsToken_
EcalClusterFunctionBaseClass * photonEcalEnergyCorrFunction_
EGEnergyCorrector * regressionCorrector_
edm::EDGetTokenT< EcalRecHitCollection > barrelEcalHitsToken_
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 69 of file PhotonEnergyCorrector.cc.

References photonUncertaintyCalculator_, and regressionCorrector_.

69  {
70  delete regressionCorrector_;
72 }
EGEnergyCorrector * regressionCorrector_
EnergyUncertaintyPhotonSpecific * photonUncertaintyCalculator_

Member Function Documentation

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

Definition at line 193 of file PhotonEnergyCorrector.cc.

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

Referenced by calculate().

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

References applyCrackCorrection(), barrelEcalHitsToken_, EnergyUncertaintyPhotonSpecific::computePhotonEnergyUncertainty_highR9(), EnergyUncertaintyPhotonSpecific::computePhotonEnergyUncertainty_lowR9(), EGEnergyCorrector::CorrectedEnergyWithError(), reco::Photon::e5x5(), DetId::Ecal, reco::Photon::ecal_photons, reco::Photon::ecal_standard, EcalBarrel, EcalEndcap, endcapEcalHitsToken_, 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().

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

76  {
77  theEventSetup.get<CaloGeometryRecord>().get(theCaloGeom_);
78 
79 
80  scEnergyFunction_->init(theEventSetup);
81  scCrackEnergyFunction_->init(theEventSetup);
82  scEnergyErrorFunction_->init(theEventSetup);
83  photonEcalEnergyCorrFunction_->init(theEventSetup);
84 
85  if ( weightsfromDB_ ) {
87  }
88  if ( !weightsfromDB_ && !(w_file_ == "none") ) {
90  }
91 
92 
93  photonUncertaintyCalculator_->init(theEventSetup);
94 
95  if( gedRegression_ )
96  gedRegression_->update(theEventSetup);
97 
98 
99 }
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 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<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().