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 
66 DetId CaloGeometryHelper::getClosestCell(const XYZPoint& point, bool ecal, bool central) const
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=HcalGeometry_->getClosestCell(GlobalPoint(point.X(),point.Y(),point.Z()));
103  HcalDetId myDetId(result);
104 
105  // special patch for HF
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 
119  if(result.subdetId()!=HcalEndcap) return result;
120  // Special patch to correct the HCAL geometry
121  if(myDetId.depth()==3) return result;
122 
123  int ieta=myDetId.ietaAbs();
124  float azmin=400.458;
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 #ifdef DEBUGGCC
145  if(result.null())
146  {
147  return result;
148  }
149  GlobalPoint ip=GlobalPoint(point.x(),point.y(),point.z());
150  GlobalPoint cc=HcalGeometry_->getPosition(result);
151  float deltaeta2 = ip.eta()-cc.eta();
152  deltaeta2 *= deltaeta2;
153  float deltaphi2 = acos(cos(ip.phi()-cc.phi()));
154  deltaphi2 *= deltaphi2;
155 
156  Histos::instance()->fill("h120",point.eta(),sqrt(deltaeta2+deltaphi2));
157 #endif
158 
159  }
160  return result;
161 }
162 
163 void CaloGeometryHelper::getWindow(const DetId& pivot,int s1,int s2,std::vector<DetId>& vec) const
164 {
165  // currently the getWindow method is the same for EcalBarrelTopology and EndcapTopology
166  // (implemented in CaloSubDetectorTopology)
167  // optimized versions are foreseen
168  vec=getEcalTopology(pivot.subdetId())->getWindow(pivot,s1,s2);
170  sort(vec.begin(),vec.end(),distance);
171 }
172 
173 void CaloGeometryHelper::buildCrystal(const DetId & cell,Crystal& xtal) const
174 {
175  if(cell.subdetId()==EcalBarrel)
176  {
177  xtal=Crystal(cell,&barrelCrystals_[EBDetId(cell).hashedIndex()]);
178  return;
179  }
180  if(cell.subdetId()==EcalEndcap)
181  {
182  xtal=Crystal(cell,&endcapCrystals_[EEDetId(cell).hashedIndex()]);
183  return;
184  }
185 }
186 
187 // Build the array of (max)8 neighbors
189 {
190 
192  NORTHEAST};
193 
194  const unsigned nbarrel = EBDetId::kSizeForDenseIndexing;
195  // Barrel first. The hashed index runs from 0 to 61199
196  barrelNeighbours_.resize(nbarrel);
197 
198  //std::cout << " Building the array of neighbours (barrel) " ;
199 
200  const std::vector<DetId>& vec(EcalBarrelGeometry_->getValidDetIds(DetId::Ecal,EcalBarrel));
201  unsigned size=vec.size();
202  for(unsigned ic=0; ic<size; ++ic)
203  {
204  // We get the 9 cells in a square.
205  std::vector<DetId> neighbours(EcalBarrelTopology_->getWindow(vec[ic],3,3));
206  // std::cout << " Cell " << EBDetId(vec[ic]) << std::endl;
207  unsigned nneighbours=neighbours.size();
208 
209  unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
210  if(hashedindex>=nbarrel)
211  {
212  LogDebug("CaloGeometryTools") << " Array overflow " << std::endl;
213  }
214 
215 
216  // If there are 9 cells, it is easy, and this order is know:
217 // 6 7 8
218 // 3 4 5
219 // 0 1 2 (0 = SOUTHWEST)
220 
221  if(nneighbours==9)
222  {
223  //barrelNeighbours_[hashedindex].reserve(8);
224  unsigned int nn=0;
225  for(unsigned in=0;in<nneighbours;++in)
226  {
227  // remove the centre
228  if(neighbours[in]!=vec[ic])
229  {
230  barrelNeighbours_[hashedindex][nn]=(neighbours[in]);
231  nn++;
232  // std::cout << " Neighbour " << in << " " << EBDetId(neighbours[in]) << std::endl;
233  }
234  }
235  }
236  else
237  {
238  DetId central(vec[ic]);
239  //barrelNeighbours_[hashedindex].resize(8,DetId(0));
240  for(unsigned idir=0;idir<8;++idir)
241  {
242  DetId testid=central;
243  bool status=move(testid,orderedDir[idir],false);
244  if(status) barrelNeighbours_[hashedindex][idir]=testid;
245  }
246 
247  }
248  }
249 
250  // Moved to the endcap
251 
252  // std::cout << " done " << size << std::endl;
253  // std::cout << " Building the array of neighbours (endcap) " ;
254 
255 
256  const std::vector<DetId> & vece(EcalEndcapGeometry_->getValidDetIds(DetId::Ecal,EcalEndcap));
257  size=vece.size();
258  // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
259  // of crystals
260  const unsigned nendcap=EEDetId::kSizeForDenseIndexing;
261 
262  endcapNeighbours_.resize(nendcap);
263  for(unsigned ic=0; ic<size; ++ic)
264  {
265  // We get the 9 cells in a square.
266  std::vector<DetId> neighbours(EcalEndcapTopology_->getWindow(vece[ic],3,3));
267  unsigned nneighbours=neighbours.size();
268  // remove the centre
269  unsigned hashedindex=EEDetId(vece[ic]).hashedIndex();
270 
271  if(hashedindex>=nendcap)
272  {
273  LogDebug("CaloGeometryTools") << " Array overflow " << std::endl;
274  }
275 
276  if(nneighbours==9)
277  {
278  //endcapNeighbours_[hashedindex].reserve(8);
279  unsigned int nn=0;
280  for(unsigned in=0;in<nneighbours;++in)
281  {
282  // remove the centre
283  if(neighbours[in]!=vece[ic])
284  {
285  endcapNeighbours_[hashedindex][nn]=(neighbours[in]);
286  nn++;
287  }
288  }
289  }
290  else
291  {
292  DetId central(vece[ic]);
293  //endcapNeighbours_[hashedindex].resize(8,DetId(0));
294  for(unsigned idir=0;idir<8;++idir)
295  {
296  DetId testid=central;
297  bool status=move(testid,orderedDir[idir],false);
298  if(status) endcapNeighbours_[hashedindex][idir]=testid;
299  }
300 
301  }
302  }
303  // std::cout << " done " << size <<std::endl;
305 }
306 
308 {
309  return (detid.subdetId()==EcalBarrel)?barrelNeighbours_[EBDetId(detid).hashedIndex()]:
311 }
312 
313 bool CaloGeometryHelper::move(DetId& cell, const CaloDirection&dir,bool fast) const
314 {
315  DetId originalcell = cell;
316  if(dir==NONE || cell==DetId(0)) return false;
317 
318  // Conversion CaloDirection and index in the table
319  // CaloDirection :NONE,SOUTH,SOUTHEAST,SOUTHWEST,EAST,WEST, NORTHEAST,NORTHWEST,NORTH
320  // Table : SOUTHWEST,SOUTH,SOUTHEAST,WEST,EAST,NORTHWEST,NORTH, NORTHEAST
321  static const int calodirections[9]={-1,1,2,0,4,3,7,5,6};
322 
323  if(fast&&neighbourmapcalculated_)
324  {
325  DetId result = (originalcell.subdetId()==EcalBarrel) ?
326  barrelNeighbours_[EBDetId(originalcell).hashedIndex()][calodirections[dir]]:
327  endcapNeighbours_[EEDetId(originalcell).hashedIndex()][calodirections[dir]];
328  bool status = !result.null();
329  cell = result;
330  return status;
331  }
332 
333  if(dir==NORTH || dir ==SOUTH || dir==EAST || dir==WEST)
334  {
335  return simplemove(cell,dir);
336  }
337  else
338  {
339  if(dir == NORTHEAST || dir==NORTHWEST || dir==SOUTHEAST || dir==SOUTHWEST)
340  return diagonalmove(cell,dir);
341  }
342 
343  cell = DetId(0);
344  return false;
345 }
346 
347 
349 {
350  std::vector<DetId> neighbours;
351  if(cell.subdetId()==EcalBarrel)
352  neighbours = EcalBarrelTopology_->getNeighbours(cell,dir);
353  else if(cell.subdetId()==EcalEndcap)
354  neighbours= EcalEndcapTopology_->getNeighbours(cell,dir);
355 
356  if ((!neighbours.empty()) && (!neighbours[0].null()))
357  {
358  cell = neighbours[0];
359  return true;
360  }
361  else
362  {
363  cell = DetId(0);
364  return false;
365  }
366 }
367 
369 {
370  bool result;
371  // One has to try both paths
372  if(dir==NORTHEAST)
373  {
374  result = simplemove(cell,NORTH);
375  if(result)
376  return simplemove(cell,EAST);
377  else
378  {
379  result = simplemove(cell,EAST);
380  if(result)
381  return simplemove(cell,NORTH);
382  else
383  return false;
384  }
385  }
386  else if(dir==NORTHWEST)
387  {
388  result = simplemove(cell,NORTH);
389  if(result)
390  return simplemove(cell,WEST);
391  else
392  {
393  result = simplemove(cell,WEST);
394  if(result)
395  return simplemove(cell,NORTH);
396  else
397  return false;
398  }
399  }
400  else if(dir == SOUTHEAST)
401  {
402  result = simplemove(cell,SOUTH);
403  if(result)
404  return simplemove(cell,EAST);
405  else
406  {
407  result = simplemove(cell,EAST);
408  if(result)
409  return simplemove(cell,SOUTH);
410  else
411  return false;
412  }
413  }
414  else if(dir == SOUTHWEST)
415  {
416  result = simplemove(cell,SOUTH);
417  if(result)
418  return simplemove(cell,WEST);
419  else
420  {
421  result = simplemove(cell,SOUTH);
422  if(result)
423  return simplemove(cell,WEST);
424  else
425  return false;
426  }
427  }
428  cell = DetId(0);
429  return false;
430 }
431 
432 bool CaloGeometryHelper::borderCrossing(const DetId& c1, const DetId& c2) const
433 {
434  if(c1.subdetId()!=c2.subdetId()) return false;
435 
436  if(c1.subdetId()==EcalBarrel)
437  {
438  // there is a crack if the two cells don't belong to the same
439  // module
440  EBDetId cc1(c1);
441  EBDetId cc2(c2);
442  return (cc1.im()!=cc2.im()||cc1.ism()!=cc2.ism() );
443  }
444 
445 if(c1.subdetId()==EcalEndcap)
446  {
447  // there is a crack if the two cells don't belong to the same
448  // module
449  return (EEDetId(c1).isc()!=EEDetId(c2).isc());
450  }
451  return false;
452 }
453 
455 {
456  const unsigned nbarrel = EBDetId::kSizeForDenseIndexing;
457  // Barrel first. The hashed index runs from 0 to 61199
458  barrelCrystals_.resize(nbarrel,BaseCrystal());
459 
460  //std::cout << " Building the array of crystals (barrel) " ;
461  const std::vector<DetId>& vec(EcalBarrelGeometry_->getValidDetIds(DetId::Ecal,EcalBarrel));
462  unsigned size=vec.size();
463  for(unsigned ic=0; ic<size; ++ic)
464  {
465  unsigned hashedindex=EBDetId(vec[ic]).hashedIndex();
466  auto geom = EcalBarrelGeometry_->getGeometry(vec[ic]);
467  BaseCrystal xtal(vec[ic]);
468  xtal.setCorners(geom->getCorners(),geom->getPosition());
469  barrelCrystals_[hashedindex]=xtal;
470  }
471 
472  // std::cout << " done " << size << std::endl;
473  // std::cout << " Building the array of crystals (endcap) " ;
474 
475 
476  const std::vector<DetId>& vece(EcalEndcapGeometry_->getValidDetIds(DetId::Ecal,EcalEndcap));
477  size=vece.size();
478  // There are some holes in the hashedIndex for the EE. Hence the array is bigger than the number
479  // of crystals
480  const unsigned nendcap=EEDetId::kSizeForDenseIndexing;
481 
482  endcapCrystals_.resize(nendcap,BaseCrystal());
483  for(unsigned ic=0; ic<size; ++ic)
484  {
485  unsigned hashedindex=EEDetId(vece[ic]).hashedIndex();
486  auto geom = EcalEndcapGeometry_->getGeometry(vece[ic]);
487  BaseCrystal xtal(vece[ic]);
488  xtal.setCorners(geom->getCorners(),geom->getPosition());
489  endcapCrystals_[hashedindex]=xtal;
490  }
491  // std::cout << " done " << size << std::endl;
492 }
#define LogDebug(id)
size
Write out results.
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:86
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
const EcalEndcapGeometry * EcalEndcapGeometry_
Definition: Calorimeter.h:78
int isc() const
Definition: EEDetId.cc:285
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:61
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
U second(std::pair< T, U > const &p)
static const CaloDirection orderedDir[8]
int im() const
get the number of module inside the SM (1-4)
Definition: EBDetId.h:66
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
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
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:56
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
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::vector< NeiVect > barrelNeighbours_
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:98
virtual DetId getClosestCell(const GlobalPoint &r) const
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
Definition: DetId.h:18
bool diagonalmove(DetId &cell, const CaloDirection &dir) const
int hashedIndex() const
Definition: EEDetId.h:182
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
bool null() const
is this a null id ?
Definition: DetId.h:45
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