CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DDHGCalHEAlgo.cc
Go to the documentation of this file.
1 // File: DDHGCalHEAlgo.cc
3 // Description: Geometry factory class for HGCal (Mix)
5 
22 
23 #include <cmath>
24 #include <memory>
25 #include <string>
26 #include <unordered_set>
27 #include <vector>
28 
29 //#define EDM_ML_DEBUG
30 using namespace angle_units::operators;
31 
32 class DDHGCalHEAlgo : public DDAlgorithm {
33 public:
34  // Constructor and Destructor
35  DDHGCalHEAlgo(); // const std::string & name);
36  ~DDHGCalHEAlgo() override;
37 
38  void initialize(const DDNumericArguments& nArgs,
39  const DDVectorArguments& vArgs,
40  const DDMapArguments& mArgs,
41  const DDStringArguments& sArgs,
42  const DDStringVectorArguments& vsArgs) override;
43  void execute(DDCompactView& cpv) override;
44 
45 protected:
46  void constructLayers(const DDLogicalPart&, DDCompactView& cpv);
47  void positionMix(const DDLogicalPart& glog,
48  const std::string& name,
49  int copy,
50  double thick,
51  const DDMaterial& matter,
52  double rin,
53  double rmid,
54  double routF,
55  double zz,
56  DDCompactView& cpv);
57  void positionSensitive(const DDLogicalPart& glog,
58  double rin,
59  double rout,
60  double zpos,
61  int layertype,
62  int layercenter,
63  DDCompactView& cpv);
64 
65 private:
67  std::unique_ptr<HGCalWaferType> waferType_;
68 
69  static constexpr double tol1_ = 0.01;
70  static constexpr double tol2_ = 0.00001;
71 
72  std::vector<std::string> wafers_; // Wafers
73  std::vector<std::string> materials_; // Materials
74  std::vector<std::string> names_; // Names
75  std::vector<double> thick_; // Thickness of the material
76  std::vector<int> copyNumber_; // Initial copy numbers
77  std::vector<int> layers_; // Number of layers in a section
78  std::vector<double> layerThick_; // Thickness of each section
79  std::vector<double> rMixLayer_; // Partition between Si/Sci part
80  std::vector<int> layerType_; // Type of the layer
81  std::vector<int> layerSense_; // Content of a layer (sensitive?)
82  int firstLayer_; // Copy # of the first sensitive layer
83  int absorbMode_; // Absorber mode
84  int sensitiveMode_; // Sensitive mode
85  std::vector<std::string> materialsTop_; // Materials of top layers
86  std::vector<std::string> namesTop_; // Names of top layers
87  std::vector<double> layerThickTop_; // Thickness of the top sections
88  std::vector<int> layerTypeTop_; // Type of the Top layer
89  std::vector<int> copyNumberTop_; // Initial copy numbers (top section)
90  std::vector<std::string> materialsBot_; // Materials of bottom layers
91  std::vector<std::string> namesBot_; // Names of bottom layers
92  std::vector<double> layerThickBot_; // Thickness of the bottom sections
93  std::vector<int> layerTypeBot_; // Type of the bottom layers
94  std::vector<int> copyNumberBot_; // Initial copy numbers (bot section)
95  std::vector<int> layerSenseBot_; // Content of bottom layer (sensitive?)
96  std::vector<int> layerCenter_; // Centering of the wafers
97 
98  double zMinBlock_; // Starting z-value of the block
99  std::vector<double> rad100to200_; // Parameters for 120-200mum trans.
100  std::vector<double> rad200to300_; // Parameters for 200-300mum trans.
101  double zMinRadPar_; // Minimum z for radius parametriz.
102  int choiceType_; // Type of parametrization to be used
103  int nCutRadPar_; // Cut off threshold for corners
104  double fracAreaMin_; // Minimum fractional conatined area
105  double waferSize_; // Width of the wafer
106  double waferSepar_; // Sensor separation
107  int sectors_; // Sectors
108  std::vector<double> slopeB_; // Slope at the lower R
109  std::vector<double> zFrontB_; // Starting Z values for the slopes
110  std::vector<double> rMinFront_; // Corresponding rMin's
111  std::vector<double> slopeT_; // Slopes at the larger R
112  std::vector<double> zFrontT_; // Starting Z values for the slopes
113  std::vector<double> rMaxFront_; // Corresponding rMax's
114  std::string nameSpace_; // Namespace of this and ALL sub-parts
115  std::unordered_set<int> copies_; // List of copy #'s
116  double alpha_, cosAlpha_;
117 };
118 
120 #ifdef EDM_ML_DEBUG
121  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Creating an instance";
122 #endif
123 }
124 
126 
128  const DDVectorArguments& vArgs,
129  const DDMapArguments&,
130  const DDStringArguments& sArgs,
131  const DDStringVectorArguments& vsArgs) {
132  wafers_ = vsArgs["WaferNames"];
133 #ifdef EDM_ML_DEBUG
134  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << wafers_.size() << " wafers";
135  for (unsigned int i = 0; i < wafers_.size(); ++i)
136  edm::LogVerbatim("HGCalGeom") << "Wafer[" << i << "] " << wafers_[i];
137 #endif
138  materials_ = vsArgs["MaterialNames"];
139  names_ = vsArgs["VolumeNames"];
140  thick_ = vArgs["Thickness"];
141  copyNumber_.resize(materials_.size(), 1);
142 #ifdef EDM_ML_DEBUG
143  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << materials_.size() << " types of volumes";
144  for (unsigned int i = 0; i < names_.size(); ++i)
145  edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << names_[i] << " of thickness " << thick_[i]
146  << " filled with " << materials_[i] << " first copy number " << copyNumber_[i];
147 #endif
148  layers_ = dbl_to_int(vArgs["Layers"]);
149  layerThick_ = vArgs["LayerThick"];
150  rMixLayer_ = vArgs["LayerRmix"];
151 #ifdef EDM_ML_DEBUG
152  edm::LogVerbatim("HGCalGeom") << "There are " << layers_.size() << " blocks";
153  for (unsigned int i = 0; i < layers_.size(); ++i)
154  edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] of thickness " << layerThick_[i] << " Rmid " << rMixLayer_[i]
155  << " with " << layers_[i] << " layers";
156 #endif
157  layerType_ = dbl_to_int(vArgs["LayerType"]);
158  layerSense_ = dbl_to_int(vArgs["LayerSense"]);
159  firstLayer_ = (int)(nArgs["FirstLayer"]);
160  absorbMode_ = (int)(nArgs["AbsorberMode"]);
161  sensitiveMode_ = (int)(nArgs["SensitiveMode"]);
162 #ifdef EDM_ML_DEBUG
163  edm::LogVerbatim("HGCalGeom") << "First Layer " << firstLayer_ << " and "
164  << "Absober:Sensitive mode " << absorbMode_ << ":" << sensitiveMode_;
165 #endif
166  layerCenter_ = dbl_to_int(vArgs["LayerCenter"]);
167 #ifdef EDM_ML_DEBUG
168  for (unsigned int i = 0; i < layerCenter_.size(); ++i)
169  edm::LogVerbatim("HGCalGeom") << "LayerCenter [" << i << "] " << layerCenter_[i];
170 #endif
171  if (firstLayer_ > 0) {
172  for (unsigned int i = 0; i < layerType_.size(); ++i) {
173  if (layerSense_[i] > 0) {
174  int ii = layerType_[i];
175  copyNumber_[ii] = firstLayer_;
176 #ifdef EDM_ML_DEBUG
177  edm::LogVerbatim("HGCalGeom") << "First copy number for layer type " << i << ":" << ii << " with "
178  << materials_[ii] << " changed to " << copyNumber_[ii];
179 #endif
180  break;
181  }
182  }
183  }
184 #ifdef EDM_ML_DEBUG
185  edm::LogVerbatim("HGCalGeom") << "There are " << layerType_.size() << " layers";
186  for (unsigned int i = 0; i < layerType_.size(); ++i)
187  edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerType_[i] << " sensitive class "
188  << layerSense_[i];
189 #endif
190  materialsTop_ = vsArgs["TopMaterialNames"];
191  namesTop_ = vsArgs["TopVolumeNames"];
192  layerThickTop_ = vArgs["TopLayerThickness"];
193  layerTypeTop_ = dbl_to_int(vArgs["TopLayerType"]);
194  copyNumberTop_.resize(materialsTop_.size(), 1);
195 #ifdef EDM_ML_DEBUG
196  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << materialsTop_.size() << " types of volumes in the top part";
197  for (unsigned int i = 0; i < materialsTop_.size(); ++i)
198  edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << namesTop_[i] << " of thickness " << layerThickTop_[i]
199  << " filled with " << materialsTop_[i] << " first copy number " << copyNumberTop_[i];
200  edm::LogVerbatim("HGCalGeom") << "There are " << layerTypeTop_.size() << " layers in the top part";
201  for (unsigned int i = 0; i < layerTypeTop_.size(); ++i)
202  edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerTypeTop_[i];
203 #endif
204  materialsBot_ = vsArgs["BottomMaterialNames"];
205  namesBot_ = vsArgs["BottomVolumeNames"];
206  layerTypeBot_ = dbl_to_int(vArgs["BottomLayerType"]);
207  layerSenseBot_ = dbl_to_int(vArgs["BottomLayerSense"]);
208  layerThickBot_ = vArgs["BottomLayerThickness"];
209  copyNumberBot_.resize(materialsBot_.size(), 1);
210 #ifdef EDM_ML_DEBUG
211  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << materialsBot_.size() << " types of volumes in the bottom part";
212  for (unsigned int i = 0; i < materialsBot_.size(); ++i)
213  edm::LogVerbatim("HGCalGeom") << "Volume [" << i << "] " << namesBot_[i] << " of thickness " << layerThickBot_[i]
214  << " filled with " << materialsBot_[i] << " first copy number " << copyNumberBot_[i];
215  edm::LogVerbatim("HGCalGeom") << "There are " << layerTypeBot_.size() << " layers in the bottom part";
216  for (unsigned int i = 0; i < layerTypeBot_.size(); ++i)
217  edm::LogVerbatim("HGCalGeom") << "Layer [" << i << "] with material type " << layerTypeBot_[i]
218  << " sensitive class " << layerSenseBot_[i];
219 #endif
220  zMinBlock_ = nArgs["zMinBlock"];
221  rad100to200_ = vArgs["rad100to200"];
222  rad200to300_ = vArgs["rad200to300"];
223  zMinRadPar_ = nArgs["zMinForRadPar"];
224  choiceType_ = (int)(nArgs["choiceType"]);
225  nCutRadPar_ = (int)(nArgs["nCornerCut"]);
226  fracAreaMin_ = nArgs["fracAreaMin"];
227  waferSize_ = nArgs["waferSize"];
228  waferSepar_ = nArgs["SensorSeparation"];
229  sectors_ = (int)(nArgs["Sectors"]);
230  alpha_ = (1._pi) / sectors_;
231  cosAlpha_ = cos(alpha_);
232 #ifdef EDM_ML_DEBUG
233  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: zStart " << zMinBlock_ << " radius for wafer type separation uses "
234  << rad100to200_.size() << " parameters; zmin " << zMinRadPar_ << " cutoff "
235  << choiceType_ << ":" << nCutRadPar_ << ":" << fracAreaMin_ << " wafer width "
236  << waferSize_ << " separations " << waferSepar_ << " sectors " << sectors_ << ":"
237  << convertRadToDeg(alpha_) << ":" << cosAlpha_;
238  for (unsigned int k = 0; k < rad100to200_.size(); ++k)
239  edm::LogVerbatim("HGCalGeom") << "[" << k << "] 100-200 " << rad100to200_[k] << " 200-300 " << rad200to300_[k];
240 #endif
241  slopeB_ = vArgs["SlopeBottom"];
242  zFrontB_ = vArgs["ZFrontBottom"];
243  rMinFront_ = vArgs["RMinFront"];
244  slopeT_ = vArgs["SlopeTop"];
245  zFrontT_ = vArgs["ZFrontTop"];
246  rMaxFront_ = vArgs["RMaxFront"];
247 #ifdef EDM_ML_DEBUG
248  for (unsigned int i = 0; i < slopeB_.size(); ++i)
249  edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] Zmin " << zFrontB_[i] << " Rmin " << rMinFront_[i]
250  << " Slope " << slopeB_[i];
251  for (unsigned int i = 0; i < slopeT_.size(); ++i)
252  edm::LogVerbatim("HGCalGeom") << "Block [" << i << "] Zmin " << zFrontT_[i] << " Rmax " << rMaxFront_[i]
253  << " Slope " << slopeT_[i];
254 #endif
255  nameSpace_ = DDCurrentNamespace::ns();
256 #ifdef EDM_ML_DEBUG
257  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: NameSpace " << nameSpace_;
258 #endif
259 
260  waferType_ = std::make_unique<HGCalWaferType>(
261  rad100to200_, rad200to300_, (waferSize_ + waferSepar_), zMinRadPar_, choiceType_, nCutRadPar_, fracAreaMin_);
262 }
263 
265 // DDHGCalHEAlgo methods...
267 
269 #ifdef EDM_ML_DEBUG
270  edm::LogVerbatim("HGCalGeom") << "==>> Constructing DDHGCalHEAlgo...";
271  copies_.clear();
272 #endif
273  constructLayers(parent(), cpv);
274 #ifdef EDM_ML_DEBUG
275  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << copies_.size() << " different wafer copy numbers";
276  int k(0);
277  for (std::unordered_set<int>::const_iterator itr = copies_.begin(); itr != copies_.end(); ++itr, ++k) {
278  edm::LogVerbatim("HGCalGeom") << "Copy [" << k << "] : " << (*itr);
279  }
280  copies_.clear();
281  edm::LogVerbatim("HGCalGeom") << "<<== End of DDHGCalHEAlgo construction...";
282 #endif
283 }
284 
286 #ifdef EDM_ML_DEBUG
287  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: \t\tInside Layers";
288 #endif
289  double zi(zMinBlock_);
290  int laymin(0);
291  for (unsigned int i = 0; i < layers_.size(); i++) {
292  double zo = zi + layerThick_[i];
293  double routF = HGCalGeomTools::radius(zi, zFrontT_, rMaxFront_, slopeT_);
294  int laymax = laymin + layers_[i];
295  double zz = zi;
296  double thickTot(0);
297  for (int ly = laymin; ly < laymax; ++ly) {
298  int ii = layerType_[ly];
299  int copy = copyNumber_[ii];
300  double hthick = 0.5 * thick_[ii];
301  double rinB = HGCalGeomTools::radius(zo, zFrontB_, rMinFront_, slopeB_);
302  zz += hthick;
303  thickTot += thick_[ii];
304 
305  std::string name = names_[ii] + std::to_string(copy);
306 #ifdef EDM_ML_DEBUG
307  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " Front " << zi << ", " << routF
308  << " Back " << zo << ", " << rinB << " superlayer thickness " << layerThick_[i];
309 #endif
310  DDName matName(DDSplit(materials_[ii]).first, DDSplit(materials_[ii]).second);
311  DDMaterial matter(matName);
312  DDLogicalPart glog;
313  if (layerSense_[ly] < 1) {
314  std::vector<double> pgonZ, pgonRin, pgonRout;
315  if (layerSense_[ly] == 0 || absorbMode_ == 0) {
316  double rmax =
317  (std::min(routF, HGCalGeomTools::radius(zz + hthick, zFrontT_, rMaxFront_, slopeT_)) * cosAlpha_) - tol1_;
318  pgonZ.emplace_back(-hthick);
319  pgonZ.emplace_back(hthick);
320  pgonRin.emplace_back(rinB);
321  pgonRin.emplace_back(rinB);
322  pgonRout.emplace_back(rmax);
323  pgonRout.emplace_back(rmax);
324  } else {
325  HGCalGeomTools::radius(zz - hthick,
326  zz + hthick,
327  zFrontB_,
328  rMinFront_,
329  slopeB_,
330  zFrontT_,
331  rMaxFront_,
332  slopeT_,
333  -layerSense_[ly],
334  pgonZ,
335  pgonRin,
336  pgonRout);
337  for (unsigned int isec = 0; isec < pgonZ.size(); ++isec) {
338  pgonZ[isec] -= zz;
339  pgonRout[isec] = pgonRout[isec] * cosAlpha_ - tol1_;
340  }
341  }
342  DDSolid solid =
343  DDSolidFactory::polyhedra(DDName(name, nameSpace_), sectors_, -alpha_, 2._pi, pgonZ, pgonRin, pgonRout);
344  glog = DDLogicalPart(solid.ddname(), matter, solid);
345 #ifdef EDM_ML_DEBUG
346  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " polyhedra of " << sectors_
347  << " sectors covering " << convertRadToDeg(-alpha_) << ":"
348  << convertRadToDeg(-alpha_ + 2._pi) << " with " << pgonZ.size() << " sections";
349  for (unsigned int k = 0; k < pgonZ.size(); ++k)
350  edm::LogVerbatim("HGCalGeom") << "[" << k << "] z " << pgonZ[k] << " R " << pgonRin[k] << ":" << pgonRout[k];
351 #endif
352  } else {
353  double rins = (sensitiveMode_ < 1) ? rinB : HGCalGeomTools::radius(zz + hthick, zFrontB_, rMinFront_, slopeB_);
354  double routs =
355  (sensitiveMode_ < 1) ? routF : HGCalGeomTools::radius(zz - hthick, zFrontT_, rMaxFront_, slopeT_);
356  DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rins, routs, 0.0, 2._pi);
357  glog = DDLogicalPart(solid.ddname(), matter, solid);
358 #ifdef EDM_ML_DEBUG
359  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
360  << " of dimensions " << rinB << ":" << rins << ", " << routF << ":" << routs
361  << ", " << hthick << ", 0.0, 360.0 and positioned in: " << glog.name()
362  << " number " << copy;
363 #endif
364  positionMix(glog, name, copy, thick_[ii], matter, rins, rMixLayer_[i], routs, zz, cpv);
365  }
366  DDTranslation r1(0, 0, zz);
367  DDRotation rot;
368  cpv.position(glog, module, copy, r1, rot);
369  ++copyNumber_[ii];
370 #ifdef EDM_ML_DEBUG
371  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << glog.name() << " number " << copy << " positioned in "
372  << module.name() << " at " << r1 << " with no rotation";
373 #endif
374  zz += hthick;
375  } // End of loop over layers in a block
376  zi = zo;
377  laymin = laymax;
378  if (std::abs(thickTot - layerThick_[i]) >= tol2_) {
379  if (thickTot > layerThick_[i]) {
380  edm::LogError("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " is smaller than " << thickTot
381  << ": thickness of all its components **** ERROR ****";
382  } else {
383  edm::LogWarning("HGCalGeom") << "Thickness of the partition " << layerThick_[i] << " does not match with "
384  << thickTot << " of the components";
385  }
386  }
387  } // End of loop over blocks
388 }
389 
391  const std::string& nameM,
392  int copyM,
393  double thick,
394  const DDMaterial& matter,
395  double rin,
396  double rmid,
397  double rout,
398  double zz,
399  DDCompactView& cpv) {
400  DDLogicalPart glog1;
401  DDTranslation tran;
402  DDRotation rot;
403  for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
404  int ii = layerTypeTop_[ly];
405  copyNumberTop_[ii] = copyM;
406  }
407  for (unsigned int ly = 0; ly < layerTypeBot_.size(); ++ly) {
408  int ii = layerTypeBot_[ly];
409  copyNumberBot_[ii] = copyM;
410  }
411  double hthick = 0.5 * thick;
412  // Make the top part first
413  std::string name = nameM + "Top";
414  DDSolid solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rmid, rout, 0.0, 2._pi);
415  glog1 = DDLogicalPart(solid.ddname(), matter, solid);
416 #ifdef EDM_ML_DEBUG
417  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matter.name()
418  << " of dimensions " << rmid << ", " << rout << ", " << hthick << ", 0.0, 360.0";
419 #endif
420  cpv.position(glog1, glog, 1, tran, rot);
421 #ifdef EDM_ML_DEBUG
422  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << glog1.name() << " number 1 positioned in " << glog.name()
423  << " at " << tran << " with no rotation";
424 #endif
425  double thickTot(0), zpos(-hthick);
426  for (unsigned int ly = 0; ly < layerTypeTop_.size(); ++ly) {
427  int ii = layerTypeTop_[ly];
428  int copy = copyNumberTop_[ii];
429  double hthickl = 0.5 * layerThickTop_[ii];
430  thickTot += layerThickTop_[ii];
431  name = namesTop_[ii] + std::to_string(copy);
432 #ifdef EDM_ML_DEBUG
433  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " R " << rmid << ":" << rout
434  << " Thick " << layerThickTop_[ii];
435 #endif
436  DDName matName(DDSplit(materialsTop_[ii]).first, DDSplit(materialsTop_[ii]).second);
437  DDMaterial matter1(matName);
438  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthickl, rmid, rout, 0.0, 2._pi);
439  DDLogicalPart glog2 = DDLogicalPart(solid.ddname(), matter1, solid);
440 #ifdef EDM_ML_DEBUG
441  double eta1 = -log(tan(0.5 * atan(rmid / zz)));
442  double eta2 = -log(tan(0.5 * atan(rout / zz)));
443  edm::LogVerbatim("HGCalGeom") << name << " z|rin|rout " << zz << ":" << rmid << ":" << rout << " eta " << eta1
444  << ":" << eta2;
445  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
446  << " of dimensions " << rmid << ", " << rout << ", " << hthickl << ", 0.0, 360.0";
447 #endif
448  zpos += hthickl;
449  DDTranslation r1(0, 0, zpos);
450  cpv.position(glog2, glog1, copy, r1, rot);
451 #ifdef EDM_ML_DEBUG
452  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Position " << glog2.name() << " number " << copy << " in "
453  << glog1.name() << " at " << r1 << " with no rotation";
454 #endif
455  ++copyNumberTop_[ii];
456  zpos += hthickl;
457  }
458  if (std::abs(thickTot - thick) >= tol2_) {
459  if (thickTot > thick) {
460  edm::LogError("HGCalGeom") << "Thickness of the partition " << thick << " is smaller than " << thickTot
461  << ": thickness of all its components in the top part **** ERROR ****";
462  } else {
463  edm::LogWarning("HGCalGeom") << "Thickness of the partition " << thick << " does not match with " << thickTot
464  << " of the components in top part";
465  }
466  }
467 
468  // Make the bottom part next
469  name = nameM + "Bottom";
470  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthick, rin, rmid, 0.0, 2._pi);
471  glog1 = DDLogicalPart(solid.ddname(), matter, solid);
472 #ifdef EDM_ML_DEBUG
473  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matter.name()
474  << " of dimensions " << rin << ", " << rmid << ", " << hthick << ", 0.0, 360.0";
475 #endif
476  cpv.position(glog1, glog, 1, tran, rot);
477 #ifdef EDM_ML_DEBUG
478  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << glog1.name() << " number 1 positioned in " << glog.name()
479  << " at " << tran << " with no rotation";
480 #endif
481  thickTot = 0;
482  zpos = -hthick;
483  for (unsigned int ly = 0; ly < layerTypeBot_.size(); ++ly) {
484  int ii = layerTypeBot_[ly];
485  int copy = copyNumberBot_[ii];
486  double hthickl = 0.5 * layerThickBot_[ii];
487  thickTot += layerThickBot_[ii];
488  name = namesBot_[ii] + std::to_string(copy);
489 #ifdef EDM_ML_DEBUG
490  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Layer " << ly << ":" << ii << " R " << rin << ":" << rmid
491  << " Thick " << layerThickBot_[ii];
492 #endif
493  DDName matName(DDSplit(materialsBot_[ii]).first, DDSplit(materialsBot_[ii]).second);
494  DDMaterial matter1(matName);
495  solid = DDSolidFactory::tubs(DDName(name, nameSpace_), hthickl, rin, rmid, 0.0, 2._pi);
496  DDLogicalPart glog2 = DDLogicalPart(solid.ddname(), matter1, solid);
497 #ifdef EDM_ML_DEBUG
498  double eta1 = -log(tan(0.5 * atan(rin / zz)));
499  double eta2 = -log(tan(0.5 * atan(rmid / zz)));
500  edm::LogVerbatim("HGCalGeom") << name << " z|rin|rout " << zz << ":" << rin << ":" << rmid << " eta " << eta1 << ":"
501  << eta2;
502  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << solid.name() << " Tubs made of " << matName
503  << " of dimensions " << rin << ", " << rmid << ", " << hthickl << ", 0.0, 360.0";
504 #endif
505  zpos += hthickl;
506  DDTranslation r1(0, 0, zpos);
507  cpv.position(glog2, glog1, copy, r1, rot);
508 #ifdef EDM_ML_DEBUG
509  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Position " << glog2.name() << " number " << copy << " in "
510  << glog1.name() << " at " << r1 << " with no rotation";
511 #endif
512  if (layerSenseBot_[ly] != 0) {
513 #ifdef EDM_ML_DEBUG
514  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: z " << (zz + zpos) << " Center " << copy << ":"
515  << (copy - firstLayer_) << ":" << layerCenter_[copy - firstLayer_];
516 #endif
517  positionSensitive(glog2, rin, rmid, zz + zpos, layerSenseBot_[ly], layerCenter_[copy - firstLayer_], cpv);
518  }
519  zpos += hthickl;
520  ++copyNumberBot_[ii];
521  }
522  if (std::abs(thickTot - thick) >= tol2_) {
523  if (thickTot > thick) {
524  edm::LogError("HGCalGeom") << "Thickness of the partition " << thick << " is smaller than " << thickTot
525  << ": thickness of all its components in the top part **** ERROR ****";
526  } else {
527  edm::LogWarning("HGCalGeom") << "Thickness of the partition " << thick << " does not match with " << thickTot
528  << " of the components in top part";
529  }
530  }
531 }
532 
534  double rin,
535  double rout,
536  double zpos,
537  int layertype,
538  int layercenter,
539  DDCompactView& cpv) {
540  static const double sqrt3 = std::sqrt(3.0);
541  double r = 0.5 * (waferSize_ + waferSepar_);
542  double R = 2.0 * r / sqrt3;
543  double dy = 0.75 * R;
544  int N = (int)(0.5 * rout / r) + 2;
545  std::pair<double, double> xyoff = geomTools_.shiftXY(layercenter, (waferSize_ + waferSepar_));
546 #ifdef EDM_ML_DEBUG
547  int ium(0), ivm(0), iumAll(0), ivmAll(0), kount(0), ntot(0), nin(0);
548  std::vector<int> ntype(6, 0);
549  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << glog.ddname() << " rin:rout " << rin << ":" << rout << " zpos "
550  << zpos << " N " << N << " for maximum u, v Offset; Shift " << xyoff.first << ":"
551  << xyoff.second << " WaferSize " << (waferSize_ + waferSepar_);
552 #endif
553  for (int u = -N; u <= N; ++u) {
554  for (int v = -N; v <= N; ++v) {
555  int nr = 2 * v;
556  int nc = -2 * u + v;
557  double xpos = xyoff.first + nc * r;
558  double ypos = xyoff.second + nr * dy;
559  std::pair<int, int> corner = HGCalGeomTools::waferCorner(xpos, ypos, r, R, rin, rout, false);
560 #ifdef EDM_ML_DEBUG
561  int iu = std::abs(u);
562  int iv = std::abs(v);
563  ++ntot;
564 #endif
565  if (corner.first > 0) {
566  int type = waferType_->getType(xpos, ypos, zpos);
567  int copy = HGCalTypes::packTypeUV(type, u, v);
568  if (layertype > 1)
569  type += 3;
570 #ifdef EDM_ML_DEBUG
571  if (iu > ium)
572  ium = iu;
573  if (iv > ivm)
574  ivm = iv;
575  kount++;
576  if (copies_.count(copy) == 0)
577  copies_.insert(copy);
578 #endif
579  if (corner.first == (int)(HGCalParameters::k_CornerSize)) {
580 #ifdef EDM_ML_DEBUG
581  if (iu > iumAll)
582  iumAll = iu;
583  if (iv > ivmAll)
584  ivmAll = iv;
585  ++nin;
586 #endif
587  DDTranslation tran(xpos, ypos, 0.0);
589  DDName name = DDName(DDSplit(wafers_[type]).first, DDSplit(wafers_[type]).second);
590  cpv.position(name, glog.ddname(), copy, tran, rotation);
591 #ifdef EDM_ML_DEBUG
592  ++ntype[type];
593  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: " << name << " number " << copy << " positioned in "
594  << glog.ddname() << " at " << tran << " with no rotation";
595 #endif
596  }
597  }
598  }
599  }
600 #ifdef EDM_ML_DEBUG
601  edm::LogVerbatim("HGCalGeom") << "DDHGCalHEAlgo: Maximum # of u " << ium << ":" << iumAll << " # of v " << ivm << ":"
602  << ivmAll << " and " << nin << ":" << kount << ":" << ntot << " wafers (" << ntype[0]
603  << ":" << ntype[1] << ":" << ntype[2] << ":" << ntype[3] << ":" << ntype[4] << ":"
604  << ntype[5] << ") for " << glog.ddname() << " R " << rin << ":" << rout;
605 #endif
606 }
607 
608 DEFINE_EDM_PLUGIN(DDAlgorithmFactory, DDHGCalHEAlgo, "hgcal:DDHGCalHEAlgo");
std::vector< double > rMixLayer_
Log< level::Info, true > LogVerbatim
static AlgebraicMatrix initialize()
std::vector< int > layerTypeTop_
std::vector< double > layerThickBot_
static std::vector< std::string > checklist log
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:59
int32_t *__restrict__ iv
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)
void position(const DDLogicalPart &self, const DDLogicalPart &parent, const std::string &copyno, const DDTranslation &trans, const DDRotation &rot, const DDDivision *div=nullptr)
std::vector< double > slopeB_
std::vector< std::string > materials_
DDMaterial is used to define and access material information.
Definition: DDMaterial.h:45
std::string nameSpace_
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
std::vector< double > zFrontB_
void execute(DDCompactView &cpv) override
DDName is used to identify DDD entities uniquely.
Definition: DDName.h:17
std::string to_string(const V &value)
Definition: OMSAccess.h:71
void initialize(const DDNumericArguments &nArgs, const DDVectorArguments &vArgs, const DDMapArguments &mArgs, const DDStringArguments &sArgs, const DDStringVectorArguments &vsArgs) override
static std::string & ns()
std::vector< double > thick_
Log< level::Error, false > LogError
int ii
Definition: cuy.py:589
std::vector< std::string > namesTop_
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:81
~DDHGCalHEAlgo() override
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
static constexpr uint32_t k_CornerSize
int nin
Represents a uniquely identifyable rotation matrix.
Definition: DDTransform.h:57
U second(std::pair< T, U > const &p)
static std::pair< int32_t, int32_t > waferCorner(double xpos, double ypos, double r, double R, double rMin, double rMax, bool oldBug=false)
std::unique_ptr< HGCalWaferType > waferType_
std::vector< int > layerSense_
std::vector< int > copyNumberBot_
T sqrt(T t)
Definition: SSEVec.h:19
std::vector< std::string > materialsBot_
std::vector< double > slopeT_
void positionSensitive(const DDLogicalPart &glog, double rin, double rout, double zpos, int layertype, int layercenter, DDCompactView &cpv)
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
void constructLayers(const cms::DDNamespace &ns, const std::vector< std::string > &wafers, const std::vector< std::string > &covers, const std::vector< int > &layerType, const std::vector< int > &layerSense, const std::vector< int > &maxModule, const std::vector< std::string > &names, const std::vector< std::string > &materials, std::vector< int > &copyNumber, const std::vector< double > &layerThick, const double &absorbW, const double &absorbH, const double &waferTot, const double &rMax, const double &rMaxFine, std::unordered_set< int > &copies, int firstLayer, int lastLayer, double zFront, double totalWidth, bool ignoreCenter, dd4hep::Volume &module)
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:667
std::vector< int > layers_
std::vector< double > rMinFront_
HGCalGeomTools geomTools_
std::vector< int > layerCenter_
std::vector< double > rad100to200_
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
std::vector< double > zFrontT_
#define N
Definition: blowfish.cc:9
std::vector< double > layerThickTop_
std::unordered_set< int > copies_
std::vector< int > copyNumber_
std::vector< int > layerType_
std::vector< std::string > materialsTop_
std::vector< int > copyNumberTop_
static int32_t packTypeUV(int type, int u, int v)
Definition: HGCalTypes.cc:3
std::vector< int > layerSenseBot_
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< std::string > wafers_
Log< level::Warning, false > LogWarning
std::pair< std::string, std::string > DDSplit(const std::string &n)
split into (name,namespace), separator = &#39;:&#39;
Definition: DDSplit.cc:3
std::vector< std::string > namesBot_
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DDTranslation
Definition: DDTranslation.h:7
std::vector< double > layerThick_
__shared__ uint32_t ntot
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:565
std::vector< std::string > names_
std::vector< int > layerTypeBot_
tuple module
Definition: callgraph.py:69
std::vector< double > rMaxFront_
std::vector< double > rad200to300_
const N & ddname() const
Definition: DDBase.h:61