CMS 3D CMS Logo

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