CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
EgammaSCEnergyCorrectionAlgo Class Reference

#include <EgammaSCEnergyCorrectionAlgo.h>

Public Member Functions

reco::SuperCluster applyCorrection (const reco::SuperCluster &cl, const EcalRecHitCollection &rhc, reco::CaloCluster::AlgoId theAlgo, const CaloSubdetectorGeometry *geometry, EcalClusterFunctionBaseClass *energyCorrectionFunction, std::string energyCorrectorName_, int modeEB_, int modeEE_)
 
reco::SuperCluster applyCrackCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
 
reco::SuperCluster applyLocalContCorrection (const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *localContCorrectionFunction)
 
 EgammaSCEnergyCorrectionAlgo (float noise, reco::CaloCluster::AlgoId theAlgo, const edm::ParameterSet &pset)
 
 ~EgammaSCEnergyCorrectionAlgo ()
 

Private Member Functions

float fNCrystals (int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
 
int nCrystalsGT2Sigma (reco::BasicCluster const &seed, EcalRecHitCollection const &rhc) const
 

Private Attributes

float sigmaElectronicNoise_
 
reco::CaloCluster::AlgoId theAlgo_
 

Detailed Description

Definition at line 20 of file EgammaSCEnergyCorrectionAlgo.h.

Constructor & Destructor Documentation

EgammaSCEnergyCorrectionAlgo::EgammaSCEnergyCorrectionAlgo ( float  noise,
reco::CaloCluster::AlgoId  theAlgo,
const edm::ParameterSet pset 
)

Definition at line 11 of file EgammaSCEnergyCorrectionAlgo.cc.

References sigmaElectronicNoise_.

EgammaSCEnergyCorrectionAlgo::~EgammaSCEnergyCorrectionAlgo ( )
inline

Member Function Documentation

reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCorrection ( const reco::SuperCluster cl,
const EcalRecHitCollection rhc,
reco::CaloCluster::AlgoId  theAlgo,
const CaloSubdetectorGeometry geometry,
EcalClusterFunctionBaseClass energyCorrectionFunction,
std::string  energyCorrectorName_,
int  modeEB_,
int  modeEE_ 
)

Definition at line 20 of file EgammaSCEnergyCorrectionAlgo.cc.

References SuperClusterShapeAlgo::Calculate_Covariances(), GetRecoTauVFromDQM_MC_cff::cl, reco::CaloCluster::dynamicHybrid, reco::CaloCluster::energy(), SuperClusterShapeAlgo::etaWidth(), fNCrystals(), EcalClusterFunctionBaseClass::getValue(), reco::CaloCluster::hybrid, LogTrace, reco::CaloCluster::multi5x5, nCrystalsGT2Sigma(), SuperClusterShapeAlgo::phiWidth(), reco::SuperCluster::preshowerEnergy(), reco::SuperCluster::rawEnergy(), reco::SuperCluster::seed(), reco::CaloCluster::setEnergy(), reco::SuperCluster::setEtaWidth(), reco::SuperCluster::setPhiWidth(), and tmp.

Referenced by ~EgammaSCEnergyCorrectionAlgo().

26  {
27 
28 
29  // A little bit of trivial info to be sure all is well
30 
31 
32  {
33  LogTrace("EgammaSCEnergyCorrectionAlgo")<< "::applyCorrection" << std::endl;
34  LogTrace()<< " SC has energy " << cl.energy() << std::endl;
35  LogTrace()<< " Will correct now.... " << std::endl;
36  }
37 
38  // Get the seed cluster
39  const reco::CaloClusterPtr& seedC = cl.seed();
40 
41 
42  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " Seed cluster energy... " << seedC->energy() << std::endl;
43 
44 
45 
46  // Find the algorithm used to construct the basic clusters making up the supercluster
47  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " The seed cluster used algo " << theAlgo;
48 
49 
50  // Find the detector region of the supercluster
51  // where is the seed cluster?
52  std::vector<std::pair<DetId, float> > seedHits = seedC->hitsAndFractions();
53  EcalSubdetector theBase = EcalSubdetector(seedHits.at(0).first.subdetId());
54 
55  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " seed cluster location == " << theBase << std::endl;
56 
57  // Get number of crystals 2sigma above noise in seed basiccluster
58  int nCryGT2Sigma = nCrystalsGT2Sigma(*seedC,rhc);
59 
60  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " nCryGT2Sigma " << nCryGT2Sigma << std::endl;
61 
62 
63  // Supercluster enery - seed basiccluster energy
64  float bremsEnergy = cl.energy() - seedC->energy();
65 
66  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " bremsEnergy " << bremsEnergy << std::endl;
67 
68 
69  //Create the pointer ot class SuperClusterShapeAlgo
70  //which calculates phiWidth and etaWidth
71  SuperClusterShapeAlgo SCShape(&rhc, geometry);
72 
73  double phiWidth = 0.;
74  double etaWidth = 0.;
75  //Calculate phiWidth & etaWidth for SuperClusters
76  SCShape.Calculate_Covariances(cl);
77  phiWidth = SCShape.phiWidth();
78  etaWidth = SCShape.etaWidth();
79 
80  // Calculate the new supercluster energy
81  //as a function of number of crystals in the seed basiccluster for Endcap
82  //or apply new Enegry SCale correction
83  float newEnergy = 0;
84 
86  tmp.setPhiWidth(phiWidth);
87  tmp.setEtaWidth(etaWidth);
88 
89 
90  if ( theAlgo == reco::CaloCluster::hybrid || theAlgo == reco::CaloCluster::dynamicHybrid ) {
91  if (energyCorrectorName_=="EcalClusterEnergyCorrection") newEnergy = tmp.rawEnergy() + energyCorrectionFunction->getValue(tmp, modeEB_);
92  if (energyCorrectorName_=="EcalClusterEnergyCorrectionObjectSpecific") {
93  //std::cout << "newEnergy="<<newEnergy<<std::endl;
94  newEnergy = energyCorrectionFunction->getValue(tmp, modeEB_);
95  }
96 
97  } else if ( theAlgo == reco::CaloCluster::multi5x5 ) {
98  if (energyCorrectorName_=="EcalClusterEnergyCorrection") newEnergy = tmp.rawEnergy() + tmp.preshowerEnergy() + energyCorrectionFunction->getValue(tmp, modeEE_);
99  if (energyCorrectorName_=="EcalClusterEnergyCorrectionObjectSpecific") newEnergy = energyCorrectionFunction->getValue(tmp, modeEE_);
100 
101  } else {
102  //Apply f(nCry) correction on island algo and fixedMatrix algo
103  newEnergy = seedC->energy()/fNCrystals(nCryGT2Sigma, theAlgo, theBase)+bremsEnergy;
104  }
105 
106  // Create a new supercluster with the corrected energy
107 
108  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " UNCORRECTED SC has energy... " << cl.energy() << std::endl;
109  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " CORRECTED SC has energy... " << newEnergy << std::endl;
110 
111  reco::SuperCluster corrCl =cl;
112 
113  corrCl.setEnergy(newEnergy);
114  corrCl.setPhiWidth(phiWidth);
115  corrCl.setEtaWidth(etaWidth);
116 
117  return corrCl;
118 }
float fNCrystals(int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
void setEnergy(double energy)
Definition: CaloCluster.h:113
void setPhiWidth(double pw)
Definition: SuperCluster.h:62
void setEtaWidth(double ew)
Definition: SuperCluster.h:63
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:126
#define LogTrace(id)
int nCrystalsGT2Sigma(reco::BasicCluster const &seed, EcalRecHitCollection const &rhc) const
double rawEnergy() const
raw uncorrected energy (sum of energies of component BasicClusters)
Definition: SuperCluster.h:47
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:50
EcalSubdetector
reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyCrackCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass crackCorrectionFunction 
)

Definition at line 231 of file EgammaSCEnergyCorrectionAlgo.cc.

References applyLocalContCorrection(), GetRecoTauVFromDQM_MC_cff::cl, reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), reco::CaloCluster::energy(), EcalClusterFunctionBaseClass::getValue(), reco::SuperCluster::rawEnergy(), and reco::CaloCluster::setEnergy().

