CMS 3D CMS Logo

DDHGCalHEAlgo.cc
Go to the documentation of this file.
1 // File: DDHGCalHEAlgo.cc
3 // Description: Geometry factory class for HGCal (Mix)
5 
15 
16 //#define EDM_ML_DEBUG
17 using namespace geant_units::operators;
18 
20 #ifdef EDM_ML_DEBUG
21  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Creating an instance";
22 #endif
23 }
24 
26 
28  const DDVectorArguments& vArgs,
29  const DDMapArguments&,
30  const DDStringArguments& sArgs,
31  const DDStringVectorArguments& vsArgs) {
32  wafers_ = vsArgs["WaferNames"];
33 #ifdef EDM_ML_DEBUG
34  edm::LogVerbatim("HGCalGeom")
35  << "DDHGCalHEAlgo: " << wafers_.size() << " wafers";
36  for (unsigned int i = 0; i < wafers_.size(); ++i)
37  edm::LogVerbatim("HGCalGeom") << "Wafer[" << i << "] " << wafers_[i];
38 #endif
39  materials_ = vsArgs["MaterialNames"];
40  names_ = vsArgs["VolumeNames"];
41  thick_ = vArgs["Thickness"];
42  for (unsigned int i = 0; i < materials_.size(); ++i) {
43  copyNumber_.emplace_back(1);
44  }
45 #ifdef EDM_ML_DEBUG
46  edm::LogVerbatim("HGCalGeom")
47  << "DDHGCalHEAlgo: " << materials_.size() << " types of volumes";
48  for (unsigned int i = 0; i < names_.size(); ++i)
49  edm::LogVerbatim("HGCalGeom")
50  << "Volume [" << i << "] " << names_[i] << " of thickness " << thick_[i]
51  << " filled with " << materials_[i] << " first copy number "
52  << copyNumber_[i];
53 #endif
54  layers_ = dbl_to_int(vArgs["Layers"]);
55  layerThick_ = vArgs["LayerThick"];
56  rMixLayer_ = vArgs["LayerRmix"];
57 #ifdef EDM_ML_DEBUG
58  edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
59  for (unsigned int i = 0; i < layers_.size(); ++i)
60  edm::LogVerbatim("HGCalGeom")
61  << "Block [" << i << "] of thickness " << layerThick_[i] << " Rmid "
62  << rMixLayer_[i] << " with " << layers_[i] << " layers";
63 #endif
64  layerType_ = dbl_to_int(vArgs["LayerType"]);
65  layerSense_ = dbl_to_int(vArgs["LayerSense"]);
66  firstLayer_ = (int)(nArgs["FirstLayer"]);
67  absorbMode_ = (int)(nArgs["AbsorberMode"]);
68 #ifdef EDM_ML_DEBUG
69  edm::LogVerbatim("HGCalGeom") << "First Layere " << firstLayer_ << " and "
70  << "Absober mode " << absorbMode_;
71 #endif
72  if (firstLayer_ > 0) {
73  for (unsigned int i = 0; i < layerType_.size(); ++i) {
74  if (layerSense_[i] > 0) {
75  int ii = layerType_[i];
76  copyNumber_[ii] = firstLayer_;
77 #ifdef EDM_ML_DEBUG
78  edm::LogVerbatim("HGCalGeom")
79  << "First copy number for layer type " << i << ":" << ii << " with "
80  << materials_[ii] << " changed to " << copyNumber_[ii];
81 #endif
82  break;
83  }
84  }
85  }
86 #ifdef EDM_ML_DEBUG
87  edm::LogVerbatim("HGCalGeom")
88  << "There are " << layerType_.size() << " layers";
89  for (unsigned int i = 0; i < layerType_.size(); ++i)
90  edm::LogVerbatim("HGCalGeom")
91  << "Layer [" << i << "] with material type " << layerType_[i]
92  << " sensitive class " << layerSense_[i];
93 #endif
94  materialsTop_ = vsArgs["TopMaterialNames"];
95  namesTop_ = vsArgs["TopVolumeNames"];
96  layerThickTop_ = vArgs["TopLayerThickness"];
97  layerTypeTop_ = dbl_to_int(vArgs["TopLayerType"]);
98  for (unsigned int i = 0; i < materialsTop_.size(); ++i) {
99  copyNumberTop_.emplace_back(1);
100  }
101 #ifdef EDM_ML_DEBUG
102  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << materialsTop_.size()
103  << " types of volumes in the top part";
104  for (unsigned int i = 0; i < materialsTop_.size(); ++i)
105  edm::LogVerbatim("HGCalGeom")
106  << "Volume [" << i << "] " << namesTop_[i] << " of thickness "
107  << layerThickTop_[i] << " filled with " << materialsTop_[i]
108  << " first copy number " << copyNumberTop_[i];
109  edm::LogVerbatim("HGCalGeom")
110  << "There are " << layerTypeTop_.size() << " layers in the top part";
111  for (unsigned int i = 0; i < layerTypeTop_.size(); ++i)
112  edm::LogVerbatim("HGCalGeom")
113  << "Layer [" << i << "] with material type " << layerTypeTop_[i];
114 #endif
115  materialsBot_ = vsArgs["BottomMaterialNames"];
116  namesBot_ = vsArgs["BottomVolumeNames"];
117  layerTypeBot_ = dbl_to_int(vArgs["BottomLayerType"]);
118  layerSenseBot_ = dbl_to_int(vArgs["BottomLayerSense"]);
119  layerThickBot_ = vArgs["BottomLayerThickness"];
120  for (unsigned int i = 0; i < materialsBot_.size(); ++i) {
121  copyNumberBot_.emplace_back(1);
122  }
123 #ifdef EDM_ML_DEBUG
124  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << materialsBot_.size()
125  << " types of volumes in the bottom part";
126  for (unsigned int i = 0; i < materialsBot_.size(); ++i)
127  edm::LogVerbatim("HGCalGeom")
128  << "Volume [" << i << "] " << namesBot_[i] << " of thickness "
129  << layerThickBot_[i] << " filled with " << materialsBot_[i]
130  << " first copy number " << copyNumberBot_[i];
131  edm::LogVerbatim("HGCalGeom")
132  << "There are " << layerTypeBot_.size() << " layers in the top part";
133  for (unsigned int i = 0; i < layerTypeBot_.size(); ++i)
134  edm::LogVerbatim("HGCalGeom")
135  << "Layer [" << i << "] with material type " << layerTypeBot_[i]
136  << " sensitive class " << layerSenseBot_[i];
137 #endif
138  zMinBlock_ = nArgs["zMinBlock"];
139  rad100to200_ = vArgs["rad100to200"];
140  rad200to300_ = vArgs["rad200to300"];
141  zMinRadPar_ = nArgs["zMinForRadPar"];
142  choiceType_ = (int)(nArgs["choiceType"]);
143  nCutRadPar_ = (int)(nArgs["nCornerCut"]);
144  fracAreaMin_ = nArgs["fracAreaMin"];
145  waferSize_ = nArgs["waferSize"];
146  waferSepar_ = nArgs["SensorSeparation"];
147  sectors_ = (int)(nArgs["Sectors"]);
148  alpha_ = (1._pi) / sectors_;
149  cosAlpha_ = cos(alpha_);
150 #ifdef EDM_ML_DEBUG
151  edm::LogVerbatim("HGCalGeom")
152  << "DDHGCalHEAlgo: zStart " << zMinBlock_
153  << " radius for wafer type separation uses " << rad100to200_.size()
154  << " parameters; zmin " << zMinRadPar_ << " cutoff " << choiceType_ << ":"
155  << nCutRadPar_ << ":" << fracAreaMin_ << " wafer width " << waferSize_
156  << " separations " << waferSepar_ << " sectors " << sectors_ << ":"
157  << convertRadToDeg(alpha_) << ":" << cosAlpha_;
158  for (unsigned int k = 0; k < rad100to200_.size(); ++k)
159  edm::LogVerbatim("HGCalGeom") << "[" << k << "] 100-200 " << rad100to200_[k]
160  << " 200-300 " << rad200to300_[k];
161 #endif
162  slopeB_ = vArgs["SlopeBottom"];
163  zFrontB_ = vArgs["ZFrontBottom"];
164  rMinFront_ = vArgs["RMinFront"];
165  slopeT_ = vArgs["SlopeTop"];
166  zFrontT_ = vArgs["ZFrontTop"];
167  rMaxFront_ = vArgs["RMaxFront"];
168 #ifdef EDM_ML_DEBUG
169  for (unsigned int i = 0; i < slopeB_.size(); ++i)
170  edm::LogVerbatim("HGCalGeom")
171  << "Block [" << i << "] Zmin " << zFrontB_[i] << " Rmin "
172  << rMinFront_[i] << " Slope " << slopeB_[i];
173  for (unsigned int i = 0; i < slopeT_.size(); ++i)
174  edm::LogVerbatim("HGCalGeom")
175  << "Block [" << i << "] Zmin " << zFrontT_[i] << " Rmax "
176  << rMaxFront_[i] << " Slope " << slopeT_[i];
177 #endif
178  nameSpace_ = DDCurrentNamespace::ns();
179 #ifdef EDM_ML_DEBUG
180  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: NameSpace " << nameSpace_;
181 #endif
182 
183  waferType_ = std::make_unique<HGCalWaferType>(
184  rad100to200_, rad200to300_, (waferSize_ + waferSepar_), zMinRadPar_,
185  choiceType_, nCutRadPar_, fracAreaMin_);
186 }
187 
189 // DDHGCalHEAlgo methods...
191 
193 #ifdef EDM_ML_DEBUG
194  edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalHEAlgo...";
195  copies_.clear();
196 #endif
197  constructLayers(parent(), cpv);
198 #ifdef EDM_ML_DEBUG
199  edm::LogVerbatim("HGCalGeom")
200  << "DDHGCalHEAlgo: " << copies_.size() << " different wafer copy numbers";
201  int k(0);
202  for (std::unordered_set<int>::const_iterator itr = copies_.begin();
203  itr != copies_.end(); ++itr, ++k) {
204  edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
205  }
206  copies_.clear();
207  edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalHEAlgo construction...";
208 #endif
209 }
210 
212  DDCompactView& cpv) {
213 #ifdef EDM_ML_DEBUG
214  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: \t\tInside Layers";
215 #endif
216  double zi(zMinBlock_);
217  int laymin(0);
218  const double tol(0.01);
219  for (unsigned int i = 0; i < layers_.size(); i++) {
220  double zo = zi + layerThick_[i];
221  double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
222  int laymax = laymin + layers_[i];
223  double zz = zi;
224  double thickTot(0);
225  for (int ly = laymin; ly < laymax; ++ly) {
226  int ii = layerType_[ly];
227  int copy = copyNumber_[ii];
228  double hthick = 0.5 * thick_[ii];
229  double rinB = HGCalGeomTools::radius(zo, zFrontB_, rMinFront_, slopeB_);
230  zz += hthick;
231  thickTot += thick_[ii];
232 
233  std::string name = names_[ii] + std::to_string(copy);
234 #ifdef EDM_ML_DEBUG
235  edm::LogVerbatim("HGCalGeom")
236  << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " Front " << zi
237  << ", " << routF << " Back " << zo << ", " << rinB
238  << " superlayer thickness " << layerThick_[i];
239 #endif
240  DDName matName(DDSplit(materials_[ii]).first,
241  DDSplit(materials_[ii]).second);
242  DDMaterial matter(matName);
243  DDLogicalPart glog;
244  if (layerSense_[ly] < 1) {
245  std::vector<double> pgonZ, pgonRin, pgonRout;
246  if (layerSense_[ly] == 0 || absorbMode_ == 0) {
247  double rmax =
248  (std::min(routF, HGCalGeomTools::radius(zz + hthick, zFrontT_,
249  rMaxFront_, slopeT_)) *
250  cosAlpha_) -
251  tol;
252  pgonZ.emplace_back(-hthick);
253  pgonZ.emplace_back(hthick);
254  pgonRin.emplace_back(rinB);
255  pgonRin.emplace_back(rinB);
256  pgonRout.emplace_back(rmax);
257  pgonRout.emplace_back(rmax);
258  } else {
259  HGCalGeomTools::radius(zz-hthick, zz+hthick, zFrontB_, rMinFront_,
260  slopeB_, zFrontT_, rMaxFront_, slopeT_,
261  -layerSense_[ly], pgonZ, pgonRin, pgonRout);
262  for (unsigned int isec=0; isec < pgonZ.size(); ++isec) {
263  pgonZ[isec] -= zz;
264  pgonRout[isec] = pgonRout[isec]*cosAlpha_ - tol;
265  }
266  }
267  DDSolid solid =
268  DDSolidFactory::polyhedra(DDName(name, nameSpace_), sectors_,
269  -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
270  glog = DDLogicalPart(solid.ddname(), matter, solid);
271 #ifdef EDM_ML_DEBUG
272  edm::LogVerbatim("HGCalGeom")
273  << "DDHGCalHEAlgo: " << solid.name() << " polyhedra of " << sectors_
274  << " sectors covering " << convertRadToDeg(-alpha_) << ":"
275  << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size()
276  << " sections";
277  for (unsigned int k = 0; k < pgonZ.size(); ++k)
278  edm::LogVerbatim("HGCalGeom")
279  << "[" << k << "] z " << pgonZ[k] << " R " << pgonRin[k] << ":"
280  << pgonRout[k];
281 #endif
282  } else {
283  DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick,
284  rinB, routF, 0.0, 2._pi);
285  glog = DDLogicalPart(solid.ddname(), matter, solid);
286 #ifdef EDM_ML_DEBUG
287  edm::LogVerbatim("HGCalGeom")
288  << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
289  << " of dimensions " << rinB << ", " << routF << ", " << hthick
290  << ", 0.0, 360.0 and positioned in: " << glog.name() << " number "
291  << copy;
292 #endif
293  positionMix(glog, name, copy, thick_[ii], matter, rinB, rMixLayer_[i],
294  routF, zz, cpv);
295  }
296  DDTranslation r1(0, 0, zz);
297  DDRotation rot;
298  cpv.position(glog, module, copy, r1, rot);
299  ++copyNumber_[ii];
300 #ifdef EDM_ML_DEBUG
301  edm::LogVerbatim("HGCalGeom")
302  << "DDHGCalHEAlgo: " << glog.name() << " number " << copy
303  << " positioned in " << module.name() << " at " << r1 << " with "
304  << rot;
305 #endif
306  zz += hthick;
307  } // End of loop over layers in a block
308  zi = zo;
309  laymin = laymax;
310  if (std::abs(thickTot - layerThick_[i]) < 0.00001) {
311  } else if (thickTot > layerThick_[i]) {
312  edm::LogError("HGCalGeom")
313  << "Thickness of the partition " << layerThick_[i]
314  << " is smaller than " << thickTot << ": thickness of all its "
315  << "components **** ERROR ****";
316  } else if (thickTot < layerThick_[i]) {
317  edm::LogWarning("HGCalGeom")
318  << "Thickness of the partition " << layerThick_[i]
319  << " does not match with " << thickTot << " of the components";
320  }
321  } // End of loop over blocks
322 }
323 
325  const std::string& nameM, int copyM,
326  double thick, const DDMaterial& matter,
327  double rin, double rmid, double rout, double zz,
328  DDCompactView& cpv) {
329  DDLogicalPart glog1;
330  DDTranslation tran;
331  DDRotation rot;
332  for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
333  int ii = layerTypeTop_[ly];
334  copyNumberTop_[ii] = copyM;
335  }
336  for (unsigned int ly = 0; ly < layerTypeBot_.size(); ++ly) {
337  int ii = layerTypeBot_[ly];
338  copyNumberBot_[ii] = copyM;
339  }
340  double hthick = 0.5 * thick;
341  // Make the top part first
342  std::string name = nameM + "Top";
343  DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rmid,
344  rout, 0.0, 2._pi);
345  glog1 = DDLogicalPart(solid.ddname(), matter, solid);
346 #ifdef EDM_ML_DEBUG
347  edm::LogVerbatim("HGCalGeom")
348  << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matter.name()
349  << " of dimensions " << rmid << ", " << rout << ", " << hthick
350  << ", 0.0, 360.0";
351 #endif
352  cpv.position(glog1, glog, 1, tran, rot);
353 #ifdef EDM_ML_DEBUG
354  edm::LogVerbatim("HGCalGeom")
355  << "DDHGCalHEAlgo: " << glog1.name() << " number 1 positioned in "
356  << glog.name() << " at " << tran << " with " << rot;
357 #endif
358  double thickTot(0), zpos(-hthick);
359  for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
360  int ii = layerTypeTop_[ly];
361  int copy = copyNumberTop_[ii];
362  double hthickl = 0.5 * layerThickTop_[ii];
363  thickTot += layerThickTop_[ii];
364  name = namesTop_[ii] + std::to_string(copy);
365 #ifdef EDM_ML_DEBUG
366  edm::LogVerbatim("HGCalGeom")
367  << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " R " << rmid << ":"
368  << rout << " Thick " << layerThickTop_[ii];
369 #endif
370  DDName matName(DDSplit(materialsTop_[ii]).first,
371  DDSplit(materialsTop_[ii]).second);
372  DDMaterial matter1(matName);
373  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthickl, rmid, rout,
374  0.0, 2._pi);
375  DDLogicalPart glog2 = DDLogicalPart(solid.ddname(), matter1, solid);
376 #ifdef EDM_ML_DEBUG
377  double eta1 = -log(tan(0.5 * atan(rmid / zz)));
378  double eta2 = -log(tan(0.5 * atan(rout / zz)));
379  edm::LogVerbatim("HGCalGeom")
380  << name << " z|rin|rout " << zz << ":" << rmid << ":" << rout << " eta "
381  << eta1 << ":" << eta2;
382  edm::LogVerbatim("HGCalGeom")
383  << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
384  << " of dimensions " << rmid << ", " << rout << ", " << hthickl
385  << ", 0.0, 360.0";
386 #endif
387  zpos += hthickl;
388  DDTranslation r1(0, 0, zpos);
389  cpv.position(glog2, glog1, copy, r1, rot);
390 #ifdef EDM_ML_DEBUG
391  edm::LogVerbatim("HGCalGeom")
392  << "DDHGCalHEAlgo: Position " << glog2.name() << " number " << copy
393  << " in " << glog1.name() << " at " << r1 << " with " << rot;
394 #endif
395  ++copyNumberTop_[ii];
396  zpos += hthickl;
397  }
398  if (std::abs(thickTot - thick) < 0.00001) {
399  } else if (thickTot > thick) {
400  edm::LogError("HGCalGeom")
401  << "Thickness of the partition " << thick << " is smaller than "
402  << thickTot << ": thickness of all its components in "
403  << "the top part **** ERROR ****";
404  } else if (thickTot < thick) {
405  edm::LogWarning("HGCalGeom")
406  << "Thickness of the partition " << thick << " does not match with "
407  << thickTot << " of the components in top part";
408  }
409 
410  // Make the bottom part next
411  name = nameM + "Bottom";
412  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rin, rmid, 0.0,
413  2._pi);
414  glog1 = DDLogicalPart(solid.ddname(), matter, solid);
415 #ifdef EDM_ML_DEBUG
416  edm::LogVerbatim("HGCalGeom")
417  << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matter.name()
418  << " of dimensions " << rin << ", " << rmid << ", " << hthick
419  << ", 0.0, 360.0";
420 #endif
421  cpv.position(glog1, glog, 1, tran, rot);
422 #ifdef EDM_ML_DEBUG
423  edm::LogVerbatim("HGCalGeom")
424  << "DDHGCalHEAlgo: " << glog1.name() << " number 1 positioned in "
425  << glog.name() << " at " << tran << " with " << rot;
426 #endif
427  thickTot = 0;
428  zpos = -hthick;
429  for (unsigned int ly = 0; ly < layerTypeBot_.size(); ++ly) {
430  int ii = layerTypeBot_[ly];
431  int copy = copyNumberBot_[ii];
432  double hthickl = 0.5 * layerThickBot_[ii];
433  thickTot += layerThickBot_[ii];
434  name = namesBot_[ii] + std::to_string(copy);
435 #ifdef EDM_ML_DEBUG
436  edm::LogVerbatim("HGCalGeom")
437  << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " R " << rin << ":"
438  << rmid << " Thick " << layerThickBot_[ii];
439 #endif
440  DDName matName(DDSplit(materialsBot_[ii]).first,
441  DDSplit(materialsBot_[ii]).second);
442  DDMaterial matter1(matName);
443  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthickl, rin, rmid,
444  0.0, 2._pi);
445  DDLogicalPart glog2 = DDLogicalPart(solid.ddname(), matter1, solid);
446 #ifdef EDM_ML_DEBUG
447  double eta1 = -log(tan(0.5 * atan(rin / zz)));
448  double eta2 = -log(tan(0.5 * atan(rmid / zz)));
449  edm::LogVerbatim("HGCalGeom")
450  << name << " z|rin|rout " << zz << ":" << rin << ":" << rmid << " eta "
451  << eta1 << ":" << eta2;
452  edm::LogVerbatim("HGCalGeom")
453  << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
454  << " of dimensions " << rin << ", " << rmid << ", " << hthickl
455  << ", 0.0, 360.0";
456 #endif
457  zpos += hthickl;
458  DDTranslation r1(0, 0, zpos);
459  cpv.position(glog2, glog1, copy, r1, rot);
460 #ifdef EDM_ML_DEBUG
461  edm::LogVerbatim("HGCalGeom")
462  << "DDHGCalHEAlgo: Position " << glog2.name() << " number " << copy
463  << " in " << glog1.name() << " at " << r1 << " with " << rot;
464 #endif
465  if (layerSenseBot_[ly] != 0)
466  positionSensitive(glog2, rin, rmid, zz + zpos, layerSenseBot_[ly], cpv);
467  zpos += hthickl;
468  ++copyNumberBot_[ii];
469  }
470  if (std::abs(thickTot - thick) < 0.00001) {
471  } else if (thickTot > thick) {
472  edm::LogError("HGCalGeom")
473  << "Thickness of the partition " << thick << " is smaller than "
474  << thickTot << ": thickness of all its components in "
475  << "the top part **** ERROR ****";
476  } else if (thickTot < thick) {
477  edm::LogWarning("HGCalGeom")
478  << "Thickness of the partition " << thick << " does not match with "
479  << thickTot << " of the components in top part";
480  }
481 }
482 
483 void DDHGCalHEAlgo::positionSensitive(const DDLogicalPart& glog, double rin,
484  double rout, double zpos, int layertype,
485  DDCompactView& cpv) {
486  static const double sqrt3 = std::sqrt(3.0);
487  double r = 0.5 * (waferSize_ + waferSepar_);
488  double R = 2.0 * r / sqrt3;
489  double dy = 0.75 * R;
490  int N = (int)(0.5 * rout / r) + 2;
491 #ifdef EDM_ML_DEBUG
492  int ium(0), ivm(0), iumAll(0), ivmAll(0), kount(0), ntot(0), nin(0);
493  std::vector<int> ntype(6, 0);
494  edm::LogVerbatim("HGCalGeom")
495  << "DDHGCalHEAlgo: " << glog.ddname() << " rout " << rout << " N " << N
496  << " for maximum u, v";
497 #endif
498  for (int u = -N; u <= N; ++u) {
499  int iu = std::abs(u);
500  for (int v = -N; v <= N; ++v) {
501  int iv = std::abs(v);
502  int nr = 2 * v;
503  int nc = -2 * u + v;
504  double xpos = nc * r;
505  double ypos = nr * dy;
506  std::pair<int, int> corner =
507  HGCalGeomTools::waferCorner(xpos, ypos, r, R, rin, rout, false);
508 #ifdef EDM_ML_DEBUG
509  ++ntot;
510 #endif
511  if (corner.first > 0) {
512  int type = waferType_->getType(xpos, ypos, zpos);
513  int copy = type * 1000000 + iv * 100 + iu;
514  if (u < 0) copy += 10000;
515  if (v < 0) copy += 100000;
516 #ifdef EDM_ML_DEBUG
517  if (iu > ium) ium = iu;
518  if (iv > ivm) ivm = iv;
519  kount++;
520  if (copies_.count(copy) == 0) copies_.insert(copy);
521 #endif
522  if (corner.first == (int)(HGCalParameters::k_CornerSize)) {
523 #ifdef EDM_ML_DEBUG
524  if (iu > iumAll) iumAll = iu;
525  if (iv > ivmAll) ivmAll = iv;
526  ++nin;
527 #endif
528  DDTranslation tran(xpos, ypos, 0.0);
530  if (layertype > 1) type += 3;
531  DDName name = DDName(DDSplit(wafers_[type]).first,
532  DDSplit(wafers_[type]).second);
533  cpv.position(name, glog.ddname(), copy, tran, rotation);
534 #ifdef EDM_ML_DEBUG
535  ++ntype[type];
536  edm::LogVerbatim("HGCalGeom")
537  << "DDHGCalHEAlgo: " << name << " number " << copy
538  << " positioned in " << glog.ddname() << " at " << tran
539  << " with " << rotation;
540 #endif
541  }
542  }
543  }
544  }
545 #ifdef EDM_ML_DEBUG
546  edm::LogVerbatim("HGCalGeom")
547  << "DDHGCalHEAlgo: Maximum # of u " << ium << ":" << iumAll << " # of v "
548  << ivm << ":" << ivmAll << " and " << nin << ":" << kount << ":" << ntot
549  << " wafers (" << ntype[0] << ":" << ntype[1] << ":" << ntype[2] << ":"
550  << ntype[3] << ":" << ntype[4] << ":" << ntype[5] << ") for "
551  << glog.ddname() << " R " << rin << ":" << rout;
552 #endif
553 }
type
Definition: HCALResponse.h:21
static void radius(double zf, double zb, std::vector< double > const &zFront1, std::vector< double > const &rFront1, std::vector< double > const &slope1, std::vector< double > const &zFront2, std::vector< double > const &rFront2, std::vector< double > const &slope2, int flag, std::vector< double > &zz, std::vector< double > &rin, std::vector< double > &rout)
const N & name() const
Definition: DDBase.h:74
void positionMix(const DDLogicalPart &glog, const std::string &name, int copy, double thick, const DDMaterial &matter, double rin, double rmid, double routF, double zz, DDCompactView &cpv)
def copy(args, dbName)
DDMaterial is used to define and access material information.
Definition: DDMaterial.h:43
void execute(DDCompactView &cpv) override
DDName is used to identify DDD entities uniquely.
Definition: DDName.h:15
constexpr NumType convertRadToDeg(NumType radians)
Definition: GeantUnits.h:98
void initialize(const DDNumericArguments &nArgs, const DDVectorArguments &vArgs, const DDMapArguments &mArgs, const DDStringArguments &sArgs, const DDStringVectorArguments &vsArgs) override
static std::string & ns()
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
~DDHGCalHEAlgo() override
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DDTranslation
Definition: DDTranslation.h:7
int nin
Represents a uniquely identifyable rotation matrix.
Definition: DDTransform.h:68
U second(std::pair< T, U > const &p)
std::vector< int > dbl_to_int(const std::vector< double > &vecdbl)
Converts a std::vector of doubles to a std::vector of int.
Definition: DDutils.h:7
static std::pair< int32_t, int32_t > waferCorner(double xpos, double ypos, double r, double R, double rMin, double rMax, bool oldBug=false)
T sqrt(T t)
Definition: SSEVec.h:18
void constructLayers(const DDLogicalPart &, DDCompactView &cpv)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:93
static DDSolid tubs(const DDName &name, double zhalf, double rIn, double rOut, double startPhi, double deltaPhi)
Definition: DDSolid.cc:865
T min(T a, T b)
Definition: MathUtil.h:58
ii
Definition: cuy.py:590
int k[5][pyjets_maxn]
static uint32_t k_CornerSize
#define N
Definition: blowfish.cc:9
void position(const DDLogicalPart &self, const DDLogicalPart &parent, const std::string &copyno, const DDTranslation &trans, const DDRotation &rot, const DDDivision *div=0)
void positionSensitive(const DDLogicalPart &glog, double rin, double rout, double zpos, int layertype, DDCompactView &cpv)
std::pair< std::string, std::string > DDSplit(const std::string &n)
split into (name,namespace), separator = &#39;:&#39;
Definition: DDSplit.cc:3
Definition: vlib.h:208
static DDSolid polyhedra(const DDName &name, int sides, double startPhi, double deltaPhi, const std::vector< double > &z, const std::vector< double > &rmin, const std::vector< double > &rmax)
Creates a polyhedra (refere to Geant3 or Geant4 documentation)
Definition: DDSolid.cc:730
const N & ddname() const
Definition: DDBase.h:76