CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
DDHCalBarrelAlgo.cc
Go to the documentation of this file.
1 // File: DDHCalBarrelAlgo.cc
3 // adapted from CCal(G4)HcalBarrel.cc
4 // Description: Geometry factory class for Hcal Barrel
6 
7 #include <cmath>
8 #include <algorithm>
9 #include <map>
10 #include <string>
11 #include <vector>
12 
25 
26 //#define EDM_ML_DEBUG
27 using namespace angle_units::operators;
28 
29 class DDHCalBarrelAlgo : public DDAlgorithm {
30 public:
31  //Constructor and Destructor
33  ~DDHCalBarrelAlgo() override;
34 
35  //Get Methods
36  const std::string& getGenMaterial() const { return genMaterial; }
37  int getNsectors() const { return nsectors; }
38  int getNsectortot() const { return nsectortot; }
39  int getNhalf() const { return nhalf; }
40  double getRin() const { return rin; }
41  double getRout() const { return rout; }
42  int getRzones() const { return rzones; }
43  double getTanTheta(unsigned int i) const { return ttheta[i]; }
44  double getTheta(unsigned int i) const { return theta[i]; }
45  double getRmax(unsigned int i) const { return rmax[i]; }
46  double getZoff(unsigned int i) const { return zoff[i]; }
47 
48  int getNLayers() const { return nLayers; }
49  int getLayerId(unsigned i) const { return layerId[i]; }
50  const std::string& getLayerLabel(unsigned i) const { return layerLabel[i]; }
51  const std::string& getLayerMaterial(unsigned i) const { return layerMat[i]; }
52  double getLayerWidth(unsigned i) const { return layerWidth[i]; }
53  double getLayerD1(unsigned i) const { return layerD1[i]; }
54  double getLayerD2(unsigned i) const { return layerD2[i]; }
55  double getLayerAlpha(unsigned i) const { return layerAlpha[i]; }
56  double getLayerT1(unsigned i) const { return layerT1[i]; }
57  double getLayerT2(unsigned i) const { return layerT2[i]; }
58  int getLayerAbsorb(unsigned int i) const { return layerAbsorb[i]; }
59  double getLayerGap(unsigned int i) const { return layerGap[i]; }
60 
61  const std::string& getSideMat(unsigned int i) const { return sideMat[i]; }
62  double getSideD(unsigned int i) const { return sideD[i]; }
63  double getSideT(unsigned int i) const { return sideT[i]; }
64  int getSideAbsorber() const { return nSideAbs; }
65  const std::string& getSideAbsName(unsigned int i) const { return sideAbsName[i]; }
66  const std::string& getSideAbsMat(unsigned int i) const { return sideAbsMat[i]; }
67  double getSideAbsW(unsigned int i) const { return sideAbsW[i]; }
68 
69  int getAbsorberN() const { return nAbsorber; }
70  const std::string& getAbsorbName(unsigned int i) const { return absorbName[i]; }
71  const std::string& getAbsorbMat(unsigned int i) const { return absorbMat[i]; }
72  double getAbsorbD(unsigned int i) const { return absorbD[i]; }
73  double getAbsorbT(unsigned int i) const { return absorbT[i]; }
74  const std::string& getMiddleMat() const { return middleMat; }
75  double getMiddleD() const { return middleD; }
76  double getMiddleW() const { return middleW; }
77  int getMidAbsorber() const { return nMidAbs; }
78  const std::string& getMidAbsName(unsigned int i) const { return midName[i]; }
79  const std::string& getMidAbsMat(unsigned int i) const { return midMat[i]; }
80  double getMidAbsW(unsigned int i) const { return midW[i]; }
81  double getMidAbsT(unsigned int i) const { return midT[i]; }
82 
83  const std::string& getDetMat() const { return detMat; }
84  const std::string& getDetMatPl() const { return detMatPl; }
85  const std::string& getDetMatSc() const { return detMatSc; }
86  int getDetType(unsigned int i) const { return detType[i]; }
87  double getDetdP1(unsigned int i) const { return detdP1[i]; }
88  double getDetdP2(unsigned int i) const { return detdP2[i]; }
89  double getDetT11(unsigned int i) const { return detT11[i]; }
90  double getDetT12(unsigned int i) const { return detT12[i]; }
91  double getDetTsc(unsigned int i) const { return detTsc[i]; }
92  double getDetT21(unsigned int i) const { return detT21[i]; }
93  double getDetT22(unsigned int i) const { return detT22[i]; }
94  double getDetWidth1(unsigned int i) const { return detWidth1[i]; }
95  double getDetWidth2(unsigned int i) const { return detWidth2[i]; }
96  int getDetPosY(unsigned int i) const { return detPosY[i]; }
97 
98  void initialize(const DDNumericArguments& nArgs,
99  const DDVectorArguments& vArgs,
100  const DDMapArguments& mArgs,
101  const DDStringArguments& sArgs,
102  const DDStringVectorArguments& vsArgs) override;
103 
104  void execute(DDCompactView& cpv) override;
105 
106 private:
107  void constructGeneralVolume(DDCompactView& cpv);
108  void constructInsideSector(const DDLogicalPart& sector, DDCompactView& cpv);
109  void constructInsideLayers(const DDLogicalPart& laylog,
110  const std::string& name,
111  int id,
112  int nAbs,
113  double rin,
114  double d1,
115  double alpha1,
116  double d2,
117  double alpha2,
118  double t1,
119  double t2,
120  DDCompactView& cpv);
121  DDLogicalPart constructSideLayer(
122  const DDLogicalPart& laylog, const std::string& nm, int nAbs, double rin, double alpha, DDCompactView& cpv);
123  DDLogicalPart constructMidLayer(
124  const DDLogicalPart& laylog, const std::string& nm, double rin, double alpha, DDCompactView& cpv);
125  void constructInsideDetectors(const DDLogicalPart& detector,
126  const std::string& name,
127  int id,
128  double dx,
129  double dy,
130  double dz,
131  int type,
132  DDCompactView& cpv);
133 
134  //General Volume
135  // <----- Zmax ------>
136  // Rout ************************-------
137  // * *Rstep2| Theta angle w.r.t. vertical
138  // * *---------------
139  // * * |
140  // * *Theta[i] Rmax[i]
141  // * *--------------- |
142  // *Theta[0] Rmax[0]| |
143  // Rin *****************----------------------
144 
145  std::string genMaterial; //General material
146  int nsectors; //Number of potenital straight edges
147  int nsectortot; //Number of straight edges (actual)
148  int nhalf; //Number of half modules
149  double rin, rout; //See picture
150  int rzones; // ....
151  std::vector<double> theta; // .... (in degrees)
152  std::vector<double> rmax; // ....
153  std::vector<double> zoff; // ....
154  std::vector<double> ttheta; //tan(theta)
155  std::string rotHalf; //Rotation matrix of the second half
156  std::string rotns; //Name space for Rotation matrices
157 
158  //Upper layers inside general volume
159  // <---- Zout ---->
160  // | **************** |
161  // | * * Wstep
162  // W * ***** |
163  // | * *
164  // | ********************
165  // <------ Zin ------->
166  // Middle layers inside general volume
167  // <------ Zout ------> Zout = Full sector Z at position
168  // | ******************** Zin = Full sector Z at position
169  // | * *
170  // W * * Angle = Theta sector
171  // | * * )
172  // | ****************--------
173  // <------ Zin ------->
174 
175  // Lower layers
176  // <------ Zout ------> Zin(i)=Zout(i-1)
177  // | ******************** Zout(i)=Zin(i)+W(i)/tan(Theta(i))
178  // | * *
179  // W * * Theta
180  // | * *
181  // | ****************--------
182  // <--- Zin ------>
183 
184  int nLayers; //Number of layers
185  std::vector<int> layerId; //Number identification
186  std::vector<std::string> layerLabel; //String identification
187  std::vector<std::string> layerMat; //Material
188  std::vector<double> layerWidth; //W in picture
189  std::vector<double> layerD1; //d1 in front picture
190  std::vector<double> layerD2; //d2 in front picture
191  std::vector<double> layerAlpha; //Angular width of the middle tiles
192  std::vector<double> layerT1; //t in front picture (side)
193  std::vector<double> layerT2; //t in front picture (front)
194  std::vector<int> layerAbsorb; //Absorber flag
195  std::vector<double> layerGap; //Gap at the edge
196 
197  int nAbsorber; //Number of absorber layers in middle
198  std::vector<std::string> absorbName; //Absorber name
199  std::vector<std::string> absorbMat; //Absorber material
200  std::vector<double> absorbD; //Distance from the bottom surface
201  std::vector<double> absorbT; //Thickness
202  std::string middleMat; //Material of the detector layer
203  double middleD; //Distance from the bottom surface
204  double middleW; //Half width
205  int nMidAbs; //Number of absorbers in front part
206  std::vector<std::string> midName; //Absorber names in the front part
207  std::vector<std::string> midMat; //Absorber material
208  std::vector<double> midW; //Half width
209  std::vector<double> midT; //Thickness
210 
211  std::vector<std::string> sideMat; //Material for special side layers
212  std::vector<double> sideD; //Depth from bottom surface
213  std::vector<double> sideT; //Thickness
214  int nSideAbs; //Number of absorbers in special side
215  std::vector<std::string> sideAbsName; //Absorber name
216  std::vector<std::string> sideAbsMat; //Absorber material
217  std::vector<double> sideAbsW; //Half width
218 
219  // Detectors. Each volume inside the layer has the shape:
220  //
221  // ******************************* |
222  // *\\\\\\\Plastic\\\\\\\\\\\\\\\* T2
223  // ******************************* |
224  // *////Scintillator/////////////* Tsc
225  // ******************************* |
226  // *\\\\\\\Plastic\\\\\\\\\\\\\\\* T1
227  // ******************************* | |
228  // * Air * dP1
229  // ******************************* |
230  //
231  std::string detMat; //fill material
232  std::string detRot; //Rotation matrix for the 2nd
233  std::string detMatPl; //Plastic material
234  std::string detMatSc; //Scintillator material
235  std::vector<int> detType;
236  std::vector<double> detdP1; //Air gap (side)
237  std::vector<double> detdP2; //Air gap (centre)
238  std::vector<double> detT11; //Back plastic thickness (side)
239  std::vector<double> detT12; //Back plastic thickness (centre)
240  std::vector<double> detTsc; //Scintillator
241  std::vector<double> detT21; //Front plastic thickness (side)
242  std::vector<double> detT22; //Front plastic thickness (centre)
243  std::vector<double> detWidth1; //Width of phi(1,4) megatiles
244  std::vector<double> detWidth2; //Width of phi(2,3) megatiles
245  std::vector<int> detPosY; //Positioning of phi(1,4) tiles - 0 centre
246 
247  std::string idName; //Name of the "parent" volume.
248  std::string idNameSpace; //Namespace of this and ALL sub-parts
249  int idOffset; // Geant4 ID's... = 3000;
250 };
251 
253  : theta(0),
254  rmax(0),
255  zoff(0),
256  ttheta(0),
257  layerId(0),
258  layerLabel(0),
259  layerMat(0),
260  layerWidth(0),
261  layerD1(0),
262  layerD2(0),
263  layerAlpha(0),
264  layerT1(0),
265  layerT2(0),
266  layerAbsorb(0),
267  layerGap(0),
268  absorbName(0),
269  absorbMat(0),
270  absorbD(0),
271  absorbT(0),
272  midName(0),
273  midMat(0),
274  midW(0),
275  midT(0),
276  sideMat(0),
277  sideD(0),
278  sideT(0),
279  sideAbsName(0),
280  sideAbsMat(0),
281  sideAbsW(0),
282  detType(0),
283  detdP1(0),
284  detdP2(0),
285  detT11(0),
286  detT12(0),
287  detTsc(0),
288  detT21(0),
289  detT22(0),
290  detWidth1(0),
291  detWidth2(0),
292  detPosY(0) {
293 #ifdef EDM_ML_DEBUG
294  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Creating an instance";
295 #endif
296 }
297 
299 
301  const DDVectorArguments& vArgs,
302  const DDMapArguments&,
303  const DDStringArguments& sArgs,
304  const DDStringVectorArguments& vsArgs) {
305  genMaterial = sArgs["MaterialName"];
306  nsectors = int(nArgs["NSector"]);
307  nsectortot = int(nArgs["NSectorTot"]);
308  nhalf = int(nArgs["NHalf"]);
309  rin = nArgs["RIn"];
310  rout = nArgs["ROut"];
311  rzones = int(nArgs["RZones"]);
312  rotHalf = sArgs["RotHalf"];
313  rotns = sArgs["RotNameSpace"];
314 
315  theta = vArgs["Theta"];
316  rmax = vArgs["RMax"];
317  zoff = vArgs["ZOff"];
318  for (int i = 0; i < rzones; i++) {
319  ttheta.emplace_back(tan(theta[i])); //*deg already done in XML
320  }
321  if (rzones > 3)
322  rmax[2] = (zoff[3] - zoff[2]) / ttheta[2];
323 
324 #ifdef EDM_ML_DEBUG
325  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: General material " << genMaterial << "\tSectors " << nsectors
326  << ", " << nsectortot << "\tHalves " << nhalf << "\tRotation matrix " << rotns << ":"
327  << rotHalf << "\n\t\t" << rin << "\t" << rout << "\t" << rzones;
328  for (int i = 0; i < rzones; i++)
329  edm::LogVerbatim("HCalGeom") << "\tTheta[" << i << "] = " << theta[i] << "\trmax[" << i << "] = " << rmax[i]
330  << "\tzoff[" << i << "] = " << zoff[i];
331 #endif
332  //Layers
334  nLayers = int(nArgs["NLayers"]);
335 #ifdef EDM_ML_DEBUG
336  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Layer\t" << nLayers;
337 #endif
338  layerId = dbl_to_int(vArgs["Id"]);
339  layerLabel = vsArgs["LayerLabel"];
340  layerMat = vsArgs["LayerMat"];
341  layerWidth = vArgs["LayerWidth"];
342  layerD1 = vArgs["D1"];
343  layerD2 = vArgs["D2"];
344  layerAlpha = vArgs["Alpha2"];
345  layerT1 = vArgs["T1"];
346  layerT2 = vArgs["T2"];
347  layerAbsorb = dbl_to_int(vArgs["AbsL"]);
348  layerGap = vArgs["Gap"];
349 #ifdef EDM_ML_DEBUG
350  for (int i = 0; i < nLayers; i++)
351  edm::LogVerbatim("HCalGeom") << layerLabel[i] << "\t" << layerId[i] << "\t" << layerMat[i] << "\t" << layerWidth[i]
352  << "\t" << layerD1[i] << "\t" << layerD2[i] << "\t" << layerAlpha[i] << "\t"
353  << layerT1[i] << "\t" << layerT2[i] << "\t" << layerAbsorb[i] << "\t" << layerGap[i];
354 #endif
355 
357  //Absorber Layers and middle part
358  absorbName = vsArgs["AbsorbName"];
359  absorbMat = vsArgs["AbsorbMat"];
360  absorbD = vArgs["AbsorbD"];
361  absorbT = vArgs["AbsorbT"];
362  nAbsorber = absorbName.size();
363 #ifdef EDM_ML_DEBUG
364  for (int i = 0; i < nAbsorber; i++)
365  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << absorbName[i] << " Material " << absorbMat[i] << " d "
366  << absorbD[i] << " t " << absorbT[i];
367 #endif
368  middleMat = sArgs["MiddleMat"];
369  middleD = nArgs["MiddleD"];
370  middleW = nArgs["MiddleW"];
371 #ifdef EDM_ML_DEBUG
372  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Middle material " << middleMat << " d " << middleD << " w "
373  << middleW;
374 #endif
375  midName = vsArgs["MidAbsName"];
376  midMat = vsArgs["MidAbsMat"];
377  midW = vArgs["MidAbsW"];
378  midT = vArgs["MidAbsT"];
379  nMidAbs = midName.size();
380 #ifdef EDM_ML_DEBUG
381  for (int i = 0; i < nMidAbs; i++)
382  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << midName[i] << " Material " << midMat[i] << " W " << midW[i]
383  << " T " << midT[i];
384 #endif
385 
386  //Absorber layers in the side part
387  sideMat = vsArgs["SideMat"];
388  sideD = vArgs["SideD"];
389  sideT = vArgs["SideT"];
390 #ifdef EDM_ML_DEBUG
391  for (unsigned int i = 0; i < sideMat.size(); i++)
392  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Side material " << sideMat[i] << " d " << sideD[i] << " t "
393  << sideT[i];
394 #endif
395  sideAbsName = vsArgs["SideAbsName"];
396  sideAbsMat = vsArgs["SideAbsMat"];
397  sideAbsW = vArgs["SideAbsW"];
398  nSideAbs = sideAbsName.size();
399 #ifdef EDM_ML_DEBUG
400  for (int i = 0; i < nSideAbs; i++)
401  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << sideAbsName[i] << " Material " << sideAbsMat[i] << " W "
402  << sideAbsW[i];
403 #endif
404 
406  // Detectors
407 
408  detMat = sArgs["DetMat"];
409  detRot = sArgs["DetRot"];
410  detMatPl = sArgs["DetMatPl"];
411  detMatSc = sArgs["DetMatSc"];
412 #ifdef EDM_ML_DEBUG
413  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Detector (" << nLayers << ") Rotation matrix " << rotns << ":"
414  << detRot << "\n\t\t" << detMat << "\t" << detMatPl << "\t" << detMatSc;
415 #endif
416  detType = dbl_to_int(vArgs["DetType"]);
417  detdP1 = vArgs["DetdP1"];
418  detdP2 = vArgs["DetdP2"];
419  detT11 = vArgs["DetT11"];
420  detT12 = vArgs["DetT12"];
421  detTsc = vArgs["DetTsc"];
422  detT21 = vArgs["DetT21"];
423  detT22 = vArgs["DetT22"];
424  detWidth1 = vArgs["DetWidth1"];
425  detWidth2 = vArgs["DetWidth2"];
426  detPosY = dbl_to_int(vArgs["DetPosY"]);
427 #ifdef EDM_ML_DEBUG
428  for (int i = 0; i < nLayers; i++)
429  edm::LogVerbatim("HCalGeom") << i + 1 << "\t" << detType[i] << "\t" << detdP1[i] << ", " << detdP2[i] << "\t"
430  << detT11[i] << ", " << detT12[i] << "\t" << detTsc[i] << "\t" << detT21[i] << ", "
431  << detT22[i] << "\t" << detWidth1[i] << "\t" << detWidth2[i] << "\t" << detPosY[i];
432 #endif
433 
434  // idName = parentName.name();
435  idName = sArgs["MotherName"];
437  idOffset = int(nArgs["IdOffset"]);
438 #ifdef EDM_ML_DEBUG
439  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Parent " << parent().name() << " idName " << idName
440  << " NameSpace " << idNameSpace << " Offset " << idOffset;
441 #endif
442 }
443 
445 // DDHCalBarrelAlgo methods...
447 
449 #ifdef EDM_ML_DEBUG
450  edm::LogVerbatim("HCalGeom") << "==>> Constructing DDHCalBarrelAlgo...";
451 #endif
453 #ifdef EDM_ML_DEBUG
454  edm::LogVerbatim("HCalGeom") << "<<== End of DDHCalBarrelAlgo construction";
455 #endif
456 }
457 
458 //----------------------start here for DDD work!!! ---------------
459 
461 #ifdef EDM_ML_DEBUG
462  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: General volume...";
463 #endif
464 
466 
467  double alpha = (1._pi) / getNsectors();
468  double dphi = getNsectortot() * (2._pi) / getNsectors();
469  int nsec, ntot = 15;
470  if (getNhalf() == 1)
471  nsec = 8;
472  else
473  nsec = 15;
474  int nf = ntot - nsec;
475 
476  //Calculate zmin... see HCalBarrel.hh picture. For polyhedra
477  //Rmin and Rmax are distances to vertex
478  double zmax = getZoff(3);
479  double zstep5 = getZoff(4);
480  double zstep4 = (getZoff(1) + getRmax(1) * getTanTheta(1));
481  if ((getZoff(2) + getRmax(1) * getTanTheta(2)) > zstep4)
482  zstep4 = (getZoff(2) + getRmax(1) * getTanTheta(2));
483  double zstep3 = (getZoff(1) + getRmax(0) * getTanTheta(1));
484  double zstep2 = (getZoff(0) + getRmax(0) * getTanTheta(0));
485  double zstep1 = (getZoff(0) + getRin() * getTanTheta(0));
486  double rout = getRout();
487  double rout1 = getRmax(3);
488  double rin = getRin();
489  double rmid1 = getRmax(0);
490  double rmid2 = getRmax(1);
491  double rmid3 = (getZoff(4) - getZoff(2)) / getTanTheta(2);
492  double rmid4 = getRmax(2);
493 
494  std::vector<double> pgonZ = {-zmax,
495  -zstep5,
496  -zstep5,
497  -zstep4,
498  -zstep3,
499  -zstep2,
500  -zstep1,
501  0,
502  zstep1,
503  zstep2,
504  zstep3,
505  zstep4,
506  zstep5,
507  zstep5,
508  zmax};
509 
510  std::vector<double> pgonRmin = {
511  rmid4, rmid3, rmid3, rmid2, rmid1, rmid1, rin, rin, rin, rmid1, rmid1, rmid2, rmid3, rmid3, rmid4};
512 
513  std::vector<double> pgonRmax = {
514  rout1, rout1, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout, rout1, rout1};
515 
516  std::vector<double> pgonZHalf = {0, zstep1, zstep2, zstep3, zstep4, zstep5, zstep5, zmax};
517 
518  std::vector<double> pgonRminHalf = {rin, rin, rmid1, rmid1, rmid2, rmid3, rmid3, rmid4};
519 
520  std::vector<double> pgonRmaxHalf = {rout, rout, rout, rout, rout, rout, rout1, rout1};
521 
522  std::string name("Null");
523  DDSolid solid;
524  if (nf == 0) {
526  DDName(idName, idNameSpace), getNsectortot(), -alpha, dphi, pgonZ, pgonRmin, pgonRmax);
527 #ifdef EDM_ML_DEBUG
528  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << DDName(idName, idNameSpace) << " Polyhedra made of "
529  << getGenMaterial() << " with " << getNsectortot() << " sectors from "
530  << convertRadToDeg(-alpha) << " to " << convertRadToDeg(-alpha + dphi) << " and with "
531  << nsec << " sections ";
532  for (unsigned int i = 0; i < pgonZ.size(); i++)
533  edm::LogVerbatim("HCalGeom") << "\t"
534  << "\tZ = " << pgonZ[i] << "\tRmin = " << pgonRmin[i] << "\tRmax = " << pgonRmax[i];
535 #endif
536  } else {
538  DDName(idName, idNameSpace), getNsectortot(), -alpha, dphi, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
539 #ifdef EDM_ML_DEBUG
540  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << DDName(idName, idNameSpace) << " Polyhedra made of "
541  << getGenMaterial() << " with " << getNsectortot() << " sectors from "
542  << convertRadToDeg(-alpha) << " to " << convertRadToDeg(-alpha + dphi) << " and with "
543  << nsec << " sections ";
544  for (unsigned int i = 0; i < pgonZHalf.size(); i++)
545  edm::LogVerbatim("HCalGeom") << "\t"
546  << "\tZ = " << pgonZHalf[i] << "\tRmin = " << pgonRminHalf[i]
547  << "\tRmax = " << pgonRmaxHalf[i];
548 #endif
549  }
550 
552  DDMaterial matter(matname);
553  DDLogicalPart genlogic(DDName(idName, idNameSpace), matter, solid);
554 
555  DDName parentName = parent().name();
556  DDTranslation r0(0, 0, 0);
557  cpv.position(DDName(idName, idNameSpace), parentName, 1, r0, rot);
558 #ifdef EDM_ML_DEBUG
559  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << DDName(idName, idNameSpace) << " number 1 positioned in "
560  << parentName << " at (0, 0, 0) with no rotation";
561 #endif
562  //Forward and backwards halfs
563  name = idName + "Half";
564  nf = (ntot + 1) / 2;
565 #ifdef EDM_ML_DEBUG
566  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << DDName(name, idNameSpace) << " Polyhedra made of "
567  << getGenMaterial() << " with " << getNsectortot() << " sectors from "
568  << convertRadToDeg(-alpha) << " to " << convertRadToDeg(-alpha + dphi) << " and with "
569  << nf << " sections ";
570  for (unsigned int i = 0; i < pgonZHalf.size(); i++)
571  edm::LogVerbatim("HCalGeom") << "\t"
572  << "\tZ = " << pgonZHalf[i] << "\tRmin = " << pgonRminHalf[i]
573  << "\tRmax = " << pgonRmaxHalf[i];
574 #endif
575 
577  DDName(name, idNameSpace), getNsectortot(), -alpha, dphi, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
578  DDLogicalPart genlogich(DDName(name, idNameSpace), matter, solid);
579 
580  cpv.position(genlogich, genlogic, 1, r0, rot);
581 #ifdef EDM_ML_DEBUG
582  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << genlogich.name() << " number 1 positioned in "
583  << genlogic.name() << " at (0, 0, 0) with no rotation";
584 #endif
585  if (getNhalf() != 1) {
586  rot = DDRotation(DDName(rotHalf, rotns));
587  cpv.position(genlogich, genlogic, 2, r0, rot);
588 #ifdef EDM_ML_DEBUG
589  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << genlogich.name() << " number 2 positioned in "
590  << genlogic.name() << " at " << r0 << " with " << rot;
591 #endif
592  } //end if (getNhalf...
593 
594  //Construct sector (from -alpha to +alpha)
595  name = idName + "Module";
596 #ifdef EDM_ML_DEBUG
597  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << DDName(name, idNameSpace) << " Polyhedra made of "
598  << getGenMaterial() << " with 1 sector from " << convertRadToDeg(-alpha) << " to "
599  << convertRadToDeg(alpha) << " and with " << nf << " sections";
600  for (unsigned int i = 0; i < pgonZHalf.size(); i++)
601  edm::LogVerbatim("HCalGeom") << "\t"
602  << "\tZ = " << pgonZHalf[i] << "\tRmin = " << pgonRminHalf[i]
603  << "\tRmax = " << pgonRmaxHalf[i];
604 #endif
605 
606  solid =
607  DDSolidFactory::polyhedra(DDName(name, idNameSpace), 1, -alpha, 2 * alpha, pgonZHalf, pgonRminHalf, pgonRmaxHalf);
608  DDLogicalPart seclogic(DDName(name, idNameSpace), matter, solid);
609 
610  double theta = 90._deg;
611  for (int ii = 0; ii < getNsectortot(); ii++) {
612  double phi = ii * 2 * alpha;
613  double phiy = phi + 90._deg;
614 
616  std::string rotstr("NULL");
617  if (phi != 0) {
618  rotstr = "R" + formatAsDegreesInInteger(phi);
619  rotation = DDRotation(DDName(rotstr, rotns));
620  if (!rotation) {
621 #ifdef EDM_ML_DEBUG
622  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Creating a new "
623  << "rotation " << rotstr << "\t 90," << convertRadToDeg(phi) << ",90,"
624  << (90 + convertRadToDeg(phi)) << ", 0, 0";
625 #endif
626  rotation = DDrot(DDName(rotstr, rotns), theta, phi, theta, phiy, 0, 0);
627  } //if !rotation
628  } //if phideg!=0
629 
630  cpv.position(seclogic, genlogich, ii + 1, r0, rotation);
631 #ifdef EDM_ML_DEBUG
632  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << seclogic.name() << " number " << ii + 1 << " positioned in "
633  << genlogich.name() << " at " << r0 << " with " << rotation;
634 #endif
635  }
636 
637  //Construct the things inside the sector
638  constructInsideSector(seclogic, cpv);
639 }
640 
642 #ifdef EDM_ML_DEBUG
643  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: Layers (" << getNLayers() << ") ...";
644 #endif
645  double alpha = (1._pi) / getNsectors();
646  double rin = getRin();
647  for (int i = 0; i < getNLayers(); i++) {
650  DDSplit(getLayerMaterial(i)).second); //idNameSpace);
651  DDMaterial matter(matname);
652 
653  double width = getLayerWidth(i);
654  double rout = rin + width;
655 
656  int in = 0, out = 0;
657  for (int j = 0; j < getRzones() - 1; j++) {
658  if (rin >= getRmax(j))
659  in = j + 1;
660  if (rout > getRmax(j))
661  out = j + 1;
662  }
663  double zout = getZoff(in) + rin * getTanTheta(in);
664 
666  //vertical walls are allowed in SolidPolyhedra
667  double deltaz = 0;
668  int nsec = 2;
669  std::vector<double> pgonZ, pgonRmin, pgonRmax;
670  // index 0
671  pgonZ.emplace_back(0);
672  pgonRmin.emplace_back(rin);
673  pgonRmax.emplace_back(rout);
674  // index 1
675  pgonZ.emplace_back(zout);
676  pgonRmin.emplace_back(rin);
677  pgonRmax.emplace_back(rout);
678  if (in == out) {
679  if (in <= 3) {
680  //index 2
681  pgonZ.emplace_back(getZoff(in) + rout * getTanTheta(in));
682  pgonRmin.emplace_back(pgonRmax[1]);
683  pgonRmax.emplace_back(pgonRmax[1]);
684  nsec++;
685  }
686  } else {
687  if (in == 3) {
688  //redo index 1, add index 2
689  pgonZ[1] = (getZoff(out) + getRmax(out) * getTanTheta(out));
690  pgonZ.emplace_back(pgonZ[1] + deltaz);
691  pgonRmin.emplace_back(pgonRmin[1]);
692  pgonRmax.emplace_back(getRmax(in));
693  //index 3
694  pgonZ.emplace_back(getZoff(in) + getRmax(in) * getTanTheta(in));
695  pgonRmin.emplace_back(pgonRmin[2]);
696  pgonRmax.emplace_back(pgonRmax[2]);
697  nsec += 2;
698  } else {
699  //index 2
700  pgonZ.emplace_back(getZoff(in) + getRmax(in) * getTanTheta(in));
701  pgonRmin.emplace_back(getRmax(in));
702  pgonRmax.emplace_back(pgonRmax[1]);
703  nsec++;
704  if (in == 0) {
705  pgonZ.emplace_back(getZoff(out) + getRmax(in) * getTanTheta(out));
706  pgonRmin.emplace_back(pgonRmin[2]);
707  pgonRmax.emplace_back(pgonRmax[2]);
708  nsec++;
709  }
710  if (in <= 1) {
711  pgonZ.emplace_back(getZoff(out) + rout * getTanTheta(out));
712  pgonRmin.emplace_back(rout);
713  pgonRmax.emplace_back(rout);
714  nsec++;
715  }
716  }
717  }
718  //Solid & volume
719  DDSolid solid;
720  double alpha1 = alpha;
721  if (getLayerGap(i) > 1.e-6) {
722  double rmid = 0.5 * (rin + rout);
723  double width = rmid * tan(alpha) - getLayerGap(i);
724  alpha1 = atan(width / rmid);
725 #ifdef EDM_ML_DEBUG
726  edm::LogVerbatim("HCalGeom") << "\t"
727  << "Alpha_1 modified from " << convertRadToDeg(alpha) << " to "
728  << convertRadToDeg(alpha1) << " Rmid " << rmid << " Reduced width " << width;
729 #endif
730  }
731 #ifdef EDM_ML_DEBUG
732  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << name << " (Layer " << i << ") Polyhedra made of "
733  << getLayerMaterial(i) << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
734  << convertRadToDeg(alpha1) << " and with " << nsec << " sections";
735  for (unsigned int k = 0; k < pgonZ.size(); k++)
736  edm::LogVerbatim("HCalGeom") << "\t"
737  << "\t" << pgonZ[k] << "\t" << pgonRmin[k] << "\t" << pgonRmax[k];
738 #endif
739  solid = DDSolidFactory::polyhedra(DDName(name, idNameSpace), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
740  DDLogicalPart glog(DDName(name, idNameSpace), matter, solid);
741 
742  cpv.position(glog, sector, getLayerId(i), DDTranslation(0.0, 0.0, 0.0), DDRotation());
743 #ifdef EDM_ML_DEBUG
744  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " number " << getLayerId(i)
745  << " positioned in " << sector.name() << " at (0,0,0) with no rotation";
746 #endif
748  getLayerLabel(i),
749  getLayerId(i),
750  getLayerAbsorb(i),
751  rin,
752  getLayerD1(i),
753  alpha1,
754  getLayerD2(i),
755  getLayerAlpha(i),
756  getLayerT1(i),
757  getLayerT2(i),
758  cpv);
759  rin = rout;
760  }
761 }
762 
764  const std::string& nm,
765  int id,
766  int nAbs,
767  double rin,
768  double d1,
769  double alpha1,
770  double d2,
771  double alpha2,
772  double t1,
773  double t2,
774  DDCompactView& cpv) {
775 #ifdef EDM_ML_DEBUG
776  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: \t\tInside layer " << id << "...";
777 #endif
778  //Pointers to the Rotation Matrices and to the Materials
781 
782  std::string nam0 = nm + "In";
783  std::string name = idName + nam0;
785  DDMaterial matter(matName);
786 
787  DDSolid solid;
788  DDLogicalPart glog, mother;
789  double rsi, dx, dy, dz, x, y;
790  int i, in;
791  //Two lower volumes
792  if (alpha1 > 0) {
793  rsi = rin + d1;
794  in = 0;
795  for (i = 0; i < getRzones() - 1; i++) {
796  if (rsi >= getRmax(i))
797  in = i + 1;
798  }
799  dx = 0.5 * t1;
800  dy = 0.5 * rsi * (tan(alpha1) - tan(alpha2));
801  dz = 0.5 * (getZoff(in) + rsi * getTanTheta(in));
802  x = rsi + dx;
803  y = 0.5 * rsi * (tan(alpha1) + tan(alpha2));
804  DDTranslation r11(x, y, dz);
805  DDTranslation r12(x, -y, dz);
806 
807  solid = DDSolidFactory::box(DDName(name + "1", idNameSpace), dx, dy, dz);
808 #ifdef EDM_ML_DEBUG
809  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << getDetMat()
810  << " of dimensions " << dx << ", " << dy << ", " << dz;
811 #endif
812  glog = DDLogicalPart(solid.ddname(), matter, solid);
813 
814  if (nAbs != 0) {
815  mother = constructSideLayer(laylog, name, nAbs, rin, alpha1, cpv);
816  } else {
817  mother = laylog;
818  }
819  cpv.position(glog, mother, idOffset + 1, r11, DDRotation());
820  cpv.position(glog, mother, idOffset + 2, r12, rot);
821 #ifdef EDM_ML_DEBUG
822  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << idOffset + 1
823  << " positioned in " << mother.name() << " at " << r11 << " with no rotation\n"
824  << "DDHCalBarrelAlgo: " << glog.name() << " Number " << idOffset + 2
825  << " positioned in " << mother.name() << " at " << r12 << " with " << rot;
826 #endif
827  //Constructin the plastics and scintillators inside
828  constructInsideDetectors(glog, nam0 + "1", id, dx, dy, dz, 1, cpv);
829  }
830 
831  //Upper volume
832  rsi = rin + d2;
833  in = 0;
834  for (i = 0; i < getRzones() - 1; i++) {
835  if (rsi >= getRmax(i))
836  in = i + 1;
837  }
838  dx = 0.5 * t2;
839  dy = 0.5 * rsi * tan(alpha2);
840  dz = 0.5 * (getZoff(in) + rsi * getTanTheta(in));
841  x = rsi + dx;
842  DDTranslation r21(x, dy, dz);
843  DDTranslation r22(x, -dy, dz);
844 
845  solid = DDSolidFactory::box(DDName(name + "2", idNameSpace), dx, dy, dz);
846 #ifdef EDM_ML_DEBUG
847  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << getDetMat()
848  << " of dimensions " << dx << ", " << dy << ", " << dz;
849 #endif
850  glog = DDLogicalPart(solid.ddname(), matter, solid);
851 
852  if (nAbs < 0) {
853  mother = constructMidLayer(laylog, name, rin, alpha1, cpv);
854  } else {
855  mother = laylog;
856  }
857  cpv.position(glog, mother, idOffset + 3, r21, DDRotation());
858  cpv.position(glog, mother, idOffset + 4, r22, rot);
859 #ifdef EDM_ML_DEBUG
860  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << idOffset + 3 << " positioned in "
861  << mother.name() << " at " << r21
862  << " with no rotation\nDDHCalBarrelAlgo: " << glog.name() << " Number " << idOffset + 4
863  << " positioned in " << mother.name() << " at " << r22 << " with " << rot;
864 #endif
865  //Constructin the plastics and scintillators inside
866  constructInsideDetectors(glog, nam0 + "2", id, dx, dy, dz, 2, cpv);
867 }
868 
870  const DDLogicalPart& laylog, const std::string& nm, int nAbs, double rin, double alpha, DDCompactView& cpv) {
871  //Extra absorber layer
872  int k = abs(nAbs) - 1;
873  std::string namek = nm + "Side";
874  double rsi = rin + getSideD(k);
875  int in = 0;
876  for (int i = 0; i < getRzones() - 1; i++) {
877  if (rsi >= getRmax(i))
878  in = i + 1;
879  }
880  std::vector<double> pgonZ, pgonRmin, pgonRmax;
881  // index 0
882  pgonZ.emplace_back(0.0);
883  pgonRmin.emplace_back(rsi);
884  pgonRmax.emplace_back(rsi + getSideT(k));
885  // index 1
886  pgonZ.emplace_back(getZoff(in) + rsi * getTanTheta(in));
887  pgonRmin.emplace_back(rsi);
888  pgonRmax.emplace_back(pgonRmax[0]);
889  // index 2
890  pgonZ.emplace_back(getZoff(in) + pgonRmax[0] * getTanTheta(in));
891  pgonRmin.emplace_back(pgonRmax[1]);
892  pgonRmax.emplace_back(pgonRmax[1]);
893  DDSolid solid =
894  DDSolidFactory::polyhedra(DDName(namek, idNameSpace), 1, -alpha, 2 * alpha, pgonZ, pgonRmin, pgonRmax);
895 #ifdef EDM_ML_DEBUG
896  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << getSideMat(k)
897  << " with 1 sector from " << convertRadToDeg(-alpha) << " to " << convertRadToDeg(alpha)
898  << " and with " << pgonZ.size() << " sections";
899  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
900  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
901  << "\tRmax = " << pgonRmax[ii];
902 #endif
903 
905  DDMaterial matter(matName);
906  DDLogicalPart glog = DDLogicalPart(solid.ddname(), matter, solid);
907 
908  cpv.position(glog, laylog, 1, DDTranslation(), DDRotation());
909 #ifdef EDM_ML_DEBUG
910  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << laylog.name()
911  << " at (0,0,0) with no rotation";
912 #endif
913  if (nAbs < 0) {
914  DDLogicalPart mother = glog;
915  double rmid = pgonRmax[0];
916  for (int i = 0; i < getSideAbsorber(); i++) {
917  double alpha1 = atan(getSideAbsW(i) / rmid);
918  if (alpha1 > 0) {
919  std::string name = namek + getSideAbsName(i);
920  solid = DDSolidFactory::polyhedra(DDName(name, idNameSpace), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
921 #ifdef EDM_ML_DEBUG
922  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of "
923  << getSideAbsMat(i) << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
924  << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
925  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
926  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
927  << "\tRmax = " << pgonRmax[ii];
928 #endif
929 
930  DDName matName(DDSplit(getSideAbsMat(i)).first, DDSplit(getSideAbsMat(i)).second);
931  DDMaterial matter(matName);
932  DDLogicalPart log = DDLogicalPart(solid.ddname(), matter, solid);
933 
934  cpv.position(log, mother, 1, DDTranslation(), DDRotation());
935 #ifdef EDM_ML_DEBUG
936  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in "
937  << mother.name() << " at (0,0,0) with no rotation";
938 #endif
939  mother = log;
940  }
941  }
942  }
943  return glog;
944 }
945 
947  const DDLogicalPart& laylog, const std::string& nm, double rin, double alpha, DDCompactView& cpv) {
948  DDSolid solid;
949  DDLogicalPart log, glog;
950  std::string name = nm + "Mid";
951  for (int k = 0; k < getAbsorberN(); k++) {
952  std::string namek = name + getAbsorbName(k);
953  double rsi = rin + getAbsorbD(k);
954  int in = 0;
955  for (int i = 0; i < getRzones() - 1; i++) {
956  if (rsi >= getRmax(i))
957  in = i + 1;
958  }
959  std::vector<double> pgonZ, pgonRmin, pgonRmax;
960  // index 0
961  pgonZ.emplace_back(0.0);
962  pgonRmin.emplace_back(rsi);
963  pgonRmax.emplace_back(rsi + getAbsorbT(k));
964  // index 1
965  pgonZ.emplace_back(getZoff(in) + rsi * getTanTheta(in));
966  pgonRmin.emplace_back(rsi);
967  pgonRmax.emplace_back(pgonRmax[0]);
968  // index 2
969  pgonZ.emplace_back(getZoff(in) + pgonRmax[0] * getTanTheta(in));
970  pgonRmin.emplace_back(pgonRmax[1]);
971  pgonRmax.emplace_back(pgonRmax[1]);
972  solid = DDSolidFactory::polyhedra(DDName(namek, idNameSpace), 1, -alpha, 2 * alpha, pgonZ, pgonRmin, pgonRmax);
973 #ifdef EDM_ML_DEBUG
974  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << getAbsorbMat(k)
975  << " with 1 sector from " << convertRadToDeg(-alpha) << " to "
976  << convertRadToDeg(alpha) << " and with " << pgonZ.size() << " sections";
977  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
978  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
979  << "\tRmax = " << pgonRmax[ii];
980 #endif
981 
983  DDMaterial matter(matName);
984  log = DDLogicalPart(solid.ddname(), matter, solid);
985 
986  cpv.position(log, laylog, 1, DDTranslation(), DDRotation());
987 #ifdef EDM_ML_DEBUG
988  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in " << laylog.name()
989  << " at (0,0,0) with no rotation";
990 #endif
991  if (k == 0) {
992  double rmin = pgonRmin[0];
993  double rmax = pgonRmax[0];
994  DDLogicalPart mother = log;
995  for (int i = 0; i < 1; i++) {
996  double alpha1 = atan(getMidAbsW(i) / rmin);
997  std::string namek = name + getMidAbsName(i);
998  solid =
999  DDSolidFactory::polyhedra(DDName(namek, idNameSpace), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
1000 #ifdef EDM_ML_DEBUG
1001  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << getMidAbsMat(i)
1002  << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
1003  << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
1004  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
1005  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
1006  << "\tRmax = " << pgonRmax[ii];
1007 #endif
1008 
1009  DDName matNam1(DDSplit(getMidAbsMat(i)).first, DDSplit(getMidAbsMat(i)).second);
1010  DDMaterial matter1(matNam1);
1011  log = DDLogicalPart(solid.ddname(), matter1, solid);
1012 
1013  cpv.position(log, mother, 1, DDTranslation(), DDRotation());
1014 #ifdef EDM_ML_DEBUG
1015  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number 1 positioned in "
1016  << mother.name() << " at (0,0,0) with no rotation";
1017 #endif
1018  mother = log;
1019  }
1020 
1021  // Now the layer with detectors
1022  double rmid = rmin + getMiddleD();
1023  pgonRmin[0] = rmid;
1024  pgonRmax[0] = rmax;
1025  pgonRmin[1] = rmid;
1026  pgonRmax[1] = rmax;
1027  pgonZ[1] = getZoff(in) + rmid * getTanTheta(in);
1028  pgonRmin[2] = rmax;
1029  pgonRmax[2] = rmax;
1030  pgonZ[2] = getZoff(in) + rmax * getTanTheta(in);
1031  double alpha1 = atan(getMiddleW() / rmin);
1032  solid = DDSolidFactory::polyhedra(DDName(name, idNameSpace), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
1033 #ifdef EDM_ML_DEBUG
1034  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << getMiddleMat()
1035  << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
1036  << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
1037  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
1038  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
1039  << "\tRmax = " << pgonRmax[ii];
1040 #endif
1041 
1042  DDName matNam1(DDSplit(getMiddleMat()).first, DDSplit(getMiddleMat()).second);
1043  DDMaterial matter1(matNam1);
1044  glog = DDLogicalPart(solid.ddname(), matter1, solid);
1045 
1046  cpv.position(glog, mother, 1, DDTranslation(), DDRotation());
1047 #ifdef EDM_ML_DEBUG
1048  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << mother.name()
1049  << " at (0,0,0) with no rotation";
1050 #endif
1051  // Now the remaining absorber layers
1052  for (int i = 1; i < getMidAbsorber(); i++) {
1053  namek = name + getMidAbsName(i);
1054  rmid = rmin + getMidAbsT(i);
1055  pgonRmin[0] = rmin;
1056  pgonRmax[0] = rmid;
1057  pgonRmin[1] = rmin;
1058  pgonRmax[1] = rmid;
1059  pgonZ[1] = getZoff(in) + rmin * getTanTheta(in);
1060  pgonRmin[2] = rmid;
1061  pgonRmax[2] = rmid;
1062  pgonZ[2] = getZoff(in) + rmid * getTanTheta(in);
1063  alpha1 = atan(getMidAbsW(i) / rmin);
1064  solid =
1065  DDSolidFactory::polyhedra(DDName(namek, idNameSpace), 1, -alpha1, 2 * alpha1, pgonZ, pgonRmin, pgonRmax);
1066 #ifdef EDM_ML_DEBUG
1067  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Polyhedra made of " << getMidAbsMat(i)
1068  << " with 1 sector from " << convertRadToDeg(-alpha1) << " to "
1069  << convertRadToDeg(alpha1) << " and with " << pgonZ.size() << " sections";
1070  for (unsigned int ii = 0; ii < pgonZ.size(); ii++)
1071  edm::LogVerbatim("HCalGeom") << "\t\tZ = " << pgonZ[ii] << "\tRmin = " << pgonRmin[ii]
1072  << "\tRmax = " << pgonRmax[ii];
1073 #endif
1074 
1075  DDName matName2(DDSplit(getMidAbsMat(i)).first, DDSplit(getMidAbsMat(i)).second);
1076  DDMaterial matter2(matName2);
1077  log = DDLogicalPart(solid.ddname(), matter2, solid);
1078 
1079  cpv.position(log, mother, i, DDTranslation(), DDRotation());
1080 #ifdef EDM_ML_DEBUG
1081  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << log.name() << " Number " << i << " positioned in "
1082  << mother.name() << " at (0,0,0) with no "
1083  << "rotation";
1084 #endif
1085  mother = log;
1086  }
1087  }
1088  }
1089  return glog;
1090 }
1091 
1093  const std::string& name,
1094  int id,
1095  double dx,
1096  double dy,
1097  double dz,
1098  int type,
1099  DDCompactView& cpv) {
1100 #ifdef EDM_ML_DEBUG
1101  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: \t\tInside detector " << id << "...";
1102 #endif
1103 
1105  DDMaterial plmatter(plmatname);
1107  DDMaterial scmatter(scmatname);
1108 
1109  std::string plname = detector.name().name() + "Plastic_";
1110  std::string scname = idName + "Scintillator" + name;
1111 
1112  id--;
1113  DDSolid solid;
1114  DDLogicalPart glog;
1115  double wid, y = 0;
1116  double dx1, dx2, shiftX;
1117 
1118  if (type == 1) {
1119  wid = 0.5 * getDetWidth1(id);
1120  if (getDetPosY(id) > 0)
1121  y = -dy + wid;
1122  dx1 = 0.5 * getDetT11(id);
1123  dx2 = 0.5 * getDetT21(id);
1124  shiftX = getDetdP1(id);
1125  } else {
1126  wid = 0.5 * getDetWidth2(id);
1127  dx1 = 0.5 * getDetT12(id);
1128  dx2 = 0.5 * getDetT22(id);
1129  shiftX = getDetdP2(id);
1130  }
1131 
1132  solid = DDSolidFactory::box(DDName(plname + "1", idNameSpace), dx1, wid, dz);
1133 #ifdef EDM_ML_DEBUG
1134  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << getDetMatPl()
1135  << " of dimensions " << dx1 << ", " << wid << ", " << dz;
1136 #endif
1137  glog = DDLogicalPart(solid.ddname(), plmatter, solid);
1138 
1139  double x = shiftX + dx1 - dx;
1140  cpv.position(glog, detector, 1, DDTranslation(x, y, 0), DDRotation());
1141 #ifdef EDM_ML_DEBUG
1142  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << detector.name()
1143  << " at (" << x << "," << y << ",0) with no rotation";
1144 #endif
1145  solid = DDSolidFactory::box(DDName(scname, idNameSpace), 0.5 * getDetTsc(id), wid, dz);
1146 #ifdef EDM_ML_DEBUG
1147  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << getDetMatSc()
1148  << " of dimensions " << 0.5 * getDetTsc(id) << ", " << wid << ", " << dz;
1149 #endif
1150  glog = DDLogicalPart(solid.ddname(), scmatter, solid);
1151 
1152  x += dx1 + 0.5 * getDetTsc(id);
1153  int copyNo = id * 10 + getDetType(id);
1154  cpv.position(glog, detector, copyNo, DDTranslation(x, y, 0), DDRotation());
1155 #ifdef EDM_ML_DEBUG
1156  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number " << copyNo << " positioned in "
1157  << detector.name() << " at (" << x << "," << y << ",0) with no rotation";
1158 #endif
1159  solid = DDSolidFactory::box(DDName(plname + "2", idNameSpace), dx2, wid, dz);
1160 #ifdef EDM_ML_DEBUG
1161  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << solid.name() << " Box made of " << getDetMatPl()
1162  << " of dimensions " << dx2 << ", " << wid << ", " << dz;
1163 #endif
1164  glog = DDLogicalPart(solid.ddname(), plmatter, solid);
1165 
1166  x += 0.5 * getDetTsc(id) + dx2;
1167  cpv.position(glog, detector, 1, DDTranslation(x, y, 0), DDRotation());
1168 #ifdef EDM_ML_DEBUG
1169  edm::LogVerbatim("HCalGeom") << "DDHCalBarrelAlgo: " << glog.name() << " Number 1 positioned in " << detector.name()
1170  << " at (" << x << "," << y << ",0) with no rotation";
1171 #endif
1172 }
1173 
1174 DEFINE_EDM_PLUGIN(DDAlgorithmFactory, DDHCalBarrelAlgo, "hcal:DDHCalBarrelAlgo");
std::vector< std::string > absorbMat
Log< level::Info, true > LogVerbatim
static AlgebraicMatrix initialize()
const std::string & getDetMat() const
int getNsectortot() const
float alpha
Definition: AMPTWrapper.h:105
static std::vector< std::string > checklist log
double getRin() const
const N & name() const
Definition: DDBase.h:59
std::vector< std::string > midMat
const std::string & getAbsorbName(unsigned int i) const
double getLayerD1(unsigned i) const
double getAbsorbT(unsigned int i) const
std::vector< double > sideAbsW
int getNhalf() const
void position(const DDLogicalPart &self, const DDLogicalPart &parent, const std::string &copyno, const DDTranslation &trans, const DDRotation &rot, const DDDivision *div=nullptr)
std::vector< std::string > layerLabel
const std::string & getSideMat(unsigned int i) const
const std::string & getAbsorbMat(unsigned int i) const
double getRmax(unsigned int i) const
DDMaterial is used to define and access material information.
Definition: DDMaterial.h:45
int getDetPosY(unsigned int i) const
double getLayerD2(unsigned i) const
double getSideAbsW(unsigned int i) const
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
std::string idNameSpace
std::vector< double > layerD2
Geom::Theta< T > theta() const
std::vector< std::string > layerMat
DDName is used to identify DDD entities uniquely.
Definition: DDName.h:17
double getLayerAlpha(unsigned i) const
std::vector< int > layerId
static std::string & ns()
std::vector< double > rmax
int getRzones() const
void constructInsideDetectors(const DDLogicalPart &detector, const std::string &name, int id, double dx, double dy, double dz, int type, DDCompactView &cpv)
int ii
Definition: cuy.py:589
double getMidAbsT(unsigned int i) const
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:81
std::string formatAsDegreesInInteger(double radianVal)
Definition: DDTypes.cc:79
A DDSolid represents the shape of a part.
Definition: DDSolid.h:39
void constructInsideSector(const DDLogicalPart &sector, DDCompactView &cpv)
const std::string & getGenMaterial() const
double getRout() const
Represents a uniquely identifyable rotation matrix.
Definition: DDTransform.h:57
double getZoff(unsigned int i) const
std::string genMaterial
U second(std::pair< T, U > const &p)
std::vector< double > layerT1
std::vector< double > detdP2
std::vector< double > detWidth2
double getDetT22(unsigned int i) const
const std::string & getDetMatSc() const
const std::string & getMidAbsName(unsigned int i) const
int getAbsorberN() const
double getLayerGap(unsigned int i) const
std::vector< std::string > sideAbsName
std::vector< double > detdP1
std::vector< double > absorbT
void execute(DDCompactView &cpv) override
std::vector< double > detWidth1
DDLogicalPart constructSideLayer(const DDLogicalPart &laylog, const std::string &nm, int nAbs, double rin, double alpha, DDCompactView &cpv)
double getDetT11(unsigned int i) const
int getNLayers() const
void constructGeneralVolume(DDCompactView &cpv)
int getLayerAbsorb(unsigned int i) const
double getLayerT2(unsigned i) const
int getDetType(unsigned int i) const
double getDetdP1(unsigned int i) const
double getLayerT1(unsigned i) const
int getSideAbsorber() const
double getSideD(unsigned int i) const
double getMidAbsW(unsigned int i) const
std::vector< double > layerGap
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< double > layerAlpha
void constructInsideLayers(const DDLogicalPart &laylog, const std::string &name, int id, int nAbs, double rin, double d1, double alpha1, double d2, double alpha2, double t1, double t2, DDCompactView &cpv)
A DDLogicalPart aggregates information concerning material, solid and sensitveness ...
Definition: DDLogicalPart.h:93
const std::string & getLayerMaterial(unsigned i) const
double getTanTheta(unsigned int i) const
std::vector< double > detT22
std::vector< int > layerAbsorb
DDRotation DDrot(const DDName &name, std::unique_ptr< DDRotationMatrix > rot)
Definition of a uniquely identifiable rotation matrix named by DDName name.
Definition: DDRotation.cc:67
int getNsectors() const
std::vector< std::string > sideAbsMat
std::vector< double > sideT
std::vector< double > detTsc
~DDHCalBarrelAlgo() override
double getDetTsc(unsigned int i) const
static DDSolid box(const DDName &name, double xHalf, double yHalf, double zHalf)
Creates a box with side length 2*xHalf, 2*yHalf, 2*zHalf.
Definition: DDSolid.cc:547
int getLayerId(unsigned i) const
double getMiddleW() const
double getTheta(unsigned int i) const
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
double getDetdP2(unsigned int i) const
std::vector< std::string > midName
const std::string & getDetMatPl() const
std::vector< std::string > sideMat
std::vector< double > detT21
std::vector< int > detPosY
const std::string & getLayerLabel(unsigned i) const
double getLayerWidth(unsigned i) const
std::vector< int > detType
double getAbsorbD(unsigned int i) const
double getDetT21(unsigned int i) const
std::vector< double > zoff
const std::string & getMidAbsMat(unsigned int i) const
const std::string & getMiddleMat() const
std::vector< double > sideD
const std::string & getSideAbsName(unsigned int i) const
std::vector< double > detT11
double getMiddleD() const
#define DEFINE_EDM_PLUGIN(factory, type, name)
std::vector< double > absorbD
std::vector< double > ttheta
double getDetWidth2(unsigned int i) const
std::vector< double > theta
int getMidAbsorber() const
std::vector< double > layerWidth
std::vector< std::string > absorbName
double getDetT12(unsigned int i) const
const std::string & getSideAbsMat(unsigned int i) const
static constexpr float d1
std::vector< double > midW
double getDetWidth1(unsigned int i) const
std::vector< double > layerD1
double getSideT(unsigned int i) const
std::pair< std::string, std::string > DDSplit(const std::string &n)
split into (name,namespace), separator = &#39;:&#39;
Definition: DDSplit.cc:3
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DDTranslation
Definition: DDTranslation.h:7
std::vector< double > detT12
const std::string & name() const
Returns the name.
Definition: DDName.cc:41
__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< double > midT
std::vector< double > layerT2
DDLogicalPart constructMidLayer(const DDLogicalPart &laylog, const std::string &nm, double rin, double alpha, DDCompactView &cpv)
void initialize(const DDNumericArguments &nArgs, const DDVectorArguments &vArgs, const DDMapArguments &mArgs, const DDStringArguments &sArgs, const DDStringVectorArguments &vsArgs) override
const N & ddname() const
Definition: DDBase.h:61