Referenced by ~EgammaSCEnergyCorrectionAlgo().

232  {
233 
234 
235  double crackcor = 1.;
236 
237  for(reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
238 
239  const reco::CaloClusterPtr cc = *cIt;
240  crackcor *= ( (cl.rawEnergy() +
241  cc->energy()*(crackCorrectionFunction->getValue(*cc)-1.)) /
242  cl.rawEnergy() );
243  }// loop on BCs
244 
245 
246  reco::SuperCluster corrCl=cl;
247  corrCl.setEnergy(cl.energy()*crackcor);
248 
249 
250  return corrCl;
251 }
void setEnergy(double energy)
Definition: CaloCluster.h:113
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:126
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
reco::SuperCluster EgammaSCEnergyCorrectionAlgo::applyLocalContCorrection ( const reco::SuperCluster cl,
EcalClusterFunctionBaseClass localContCorrectionFunction 
)

Definition at line 259 of file EgammaSCEnergyCorrectionAlgo.cc.

References GetRecoTauVFromDQM_MC_cff::cl, reco::CaloCluster::energy(), EcalClusterFunctionBaseClass::getValue(), reco::SuperCluster::seed(), and reco::CaloCluster::setEnergy().

Referenced by applyCrackCorrection(), and ~EgammaSCEnergyCorrectionAlgo().

