CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EgammaSCEnergyCorrectionAlgo.cc
Go to the documentation of this file.
1 //
2 // $Id: EgammaSCEnergyCorrectionAlgo.cc,v 1.48 2012/04/19 13:13:12 argiro Exp $
3 // Author: David Evans, Bristol
4 //
8 #include <iostream>
9 #include <string>
10 #include <vector>
11 
14  const edm::ParameterSet& pset)
15 
16 {
17  sigmaElectronicNoise_ = noise;
18 }
19 
20 
24  EcalClusterFunctionBaseClass* energyCorrectionFunction,
25  std::string energyCorrectorName_,
26  const int modeEB_,
27  const int modeEE_) {
28 
29 
30  // A little bit of trivial info to be sure all is well
31 
32 
33  {
34  LogTrace("EgammaSCEnergyCorrectionAlgo")<< "::applyCorrection" << std::endl;
35  LogTrace()<< " SC has energy " << cl.energy() << std::endl;
36  LogTrace()<< " Will correct now.... " << std::endl;
37  }
38 
39  // Get the seed cluster
40  reco::CaloClusterPtr seedC = cl.seed();
41 
42 
43  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " Seed cluster energy... " << seedC->energy() << std::endl;
44 
45 
46 
47  // Find the algorithm used to construct the basic clusters making up the supercluster
48  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " The seed cluster used algo " << theAlgo;
49 
50 
51  // Find the detector region of the supercluster
52  // where is the seed cluster?
53  std::vector<std::pair<DetId, float> > seedHits = seedC->hitsAndFractions();
54  EcalSubdetector theBase = EcalSubdetector(seedHits.at(0).first.subdetId());
55 
56  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " seed cluster location == " << theBase << std::endl;
57 
58  // Get number of crystals 2sigma above noise in seed basiccluster
59  int nCryGT2Sigma = nCrystalsGT2Sigma(*seedC,rhc);
60 
61  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " nCryGT2Sigma " << nCryGT2Sigma << std::endl;
62 
63 
64  // Supercluster enery - seed basiccluster energy
65  float bremsEnergy = cl.energy() - seedC->energy();
66 
67  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " bremsEnergy " << bremsEnergy << std::endl;
68 
69 
70  //Create the pointer ot class SuperClusterShapeAlgo
71  //which calculates phiWidth and etaWidth
72  SuperClusterShapeAlgo SCShape(&rhc, geometry);
73 
74  double phiWidth = 0.;
75  double etaWidth = 0.;
76  //Calculate phiWidth & etaWidth for SuperClusters
77  SCShape.Calculate_Covariances(cl);
78  phiWidth = SCShape.phiWidth();
79  etaWidth = SCShape.etaWidth();
80 
81  // Calculate the new supercluster energy
82  //as a function of number of crystals in the seed basiccluster for Endcap
83  //or apply new Enegry SCale correction
84  float newEnergy = 0;
85 
87  tmp.setPhiWidth(phiWidth);
88  tmp.setEtaWidth(etaWidth);
89 
90 
91  if ( theAlgo == reco::CaloCluster::hybrid || theAlgo == reco::CaloCluster::dynamicHybrid ) {
92  if (energyCorrectorName_=="EcalClusterEnergyCorrection") newEnergy = tmp.rawEnergy() + energyCorrectionFunction->getValue(tmp, modeEB_);
93  if (energyCorrectorName_=="EcalClusterEnergyCorrectionObjectSpecific") {
94  //std::cout << "newEnergy="<<newEnergy<<std::endl;
95  newEnergy = energyCorrectionFunction->getValue(tmp, modeEB_);
96  }
97 
98  } else if ( theAlgo == reco::CaloCluster::multi5x5 ) {
99  if (energyCorrectorName_=="EcalClusterEnergyCorrection") newEnergy = tmp.rawEnergy() + tmp.preshowerEnergy() + energyCorrectionFunction->getValue(tmp, modeEE_);
100  if (energyCorrectorName_=="EcalClusterEnergyCorrectionObjectSpecific") newEnergy = energyCorrectionFunction->getValue(tmp, modeEE_);
101 
102  } else {
103  //Apply f(nCry) correction on island algo and fixedMatrix algo
104  newEnergy = seedC->energy()/fNCrystals(nCryGT2Sigma, theAlgo, theBase)+bremsEnergy;
105  }
106 
107  // Create a new supercluster with the corrected energy
108 
109  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " UNCORRECTED SC has energy... " << cl.energy() << std::endl;
110  LogTrace("EgammaSCEnergyCorrectionAlgo")<< " CORRECTED SC has energy... " << newEnergy << std::endl;
111 
112  reco::SuperCluster corrCl =cl;
113 
114  corrCl.setEnergy(newEnergy);
115  corrCl.setPhiWidth(phiWidth);
116  corrCl.setEtaWidth(etaWidth);
117 
118  return corrCl;
119 }
120 
122 
123  float p0 = 0, p1 = 0, p2 = 0, p3 = 0, p4 = 0;
124  float x = nCry;
125  float result =1.f;
126 
127  if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::hybrid)) {
128  if (nCry<=10)
129  {
130  p0 = 6.32879e-01f;
131  p1 = 1.14893e-01f;
132  p2 = -2.45705e-02f;
133  p3 = 2.53074e-03f;
134  p4 = -9.29654e-05f;
135  }
136  else if (nCry>10 && nCry<=30)
137  {
138  p0 = 6.93196e-01f;
139  p1 = 4.44034e-02f;
140  p2 = -2.82229e-03f;
141  p3 = 8.19495e-05f;
142  p4 = -8.96645e-07f;
143  }
144  else
145  {
146  p0 = 5.65474e+00f;
147  p1 = -6.31640e-01f;
148  p2 = 3.14218e-02f;
149  p3 = -6.84256e-04f;
150  p4 = 5.50659e-06f;
151  }
152  if (x > 40.f) x = 40.f;
153  }
154 
155  else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::hybrid)) {
156 
157  LogTrace("EgammaSCEnergyCorrectionAlgo") << "ERROR! HybridEFRYsc called" << std::endl;
158 
159  return 1.f;
160  }
161 
162  else if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
163  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
164  p1 = 1.45900e-01f;
165  p2 = -1.61359e-02f;
166  p3 = 7.99423e-04f;
167  p4 = -1.47873e-05f;
168  if (x > 16.f) x = 16.f;
169  }
170 
171  else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
172  p0 = 4.69976e-01f; // extracted from fit to all endcap classes with Ebrem = 0.
173  p1 = 1.45900e-01f;
174  p2 = -1.61359e-02f;
175  p3 = 7.99423e-04f;
176  p4 = -1.47873e-05f;
177  if (x > 16.f) x = 16.f;
178  }
179 
180  else {
181 
182  LogTrace("EgammaSCEnergyCorrectionAlgo")<< "trying to correct unknown cluster!!!" << std::endl;
183  return 1.f;
184  }
185  result = p0 + x*(p1 + x*(p2 + x*(p3 + x*p4)));
186 
187  //Rescale energy scale correction to take into account change in calibrated
188  //RecHit definition introduced in CMSSW_1_5_0
189  float const ebfact = 1.f/0.965f;
190  float const eefact = 1.f/0.975f;
191 
192  if(theBase == EcalBarrel) {
193  result*=ebfact;
194  } else {
195  result*=eefact;
196  }
197 
198  return result;
199 }
200 
202 {
203  // return number of crystals 2Sigma above
204  // electronic noise
205 
206  std::vector<std::pair<DetId,float > > const & hits = seed.hitsAndFractions();
207 
208 
209  LogTrace("EgammaSCEnergyCorrectionAlgo") << " EgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma" << std::endl;
210  LogTrace("EgammaSCEnergyCorrectionAlgo") << " Will calculate number of crystals above 2sigma noise" << std::endl;
211  LogTrace("EgammaSCEnergyCorrectionAlgo") << " sigmaElectronicNoise = " << sigmaElectronicNoise_ << std::endl;
212  LogTrace("EgammaSCEnergyCorrectionAlgo") << " There are " << hits.size() << " recHits" << std::endl;
213 
214  int nCry = 0;
215  std::vector<std::pair<DetId,float > >::const_iterator hit;
217  for(hit = hits.begin(); hit != hits.end(); hit++) {
218  // need to get hit by DetID in order to get energy
219  aHit = rhc.find((*hit).first);
220  if((*aHit).energy()>2.f*sigmaElectronicNoise_) nCry++;
221  }
222 
223  LogTrace("EgammaSCEnergyCorrectionAlgo") << " " << nCry << " of these above 2sigma noise" << std::endl;
224 
225 
226  return nCry;
227 }
228 
229 // apply crack correction
230 
233  EcalClusterFunctionBaseClass* crackCorrectionFunction){
234 
235 
236  double crackcor = 1.;
237 
238  for(reco::CaloCluster_iterator cIt = cl.clustersBegin(); cIt != cl.clustersEnd(); ++cIt) {
239 
240  const reco::CaloClusterPtr cc = *cIt;
241  crackcor *= ( (cl.rawEnergy() +
242  cc->energy()*(crackCorrectionFunction->getValue(*cc)-1.)) /
243  cl.rawEnergy() );
244  }// loop on BCs
245 
246 
247  reco::SuperCluster corrCl=cl;
248  corrCl.setEnergy(cl.energy()*crackcor);
249 
250 
251  return corrCl;
252 }
253 
254 
255 // apply local containment correction
256 // Assume that the correction function provides correction for the seed Basic Cluster
257 
261  EcalClusterFunctionBaseClass* localContCorrectionFunction){
262 
263 
264  const EcalRecHitCollection dummy;
265 
266  const reco::CaloClusterPtr & seedBC = cl.seed();
267  float seedBCene = seedBC->energy();
268  float correctedSeedBCene = localContCorrectionFunction->getValue(*seedBC,dummy) * seedBCene;
269 
270 
271  reco::SuperCluster correctedSC = cl;
272  correctedSC.setEnergy(cl.energy() - seedBCene + correctedSeedBCene);
273 
274  return correctedSC;
275 
276 
277 
278 }
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_)
void Calculate_Covariances(const reco::SuperCluster &passedCluster)
EgammaSCEnergyCorrectionAlgo(float noise, reco::CaloCluster::AlgoId theAlgo, const edm::ParameterSet &pset)
float fNCrystals(int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
std::vector< EcalRecHit >::const_iterator const_iterator
void setEnergy(double energy)
Definition: CaloCluster.h:109
void setPhiWidth(double pw)
Definition: SuperCluster.h:58
virtual float getValue(const reco::BasicCluster &, const EcalRecHitCollection &) const =0
void setEtaWidth(double ew)
Definition: SuperCluster.h:59
double p4[4]
Definition: TauolaWrapper.h:92
tuple result
Definition: query.py:137
double f[11][100]
double energy() const
cluster energy
Definition: CaloCluster.h:120
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
int nCrystalsGT2Sigma(reco::BasicCluster const &seed, EcalRecHitCollection const &rhc) const
float cl
Definition: Combine.cc:71
reco::SuperCluster applyCrackCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *crackCorrectionFunction)
reco::SuperCluster applyLocalContCorrection(const reco::SuperCluster &cl, EcalClusterFunctionBaseClass *localContCorrectionFunction)
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
ESHandle< TrackerGeometry > geometry
iterator find(key_type k)
double p1[4]
Definition: TauolaWrapper.h:89
CaloCluster_iterator clustersBegin() const
fist iterator over BasicCluster constituents
Definition: SuperCluster.h:65
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:62
double preshowerEnergy() const
energy deposited in preshower
Definition: SuperCluster.h:50
x
Definition: VDTMath.h:216
EcalSubdetector
CaloCluster_iterator clustersEnd() const
last iterator over BasicCluster constituents
Definition: SuperCluster.h:68
double p3[4]
Definition: TauolaWrapper.h:91