CMS 3D CMS Logo

EcalClusterLazyTools.cc
Go to the documentation of this file.
2 
5 
8 
16 
20 
24 
26  const edm::EventSetup &es,
29  : geometry_(&edm::get<CaloGeometry, CaloGeometryRecord>(es)),
30  topology_(&edm::get<CaloTopology, CaloTopologyRecord>(es)),
31  ebRecHits_(&edm::get(ev, token1)),
32  eeRecHits_(&edm::get(ev, token2)) {
34  getADCToGeV(es);
36 }
37 
39  const edm::EventSetup &es,
45  ebRecHits_(&edm::get(ev, token1)),
46  eeRecHits_(&edm::get(ev, token2)) {
48  if (geometryES) {
49  ecalPS_topology_ = std::make_shared<EcalPreshowerTopology>();
50  } else {
51  edm::LogInfo("subdetector geometry not available") << "EcalPreshower geometry is missing" << std::endl;
52  }
53 
54  esRHToken_ = token3;
55  getESRecHits(ev);
56 
58  getADCToGeV(es);
60 }
61 
64  ev.getByToken(esRHToken_, pESRecHits);
65  esRecHits_ = pESRecHits.product();
66  // make the map of rechits
67  rechits_map_.clear();
68  if (pESRecHits.isValid()) {
70  for (it = pESRecHits->begin(); it != pESRecHits->end(); ++it) {
71  // remove bad ES rechits
72  std::vector<int> badf = {
73  EcalRecHit::ESFlags::kESDead, // 1
74  EcalRecHit::ESFlags::kESTwoGoodRatios,
75  EcalRecHit::ESFlags::kESBadRatioFor12, // 5
76  EcalRecHit::ESFlags::kESBadRatioFor23Upper,
77  EcalRecHit::ESFlags::kESBadRatioFor23Lower,
78  EcalRecHit::ESFlags::kESTS1Largest,
79  EcalRecHit::ESFlags::kESTS3Largest,
80  EcalRecHit::ESFlags::kESTS3Negative, // 10
81  EcalRecHit::ESFlags::kESTS13Sigmas, // 14
82  };
83 
84  if (it->checkFlags(badf))
85  continue;
86 
87  //Make the map of DetID, EcalRecHit pairs
88  rechits_map_.insert(std::make_pair(it->id(), *it));
89  }
90  }
91 }
92 
94  // get IC's
96  icalMap = &ical->getMap();
97 }
98 
100  // get ADCtoGeV
101  es.get<EcalADCToGeVConstantRcd>().get(agc);
102 }
103 
105  // transp corrections
106  es.get<EcalLaserDbRecord>().get(laser);
107 }
108 
110  if (cluster.size() == 0) {
111  throw cms::Exception("InvalidCluster") << "The cluster has no crystals!";
112  }
113  DetId id = (cluster.hitsAndFractions()[0]).first; // size is by definition > 0 -- FIXME??
114  const EcalRecHitCollection *recHits = nullptr;
115  if (id.subdetId() == EcalBarrel) {
116  recHits = ebRecHits_;
117  } else if (id.subdetId() == EcalEndcap) {
118  recHits = eeRecHits_;
119  } else {
120  throw cms::Exception("InvalidSubdetector")
121  << "The subdetId() " << id.subdetId() << " does not correspond to EcalBarrel neither EcalEndcap";
122  }
123  return recHits;
124 }
125 
126 // get time of basic cluster seed crystal
129 
130  DetId id = cluster.seed();
131  EcalRecHitCollection::const_iterator theSeedHit = recHits->find(id);
132  // std::cout << "the seed of the BC has time: "
133  //<< (*theSeedHit).time()
134  //<< "and energy: " << (*theSeedHit).energy() << " collection size: " << recHits->size()
135  //<< "\n" <<std::endl; // GF debug
136 
137  return (*theSeedHit).time();
138 }
139 
140 // error-weighted average of time from constituents of basic cluster
142  std::vector<std::pair<DetId, float> > clusterComponents = (cluster).hitsAndFractions();
143  //std::cout << "BC has this many components: " << clusterComponents.size() << std::endl; // GF debug
144 
146  //std::cout << "BasicClusterClusterTime - rechits are this many: " << recHits->size() << std::endl; // GF debug
147 
148  float weightedTsum = 0;
149  float sumOfWeights = 0;
150 
151  for (std::vector<std::pair<DetId, float> >::const_iterator detitr = clusterComponents.begin();
152  detitr != clusterComponents.end();
153  detitr++) {
154  // EcalRecHitCollection::const_iterator theSeedHit = recHits->find (id); // trash this
155  EcalRecHitCollection::const_iterator oneHit = recHits->find((detitr->first));
156 
157  // in order to get back the ADC counts from the recHit energy, three ingredients are necessary:
158  // 1) get laser correction coefficient
159  float lasercalib = 1.;
160  lasercalib = laser->getLaserCorrection(detitr->first, ev.time());
161  // 2) get intercalibration
162  EcalIntercalibConstantMap::const_iterator icalit = icalMap->find(detitr->first);
163  EcalIntercalibConstant icalconst = 1.;
164  if (icalit != icalMap->end()) {
165  icalconst = (*icalit);
166  // std::cout << "icalconst set to: " << icalconst << std::endl;
167  } else {
168  edm::LogError("EcalClusterLazyTools")
169  << "No intercalib const found for xtal " << (detitr->first).rawId() << "bailing out";
170  assert(false);
171  }
172  // 3) get adc2GeV
173  float adcToGeV = 1.;
174  if ((detitr->first).subdetId() == EcalBarrel)
175  adcToGeV = float(agc->getEBValue());
176  else if ((detitr->first).subdetId() == EcalEndcap)
177  adcToGeV = float(agc->getEEValue());
178  float adc = 2.;
179  if (icalconst > 0 && lasercalib > 0 && adcToGeV > 0)
180  adc = (*oneHit).energy() / (icalconst * lasercalib * adcToGeV);
181 
182  // don't consider recHits with too little amplitude; take sigma_noise_total into account
183  if ((detitr->first).subdetId() == EcalBarrel && adc < (1.1 * 20))
184  continue;
185  if ((detitr->first).subdetId() == EcalEndcap && adc < (2.2 * 20))
186  continue;
187 
188  // count only on rechits whose error is trusted by the method (ratio)
189  if (!(*oneHit).isTimeErrorValid())
190  continue;
191 
192  float timeError = (*oneHit).timeError();
193  // the constant used to build timeError is largely over-estimated ; remove in quadrature 0.6 and add 0.15 back.
194  // could be prettier if value of constant term was changed at recHit production level
195  if (timeError > 0.6)
196  timeError = sqrt(timeError * timeError - 0.6 * 0.6 + 0.15 * 0.15);
197  else
198  timeError = sqrt(timeError * timeError + 0.15 * 0.15);
199 
200  // do the error weighting
201  weightedTsum += (*oneHit).time() / (timeError * timeError);
202  sumOfWeights += 1. / (timeError * timeError);
203  }
204 
205  // what if no crytal is available for weighted average?
206  if (sumOfWeights == 0)
207  return -999;
208  else
209  return (weightedTsum / sumOfWeights);
210 }
211 
212 // get BasicClusterSeedTime of the seed basic cluser of the supercluster
214  return BasicClusterSeedTime((*cluster.seed()));
215 }
216 
217 // get BasicClusterTime of the seed basic cluser of the supercluster
219  return BasicClusterTime((*cluster.seed()), ev);
220 }
221 
222 // get Preshower effective sigmaIRIR
224  if (!(fabs(cluster.eta()) > 1.6 && fabs(cluster.eta()) < 3.))
225  return 0.;
226 
227  if (!ecalPS_topology_)
228  return 0.;
229 
230  std::vector<float> phoESHitsIXIX =
231  getESHits(cluster.x(), cluster.y(), cluster.z(), rechits_map_, geometry_, ecalPS_topology_.get(), 0, 1);
232  std::vector<float> phoESHitsIYIY =
233  getESHits(cluster.x(), cluster.y(), cluster.z(), rechits_map_, geometry_, ecalPS_topology_.get(), 0, 2);
234  float phoESShapeIXIX = getESShape(phoESHitsIXIX);
235  float phoESShapeIYIY = getESShape(phoESHitsIYIY);
236 
237  return sqrt(phoESShapeIXIX * phoESShapeIXIX + phoESShapeIYIY * phoESShapeIYIY);
238 }
239 
240 // get Preshower effective sigmaIXIX
242  if (!(fabs(cluster.eta()) > 1.6 && fabs(cluster.eta()) < 3.))
243  return 0.;
244 
245  if (!ecalPS_topology_)
246  return 0.;
247 
248  std::vector<float> phoESHitsIXIX =
249  getESHits(cluster.x(), cluster.y(), cluster.z(), rechits_map_, geometry_, ecalPS_topology_.get(), 0, 1);
250  float phoESShapeIXIX = getESShape(phoESHitsIXIX);
251 
252  return phoESShapeIXIX;
253 }
254 
255 // get Preshower effective sigmaIYIY
257  if (!(fabs(cluster.eta()) > 1.6 && fabs(cluster.eta()) < 3.))
258  return 0.;
259 
260  if (!ecalPS_topology_)
261  return 0.;
262 
263  std::vector<float> phoESHitsIYIY =
264  getESHits(cluster.x(), cluster.y(), cluster.z(), rechits_map_, geometry_, ecalPS_topology_.get(), 0, 2);
265  float phoESShapeIYIY = getESShape(phoESHitsIYIY);
266 
267  return phoESShapeIYIY;
268 }
269 
270 // get Preshower Rechits
271 std::vector<float> EcalClusterLazyToolsBase::getESHits(double X,
272  double Y,
273  double Z,
274  const std::map<DetId, EcalRecHit> &_rechits_map,
275  const CaloGeometry *geometry,
276  CaloSubdetectorTopology const *topology_p,
277  int row,
278  int plane) {
279  std::map<DetId, EcalRecHit> rechits_map = _rechits_map;
280  std::vector<float> esHits;
281 
282  const GlobalPoint point(X, Y, Z);
283 
285 
286  DetId esId = (dynamic_cast<const EcalPreshowerGeometry *>(geometry_p))->getClosestCellInPlane(point, plane);
287  ESDetId esDetId = (esId == DetId(0)) ? ESDetId(0) : ESDetId(esId);
288 
289  std::map<DetId, EcalRecHit>::iterator it;
290  ESDetId next;
291  ESDetId strip;
292  strip = esDetId;
293 
294  EcalPreshowerNavigator theESNav(strip, topology_p);
295  theESNav.setHome(strip);
296 
297  if (row == 1) {
298  if (plane == 1 && strip != ESDetId(0))
299  strip = theESNav.north();
300  if (plane == 2 && strip != ESDetId(0))
301  strip = theESNav.east();
302  } else if (row == -1) {
303  if (plane == 1 && strip != ESDetId(0))
304  strip = theESNav.south();
305  if (plane == 2 && strip != ESDetId(0))
306  strip = theESNav.west();
307  }
308 
309  if (strip == ESDetId(0)) {
310  for (int i = 0; i < 31; ++i)
311  esHits.push_back(0);
312  } else {
313  it = rechits_map.find(strip);
314  if (it != rechits_map.end() && it->second.energy() > 1.0e-10)
315  esHits.push_back(it->second.energy());
316  else
317  esHits.push_back(0);
318  //cout<<"center : "<<strip<<" "<<it->second.energy()<<endl;
319 
320  // Front Plane
321  if (plane == 1) {
322  // east road
323  for (int i = 0; i < 15; ++i) {
324  next = theESNav.east();
325  if (next != ESDetId(0)) {
326  it = rechits_map.find(next);
327  if (it != rechits_map.end() && it->second.energy() > 1.0e-10)
328  esHits.push_back(it->second.energy());
329  else
330  esHits.push_back(0);
331  //cout<<"east "<<i<<" : "<<next<<" "<<it->second.energy()<<endl;
332  } else {
333  for (int j = i; j < 15; j++)
334  esHits.push_back(0);
335  break;
336  //cout<<"east "<<i<<" : "<<next<<" "<<0<<endl;
337  }
338  }
339 
340  // west road
341  theESNav.setHome(strip);
342  theESNav.home();
343  for (int i = 0; i < 15; ++i) {
344  next = theESNav.west();
345  if (next != ESDetId(0)) {
346  it = rechits_map.find(next);
347  if (it != rechits_map.end() && it->second.energy() > 1.0e-10)
348  esHits.push_back(it->second.energy());
349  else
350  esHits.push_back(0);
351  //cout<<"west "<<i<<" : "<<next<<" "<<it->second.energy()<<endl;
352  } else {
353  for (int j = i; j < 15; j++)
354  esHits.push_back(0);
355  break;
356  //cout<<"west "<<i<<" : "<<next<<" "<<0<<endl;
357  }
358  }
359  } // End of Front Plane
360 
361  // Rear Plane
362  if (plane == 2) {
363  // north road
364  for (int i = 0; i < 15; ++i) {
365  next = theESNav.north();
366  if (next != ESDetId(0)) {
367  it = rechits_map.find(next);
368  if (it != rechits_map.end() && it->second.energy() > 1.0e-10)
369  esHits.push_back(it->second.energy());
370  else
371  esHits.push_back(0);
372  //cout<<"north "<<i<<" : "<<next<<" "<<it->second.energy()<<endl;
373  } else {
374  for (int j = i; j < 15; j++)
375  esHits.push_back(0);
376  break;
377  //cout<<"north "<<i<<" : "<<next<<" "<<0<<endl;
378  }
379  }
380 
381  // south road
382  theESNav.setHome(strip);
383  theESNav.home();
384  for (int i = 0; i < 15; ++i) {
385  next = theESNav.south();
386  if (next != ESDetId(0)) {
387  it = rechits_map.find(next);
388  if (it != rechits_map.end() && it->second.energy() > 1.0e-10)
389  esHits.push_back(it->second.energy());
390  else
391  esHits.push_back(0);
392  //cout<<"south "<<i<<" : "<<next<<" "<<it->second.energy()<<endl;
393  } else {
394  for (int j = i; j < 15; j++)
395  esHits.push_back(0);
396  break;
397  //cout<<"south "<<i<<" : "<<next<<" "<<0<<endl;
398  }
399  }
400  } // End of Rear Plane
401  } // Fill ES RecHits
402 
403  return esHits;
404 }
405 
406 // get Preshower hit shape
407 float EcalClusterLazyToolsBase::getESShape(const std::vector<float> &ESHits0) {
408  const int nBIN = 21;
409  float esRH[nBIN];
410  for (int idx = 0; idx < nBIN; idx++) {
411  esRH[idx] = 0.;
412  }
413 
414  for (int ibin = 0; ibin < ((nBIN + 1) / 2); ibin++) {
415  if (ibin == 0) {
416  esRH[(nBIN - 1) / 2] = ESHits0[ibin];
417  } else {
418  esRH[(nBIN - 1) / 2 + ibin] = ESHits0[ibin];
419  esRH[(nBIN - 1) / 2 - ibin] = ESHits0[ibin + 15];
420  }
421  }
422 
423  // ---- Effective Energy Deposit Width ---- //
424  double EffWidthSigmaISIS = 0.;
425  double totalEnergyISIS = 0.;
426  double EffStatsISIS = 0.;
427  for (int id_X = 0; id_X < 21; id_X++) {
428  totalEnergyISIS += esRH[id_X];
429  EffStatsISIS += esRH[id_X] * (id_X - 10) * (id_X - 10);
430  }
431  EffWidthSigmaISIS = (totalEnergyISIS > 0.) ? sqrt(fabs(EffStatsISIS / totalEnergyISIS)) : 0.;
432 
433  return EffWidthSigmaISIS;
434 }
edm::ESHandle< EcalIntercalibConstants > ical
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:34
const CaloGeometry * geometry_
const EcalRecHitCollection * esRecHits_
std::vector< float > getESHits(double X, double Y, double Z, const std::map< DetId, EcalRecHit > &rechits_map, const CaloGeometry *geometry, CaloSubdetectorTopology const *topology_p, int row=0, int plane=1)
float SuperClusterSeedTime(const reco::SuperCluster &cluster)
float BasicClusterTime(const reco::BasicCluster &cluster, const edm::Event &ev)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
float SuperClusterTime(const reco::SuperCluster &cluster, const edm::Event &ev)
float getESShape(const std::vector< float > &ESHits0)
float getLaserCorrection(DetId const &xid, edm::Timestamp const &iTime) const
std::vector< EcalRecHit >::const_iterator const_iterator
const EcalIntercalibConstantMap * icalMap
#define X(str)
Definition: MuonsGrabber.cc:38
void getADCToGeV(const edm::EventSetup &es)
bool ev
const EcalRecHitCollection * ebRecHits_
float eseffsixix(const reco::SuperCluster &cluster)
void getIntercalibConstants(const edm::EventSetup &es)
double z() const
z coordinate of cluster centroid
Definition: CaloCluster.h:177
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:180
void setHome(const T &startingPoint)
set the starting position
Definition: CaloNavigator.h:91
edm::ESHandle< EcalADCToGeVConstant > agc
T west() const
move the navigator west
Definition: CaloNavigator.h:49
T sqrt(T t)
Definition: SSEVec.h:19
EcalClusterLazyToolsBase(const edm::Event &ev, const edm::EventSetup &es, edm::EDGetTokenT< EcalRecHitCollection > token1, edm::EDGetTokenT< EcalRecHitCollection > token2)
double x() const
x coordinate of cluster centroid
Definition: CaloCluster.h:171
T south() const
move the navigator south
Definition: CaloNavigator.h:37
constexpr int adc(sample_type sample)
get the ADC sample (12 bits)
void getESRecHits(const edm::Event &ev)
bool isValid() const
Definition: HandleBase.h:70
edm::EDGetTokenT< EcalRecHitCollection > esRHToken_
const_iterator end() const
T east() const
move the navigator east
Definition: CaloNavigator.h:43
void home() const
move the navigator back to the starting point
Definition: CaloNavigator.h:97
Definition: DetId.h:17
float eseffsiyiy(const reco::SuperCluster &cluster)
std::shared_ptr< CaloSubdetectorTopology const > ecalPS_topology_
T const * product() const
Definition: Handle.h:69
EcalRecHitCollection const * getEcalRecHitCollection(const reco::BasicCluster &cluster) const
std::vector< Item >::const_iterator const_iterator
iterator find(key_type k)
HLT enums.
void getLaserDbService(const edm::EventSetup &es)
T get() const
Definition: EventSetup.h:73
const EcalRecHitCollection * eeRecHits_
const_iterator find(uint32_t rawId) const
T north() const
move the navigator north
Definition: CaloNavigator.h:31
const_iterator end() const
const CaloClusterPtr & seed() const
seed BasicCluster
Definition: SuperCluster.h:77
double y() const
y coordinate of cluster centroid
Definition: CaloCluster.h:174
std::map< DetId, EcalRecHit > rechits_map_
edm::ESHandle< EcalLaserDbService > laser
edm::Timestamp time() const
Definition: EventBase.h:60
float BasicClusterSeedTime(const reco::BasicCluster &cluster)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
const_iterator begin() const
const CaloTopology * topology_
float EcalIntercalibConstant
float eseffsirir(const reco::SuperCluster &cluster)