260  {
261 
262 
264 
265  const reco::CaloClusterPtr & seedBC = cl.seed();
266  float seedBCene = seedBC->energy();
267  float correctedSeedBCene = localContCorrectionFunction->getValue(*seedBC,dummy) * seedBCene;
268 
269 
270  reco::SuperCluster correctedSC = cl;
271  correctedSC.setEnergy(cl.energy() - seedBCene + correctedSeedBCene);
272 
273  return correctedSC;
274 
275 
276 
277 }
void setEnergy(double energy)
Definition: CaloCluster.h:113
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
double energy() const
cluster energy
Definition: CaloCluster.h:126
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:66
float EgammaSCEnergyCorrectionAlgo::fNCrystals ( int  nCry,
reco::CaloCluster::AlgoId  theAlgo,
EcalSubdetector  theBase 
) const
private

Definition at line 120 of file EgammaSCEnergyCorrectionAlgo.cc.

References EcalBarrel, EcalEndcap, f, reco::CaloCluster::hybrid, reco::CaloCluster::island, LogTrace, p1, p2, p3, p4, mps_fire::result, and x.

Referenced by applyCorrection(), and ~EgammaSCEnergyCorrectionAlgo().

120  {
121 
122  float p0 = 0, p1 = 0, p2 = 0, p3 = 0, p4 = 0;
123  float x = nCry;
124  float result =1.f;
125 
126  if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::hybrid)) {
127  if (nCry<=10)
128  {
129  p0 = 6.32879e-01f;
130  p1 = 1.14893e-01f;
131  p2 = -2.45705e-02f;
132  p3 = 2.53074e-03f;
133  p4 = -9.29654e-05f;
134  }
135  else if (nCry>10 && nCry<=30)
136  {
137  p0 = 6.93196e-01f;
138  p1 = 4.44034e-02f;
139  p2 = -2.82229e-03f;
140  p3 = 8.19495e-05f;
141  p4 = -8.96645e-07f;
142  }
143  else
144  {
145  p0 = 5.65474e+00f;
146  p1 = -6.31640e-01f;
147  p2 = 3.14218e-02f;
148  p3 = -6.84256e-04f;
149  p4 = 5.50659e-06f;
150  }
151  if (x > 40.f) x = 40.f;
152  }
153 
154  else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::hybrid)) {
155 
156  LogTrace("EgammaSCEnergyCorrectionAlgo") << "ERROR! HybridEFRYsc called" << std::endl;
157 
158  return 1.f;
159  }
160 
161  else if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
162  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
163  p1 = 1.45900e-01f;
164  p2 = -1.61359e-02f;
165  p3 = 7.99423e-04f;
166  p4 = -1.47873e-05f;
167  if (x > 16.f) x = 16.f;
168  }
169 
170  else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
171  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
172  p1 = 1.45900e-01f;
173  p2 = -1.61359e-02f;
174  p3 = 7.99423e-04f;
175  p4 = -1.47873e-05f;
176  if (x > 16.f) x = 16.f;
177  }
178 
179  else {
180 
181  LogTrace("EgammaSCEnergyCorrectionAlgo")<< "trying to correct unknown cluster!!!" << std::endl;
182  return 1.f;
183  }
184  result = p0 + x*(p1 + x*(p2 + x*(p3 + x*p4)));
185 
186  //Rescale energy scale correction to take into account change in calibrated
187  //RecHit definition introduced in CMSSW_1_5_0
188  float const ebfact = 1.f/0.965f;
189  float const eefact = 1.f/0.975f;
190 
191  if(theBase == EcalBarrel) {
192  result*=ebfact;
193  } else {
194  result*=eefact;
195  }
196 
197  return result;
198 }
double p4[4]
Definition: TauolaWrapper.h:92
double f[11][100]
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
double p1[4]
Definition: TauolaWrapper.h:89
double p3[4]
Definition: TauolaWrapper.h:91
int EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma ( reco::BasicCluster const &  seed,
EcalRecHitCollection const &  rhc 
) const
private

