CMS 3D CMS Logo

CaloGeometryHelper.cc
Go to the documentation of this file.
2 
3 //#include "FWCore/ParameterSet/interface/ParameterSet.h"
4 
5 // needed for the debugging
11 
18 
21 
23 
24 #include <algorithm>
25 
27 {
29  psLayer1Z_ = 303;
30  psLayer2Z_ = 307;
31 }
32 
34 {
35  // std::cout << " In the constructor with ParameterSet " << std::endl;
36  psLayer1Z_ = 303;
37  psLayer2Z_ = 307;
38 }
39 
41 {
44  bfield_ = bField;
46 
48  {
49  ESDetId cps1(getEcalPreshowerGeometry()->getClosestCellInPlane(GlobalPoint(80.,80.,303.),1));
50  psLayer1Z_ = getEcalPreshowerGeometry()->getGeometry(cps1)->getPosition().z();
51  ESDetId cps2(getEcalPreshowerGeometry()->getClosestCellInPlane(GlobalPoint(80.,80.,307.),2));
52  psLayer2Z_ = getEcalPreshowerGeometry()->getGeometry(cps2)->getPosition().z();
53  LogDebug("CaloGeometryTools") << " Preshower layer positions " << psLayer1Z_ << " " << psLayer2Z_ << std::endl;
54  }
55  else
56  LogDebug("CaloGeometryTools") << " No preshower present" << std::endl;
57 
58  // std::cout << " Preshower layer positions " << psLayer1Z_ << " " << psLayer2Z_ << std::endl;
59 
60 }
61 
63 {;
64 }
65 
67 {
68  DetId result;
69  if(ecal)
70  {
71  if(central)
72  {
73  // std::cout << "EcalBarrelGeometry_" << " " << EcalBarrelGeometry_ << std::endl;
74  result = EcalBarrelGeometry_->getClosestCell(GlobalPoint(point.X(),point.Y(),point.Z()));
75 #ifdef DEBUGGCC
76  if(result.null()) return result;
77  GlobalPoint ip=GlobalPoint(point.X(),point.Y(),point.Z());
78  GlobalPoint cc=EcalBarrelGeometry_->getGeometry(result)->getPosition();
79  float deltaeta2 = ip.eta()-cc.eta();
80  deltaeta2 *= deltaeta2;
81  float deltaphi2 = acos(cos(ip.phi()-cc.phi()));
82  deltaphi2 *= deltaphi2;
83  Histos::instance()->fill("h100",point.eta(),sqrt(deltaeta2+deltaphi2));
84 #endif
85  }
86  else
87  {
88  result = EcalEndcapGeometry_->getClosestCell(GlobalPoint(point.X(),point.Y(),point.Z()));
89 #ifdef DEBUGGCC
90  if(result.null())
91  {
92  return result;
93  }
94  GlobalPoint ip=GlobalPoint(point.X(),point.Y(),point.Z());
95  GlobalPoint cc=EcalEndcapGeometry_->getGeometry(result)->getPosition();
96  Histos::instance()->fill("h110",point.eta(),(ip-cc).perp());
97 #endif
98  }
99  }
100  else
101  {
102  result=static_cast<const HcalGeometry*>(HcalGeometry_)->getClosestCell(GlobalPoint(point.X(),point.Y(),point.Z()),true);
103  HcalDetId myDetId(result);
104 
105  // special patch for HF (this is already a part of HcalGeometry)
106  if ( myDetId.subdetId() == HcalForward ) {
107  int mylayer;
108  if ( fabs(point.Z()) > 1132. ) {
109  mylayer = 2;
110  } else {
111  mylayer = 1;
112  }
113  HcalDetId myDetId2((HcalSubdetector)myDetId.subdetId(),myDetId.ieta(),myDetId.iphi(),mylayer);
114  result = myDetId2;
115 // return result;
116  }
117 
118  // Special patch to correct the HCAL geometry (does not work)
119  /*
120  if(result.subdetId()!=HcalEndcap) return result;
121  if(myDetId.depth()==3) return result;
122 
123  int ieta=myDetId.ietaAbs();
124  float azmin=400.458; /// in sync with BaseParticlePropagator
125 
126  if(ieta<=17)
127  return result;
128  else if(ieta>=18 && ieta<=26)
129  azmin += 35.0; // don't consider ieta=18 nose separately
130  else if(ieta>=27)
131  azmin += 21.0;
132 
133  HcalDetId first(HcalEndcap,myDetId.ieta(),myDetId.iphi(),1);
134  bool layer2=(fabs(point.Z())>azmin);
135  if(!layer2)
136  {
137  return first;
138  }
139  else
140  {
141  HcalDetId second(HcalEndcap,myDetId.ieta(),myDetId.iphi(),2);
142  if(second!=HcalDetId()) result=second;
143  }
144  */
145 #ifdef DEBUGGCC
146  if(result.null())
147  {
148  return result;
149  }
150  GlobalPoint ip=GlobalPoint(point.x(),point.y(),point.z());
151  GlobalPoint cc=HcalGeometry_->getPosition(result);
152  float deltaeta2 = ip.eta()-cc.eta();
153  deltaeta2 *= deltaeta2;
154  float deltaphi2 = acos(cos(ip.phi()-cc.phi()));
155  deltaphi2 *= deltaphi2;
156 
157  Histos::instance()->fill("h120",point.eta(),sqrt(deltaeta2+deltaphi2));
158 #endif
159 
160  }
161  return result;
162 }
163 
164 void CaloGeometryHelper::getWindow(const DetId& pivot,int s1,int s2,std::vector<DetId>& vec) const
165 {
166  // currently the getWindow method is the same for EcalBarrelTopology and EndcapTopology
167  // (implemented in CaloSubDetectorTopology)
168  // optimized versions are foreseen
169  vec=getEcalTopology(pivot.subdetId())->getWindow(pivot,s1,s2);
171  sort(vec.begin(),vec.end(),distance);
172 }
173 
174 void CaloGeometryHelper::buildCrystal(const DetId & cell,Crystal& xtal) const
175 {
176  if(cell.subdetId()==EcalBarrel)
177  {
178  xtal=Crystal(cell,&barrelCrystals_[EBDetId(cell).hashedIndex()]);
179  return;
180  }
181  if(cell.subdetId()==EcalEndcap)
182  {
183  xtal=Crystal(cell,&endcapCrystals_[EEDetId(cell).hashedIndex()]);
184  return;
185  }
186 }
187 
188 // Build the array of (max)8 neighbors
190 {
191 
193  NORTHEAST};
194 
195  const unsigned nbarrel = EBDetId::kSizeForDenseIndexing;
196  // Barrel first. The hashed index runs from 0 to 61199
197  barrelNeighbours_.resize(nbarrel);
198 
199  //std::cout << " Building the array of neighbours (barrel) " ;
200 
201  const std::vector<DetId>& vec(EcalBarrelGeometry_->getValidDetIds(DetId::Ecal,EcalBarrel));
202  unsigned size=vec.size();
203  for(unsigned ic=0; ic<size; ++ic)
204  {
205  // We get the 9 cells in a square.
206  std::vector<DetId> neighbours(EcalBarrelTopology_->getWindow(vec[ic],3,3));
207  // std::cout << " Cell " << EBDetId(vec[ic]) << std::endl;
208  unsigned nneighbours=neighbours.size();
209 
210  unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
211  if(hashedindex>=nbarrel)
212  {
213  LogDebug("CaloGeometryTools") << " Array overflow " << std::endl;
214  }
215 
216 
217  // If there are 9 cells, it is easy, and this order is know:
218 // 6 7 8
219 // 3 4 5
220 // 0 1 2 (0 = SOUTHWEST)
221 
222  if(nneighbours==9)
223  {
224  //barrelNeighbours_[hashedindex].reserve(8);
225  unsigned int nn=0;
226  for(unsigned in=0;in<nneighbours;++in)
227  {
228  // remove the centre
229  if(neighbours[in]!=vec[ic])
230  {
231  barrelNeighbours_[hashedindex][nn]=(neighbours[in]);
232  nn++;
233  // std::cout << " Neighbour " << in << " " << EBDetId(neighbours[in]) << std::endl;
234  }
235  }
236  }
237  else
238  {
239  DetId central(vec[ic]);
240  //barrelNeighbours_[hashedindex].resize(8,DetId(0));
241  for(unsigned idir=0;idir<8;++idir)
242  {
243  DetId testid=central;
244  bool status=move(testid,orderedDir[idir],false);
245  if(status) barrelNeighbours_[hashedindex][idir]=testid;
246  }
247 
248  }
249  }
250 
251  // Moved to the endcap
252 
253  // std::cout << " done " << size << std::endl;
254  // std::cout << " Building the array of neighbours (endcap) " ;
255 
256 
257  const std::vector<DetId> & vece(EcalEndcapGeometry_->getValidDetIds(DetId::Ecal,EcalEndcap));
258  size=vece.size();
259  // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
260  // of crystals
261  const unsigned nendcap=EEDetId::kSizeForDenseIndexing;
262 
263  endcapNeighbours_.resize(nendcap);
264  for(unsigned ic=0; ic<size; ++ic)
265  {
266  // We get the 9 cells in a square.
267  std::vector<DetId> neighbours(EcalEndcapTopology_->getWindow(vece[ic],3,3));
268  unsigned nneighbours=neighbours.size();
269  // remove the centre
270  unsigned hashedindex=EEDetId(vece[ic]).hashedIndex();
271 
272  if(hashedindex>=nendcap)
273  {
274  LogDebug("CaloGeometryTools") << " Array overflow " << std::endl;
275  }
276 
277  if(nneighbours==9)
278  {
279  //endcapNeighbours_[hashedindex].reserve(8);
280  unsigned int nn=0;
281  for(unsigned in=0;in<nneighbours;++in)
282  {
283  // remove the centre
284  if(neighbours[in]!=vece[ic])
285  {
286  endcapNeighbours_[hashedindex][nn]=(neighbours[in]);
287  nn++;
288  }
289  }
290  }
291  else
292  {
293  DetId central(vece[ic]);
294  //endcapNeighbours_[hashedindex].resize(8,DetId(0));
295  for(unsigned idir=0;idir<8;++idir)
296  {
297  DetId testid=central;
298  bool status=move(testid,orderedDir[idir],false);
299  if(status) endcapNeighbours_[hashedindex][idir]=testid;
300  }
301 
302  }
303  }
304  // std::cout << " done " << size <<std::endl;
306 }
307 
309 {
310  return (detid.subdetId()==EcalBarrel)?barrelNeighbours_[EBDetId(detid).hashedIndex()]:
312 }
313 
314 bool CaloGeometryHelper::move(DetId& cell, const CaloDirection&dir,bool fast) const
315 {
316  DetId originalcell = cell;
317  if(dir==NONE || cell==DetId(0)) return false;
318 
319  // Conversion CaloDirection and index in the table
320  // CaloDirection :NONE,SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST, NORTHEAST,NORTHWEST,NORTH
321  // Table : SOUTHWEST,SOUTH,SOUTHEAST,WEST,EAST,NORTHWEST,NORTH, NORTHEAST
322  static const int calodirections[9]={-1,1,2,0,4,3,7,5,6};
323 
324  if(fast&&neighbourmapcalculated_)
325  {
326  DetId result = (originalcell.subdetId()==EcalBarrel) ?
327  barrelNeighbours_[EBDetId(originalcell).hashedIndex()][calodirections[dir]]:
328  endcapNeighbours_[EEDetId(originalcell).hashedIndex()][calodirections[dir]];
329  bool status = !result.null();
330  cell = result;
331  return status;
332  }
333 
334  if(dir==NORTH || dir ==SOUTH || dir==EAST || dir==WEST)
335  {
336  return simplemove(cell,dir);
337  }
338  else
339  {
340  if(dir == NORTHEAST || dir==NORTHWEST || dir==SOUTHEAST || dir==SOUTHWEST)
341  return diagonalmove(cell,dir);
342  }
343 
344  cell = DetId(0);
345  return false;
346 }
347 
348 
350 {
351  std::vector<DetId> neighbours;
352  if(cell.subdetId()==EcalBarrel)
353  neighbours = EcalBarrelTopology_->getNeighbours(cell,dir);
354  else if(cell.subdetId()==EcalEndcap)
355  neighbours= EcalEndcapTopology_->getNeighbours(cell,dir);
356 
357  if ((!neighbours.empty()) && (!neighbours[0].null()))
358  {
359  cell = neighbours[0];
360  return true;
361  }
362  else
363  {
364  cell = DetId(0);
365  return false;
366  }
367 }
368 
370 {
371  bool result;
372  // One has to try both paths
373  if(dir==NORTHEAST)
374  {
375  result = simplemove(cell,NORTH);
376  if(result)
377  return simplemove(cell,EAST);
378  else
379  {
380  result = simplemove(cell,EAST);
381  if(result)
382  return simplemove(cell,NORTH);
383  else
384  return false;
385  }
386  }
387  else if(dir==NORTHWEST)
388  {
389  result = simplemove(cell,NORTH);
390  if(result)
391  return simplemove(cell,WEST);
392  else
393  {
394  result = simplemove(cell,WEST);
395  if(result)
396  return simplemove(cell,NORTH);
397  else
398  return false;
399  }
400  }
401  else if(dir == SOUTHEAST)
402  {
403  result = simplemove(cell,SOUTH);
404  if(result)
405  return simplemove(cell,EAST);
406  else
407  {
408  result = simplemove(cell,EAST);
409  if(result)
410  return simplemove(cell,SOUTH);
411  else
412  return false;
413  }
414  }
415  else if(dir == SOUTHWEST)
416  {
417  result = simplemove(cell,SOUTH);
418  if(result)
419  return simplemove(cell,WEST);
420  else
421  {
422  result = simplemove(cell,SOUTH);
423  if(result)
424  return simplemove(cell,WEST);
425  else
426  return false;
427  }
428  }
429  cell = DetId(0);
430  return false;
431 }
432 
433 bool CaloGeometryHelper::borderCrossing(const DetId& c1, const DetId& c2) const
434 {
435  if(c1.subdetId()!=c2.subdetId()) return false;
436 
437  if(c1.subdetId()==EcalBarrel)
438  {
439  // there is a crack if the two cells don't belong to the same
440  // module
441  EBDetId cc1(c1);
442  EBDetId cc2(c2);
443  return (cc1.im()!=cc2.im()||cc1.ism()!=cc2.ism() );
444  }
445 
446 if(c1.subdetId()==EcalEndcap)
447  {
448  // there is a crack if the two cells don't belong to the same
449  // module
450  return (EEDetId(c1).isc()!=EEDetId(c2).isc());
451  }
452  return false;
453 }
454 
456 {
457  const unsigned nbarrel = EBDetId::kSizeForDenseIndexing;
458  // Barrel first. The hashed index runs from 0 to 61199
459  barrelCrystals_.resize(nbarrel,BaseCrystal());
460 
461  //std::cout << " Building the array of crystals (barrel) " ;
462  const std::vector<DetId>& vec(EcalBarrelGeometry_->getValidDetIds(DetId::Ecal,EcalBarrel));
463  unsigned size=vec.size();
464  for(unsigned ic=0; ic<size; ++ic)
465  {
466  unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
467  auto geom = EcalBarrelGeometry_->getGeometry(vec[ic]);
468  BaseCrystal xtal(vec[ic]);
469  xtal.setCorners(geom->getCorners(),geom->getPosition());
470  barrelCrystals_[hashedindex]=xtal;
471  }
472 
473  // std::cout << " done " << size << std::endl;
474  // std::cout << " Building the array of crystals (endcap) " ;
475 
476 
477  const std::vector<DetId>& vece(EcalEndcapGeometry_->getValidDetIds(DetId::Ecal,EcalEndcap));
478  size=vece.size();
479  // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
480  // of crystals
481  const unsigned nendcap=EEDetId::kSizeForDenseIndexing;
482 
483  endcapCrystals_.resize(nendcap,BaseCrystal());
484  for(unsigned ic=0; ic<size; ++ic)
485  {
486  unsigned hashedindex=EEDetId(vece[ic]).hashedIndex();
487  auto geom = EcalEndcapGeometry_->getGeometry(vece[ic]);
488  BaseCrystal xtal(vece[ic]);
489  xtal.setCorners(geom->getCorners(),geom->getPosition());
490  endcapCrystals_[hashedindex]=xtal;
491  }
492  // std::cout << " done " << size << std::endl;
493 }
#define LogDebug(id)
size
Write out results.
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:82
void setCorners(const CaloCellGeometry::CornersVec &vec, const GlobalPoint &pos)
Definition: BaseCrystal.cc:14
virtual std::vector< DetId > getNeighbours(const DetId &id, const CaloDirection &dir) const
void buildCrystal(const DetId &id, Crystal &) const
DetId getClosestCell(const GlobalPoint &r) const override
constexpr bool null() const
is this a null id ?
Definition: DetId.h:52
const EcalEndcapGeometry * EcalEndcapGeometry_
Definition: Calorimeter.h:78
int isc() const
Definition: EEDetId.cc:222
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
const EcalPreshowerGeometry * getEcalPreshowerGeometry() const
Definition: Calorimeter.h:55
const CaloSubdetectorTopology * EcalEndcapTopology_
Definition: Calorimeter.h:84
std::vector< NeiVect > endcapNeighbours_
bool borderCrossing(const DetId &, const DetId &) const
int ism() const
get the ECAL/SM id
Definition: EBDetId.h:59
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
const CaloSubdetectorGeometry * getEcalGeometry(int subdetn) const
Definition: Calorimeter.cc:136
math::XYZVector XYZPoint
static const CaloDirection orderedDir[8]
int im() const
get the number of module inside the SM (1-4)
Definition: EBDetId.h:64
int hashedIndex(int ieta, int iphi)
Definition: EcalPyUtils.cc:42
void fill(const std::string &name, float val1, float val2=1., float val3=1.)
Fill an histogram.
Definition: Histos.cc:168
T sqrt(T t)
Definition: SSEVec.h:18
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int ieta() const
get the cell ieta
Definition: HcalDetId.h:159
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::vector< BaseCrystal > endcapCrystals_
static Histos * instance()
Definition: Histos.cc:15
HcalSubdetector
Definition: HcalAssistant.h:31
const EcalBarrelGeometry * EcalBarrelGeometry_
Definition: Calorimeter.h:77
const CaloSubdetectorTopology * getEcalTopology(int subdetn) const
Definition: Calorimeter.cc:145
std::vector< NeiVect > barrelNeighbours_
int iphi() const
get the cell iphi
Definition: HcalDetId.h:161
Definition: DetId.h:18
bool diagonalmove(DetId &cell, const CaloDirection &dir) const
int hashedIndex() const
Definition: EEDetId.h:183
virtual std::vector< DetId > getWindow(const DetId &id, const int &northSouthSize, const int &eastWestSize) const
std::array< DetId, 8 > NeiVect
virtual std::shared_ptr< const CaloCellGeometry > getGeometry(const DetId &id) const
Get the cell geometry of a given detector id. Should return false if not found.
void getWindow(const DetId &pivot, int s1, int s2, std::vector< DetId > &) const
const CaloSubdetectorGeometry * HcalGeometry_
Definition: Calorimeter.h:79
T eta() const
Definition: PV3DBase.h:76
T perp() const
Magnitude of transverse component.
const CaloSubdetectorTopology * EcalBarrelTopology_
Definition: Calorimeter.h:83
std::vector< BaseCrystal > barrelCrystals_
dbl *** dir
Definition: mlp_gen.cc:35
CaloDirection
Codes the local directions in the cell lattice.
Definition: CaloDirection.h:9
DetId getClosestCell(const GlobalPoint &r) const override
DetId getClosestCell(const XYZPoint &point, bool ecal, bool central) const
void initialize(double bField)
bool simplemove(DetId &cell, const CaloDirection &dir) const
*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 NeiVect & getNeighbours(const DetId &det) const
bool move(DetId &cell, const CaloDirection &dir, bool fast=true) const