Definition at line 200 of file EgammaSCEnergyCorrectionAlgo.cc.

References edm::SortedCollection< T, SORT >::find(), hfClusterShapes_cfi::hits, LogTrace, and sigmaElectronicNoise_.

Referenced by applyCorrection(), and ~EgammaSCEnergyCorrectionAlgo().

201 {
202  // return number of crystals 2Sigma above
203  // electronic noise
204 
205  std::vector<std::pair<DetId,float > > const & hits = seed.hitsAndFractions();
206 
207 
208  LogTrace("EgammaSCEnergyCorrectionAlgo") << " EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma" << std::endl;
209  LogTrace("EgammaSCEnergyCorrectionAlgo") << " Will calculate number of crystals above 2sigma noise" << std::endl;
210  LogTrace("EgammaSCEnergyCorrectionAlgo") << " sigmaElectronicNoise = " << sigmaElectronicNoise_ << std::endl;
211  LogTrace("EgammaSCEnergyCorrectionAlgo") << " There are " << hits.size() << " recHits" << std::endl;
212 
213  int nCry = 0;
214  std::vector<std::pair<DetId,float > >::const_iterator hit;
216  for(hit = hits.begin(); hit != hits.end(); hit++) {
217  // need to get hit by DetID in order to get energy
218  aHit = rhc.find((*hit).first);
219  if((*aHit).energy()>2.f*sigmaElectronicNoise_) nCry++;
220  }
221 
222  LogTrace("EgammaSCEnergyCorrectionAlgo") << " " << nCry << " of these above 2sigma noise" << std::endl;
223 
224 
225  return nCry;
226 }
std::vector< EcalRecHit >::const_iterator const_iterator
#define LogTrace(id)

Member Data Documentation

float EgammaSCEnergyCorrectionAlgo::sigmaElectronicNoise_
private
reco::CaloCluster::AlgoId EgammaSCEnergyCorrectionAlgo::theAlgo_
private

Definition at line 60 of file EgammaSCEnergyCorrectionAlgo.h.