Belle II Software  release-05-01-25
backward.cc
1 #include <ecl/geometry/GeoECLCreator.h>
2 #include "ecl/geometry/BelleLathe.h"
3 #include "ecl/geometry/BelleCrystal.h"
4 #include "G4LogicalVolume.hh"
5 #include "G4PVPlacement.hh"
6 #include <G4VisAttributes.hh>
7 #include <G4Box.hh>
8 #include <G4AssemblyVolume.hh>
9 #include <G4Region.hh>
10 #include <G4TwoVector.hh>
11 #include <G4PVReplica.hh>
12 #include "G4UserLimits.hh"
13 #include <G4Trap.hh>
14 
15 #include <iostream>
16 #include "CLHEP/Matrix/Vector.h"
17 #include "ecl/geometry/shapes.h"
18 #include <geometry/Materials.h>
19 
20 using namespace std;
21 using namespace Belle2;
22 using namespace Belle2::geometry;
23 
24 void Belle2::ECL::GeoECLCreator::backward(G4LogicalVolume& _top)
25 {
26  G4LogicalVolume* top = &_top;
27 
28  bool sec = 0;
29  double phi0 = 0, dphi = (sec) ? M_PI / 16 : 2 * M_PI;
30 
31  bool b_inner_support_ring = 1;
32  bool b_outer_support_ring = 1;
33  bool b_support_wall = 1;
34  bool b_ribs = 1;
35  bool b_septum_wall = 1;
36  bool b_crystals = 1;
37  bool b_preamplifier = 1;
38  bool b_support_leg = 1;
39  int overlap = m_overlap;
40 
41  int npoints = 1000 * 1000;
42 
43  // vector<cplacement_t> bp = load_placements("/ecl/data/crystal_placement_backward.dat");
44  vector<cplacement_t> bp = load_placements(m_sap, ECLParts::backward);
45  vector<cplacement_t>::iterator fp = find_if(bp.begin(), bp.end(), [](const cplacement_t& p) {
46  const int ECL_backward_part = 1000; // lookup part
47  return p.nshape == ECL_backward_part;
48  });
49  // global transformation before placing the whole backward part in the top logical volume
50  G4Transform3D gTrans = (fp == bp.end()) ? G4Translate3D(0, 0, -1020) : get_transform(*fp);
51  if (fp != bp.end()) bp.erase(fp); // now not needed
52 
53  // cppcheck-suppress knownConditionTrueFalse
54  if (b_inner_support_ring) {
55  zr_t vc1[] = {{0., 452.3 + 3}, {0., 452.3}, {3., 474.9 - 20 / cosd(27.81)}, {434., 702.27 - 20 / cosd(27.81)}, {434., 702.27}, {3., 474.9}, {3., 452.3 + 3}};
56  std::vector<zr_t> contour1(vc1, vc1 + sizeof(vc1) / sizeof(zr_t));
57  G4VSolid* part1solid = new BelleLathe("part1solid", phi0, dphi, contour1);
58  G4LogicalVolume* part1logical = new G4LogicalVolume(part1solid, Materials::get("SUS304"), "part1logical", 0, 0, 0);
59  part1logical->SetVisAttributes(att("iron"));
60  auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part1logical, "part1physical", top, false, 0, 0);
61  if (overlap) pv->CheckOverlaps(npoints);
62  }
63 
64  // cppcheck-suppress knownConditionTrueFalse
65  if (b_support_wall) {
66  // solveing equation to get L : 3+L*cosd(52.90)+1.6*cosd(52.90+90) = 435 - 202.67
67  double L = (435 - 202.67 - 3 - 1.6 * cosd(52.90 + 90)) / cosd(52.90);
68  zr_t vc23[] = {{0, 452.3 + 3}, {3, 452.3 + 3}, {3, 1190.2}, {3 + L * cosd(52.90), 1190.2 + L * sind(52.90)},
69  {3 + L * cosd(52.90) + 1.6 * cosd(52.90 + 90), 1190.2 + L * sind(52.90) + 1.6 * sind(52.90 + 90)}, {3 + 1.6 * cosd(52.90 + 90), 1190.2 + 1.6 * sind(52.90 + 90)}, {0, 1190.2}
70  };
71  std::vector<zr_t> contour23(vc23, vc23 + sizeof(vc23) / sizeof(zr_t));
72  G4VSolid* part23solid = new BelleLathe("part23solid", phi0, dphi, contour23);
73  G4LogicalVolume* part23logical = new G4LogicalVolume(part23solid, Materials::get("A5052"), "part23logical", 0, 0, 0);
74  part23logical->SetVisAttributes(att("alum"));
75  auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part23logical, "part23physical", top, false, 0, 0);
76  if (overlap) pv->CheckOverlaps(npoints);
77  }
78 
79  // cppcheck-suppress knownConditionTrueFalse
80  if (b_outer_support_ring) {
81  zr_t vc4[] = {{434 - 214.8, 1496 - 20}, {434, 1496 - 20}, {434, 1496 - 5}, {434 + 5, 1496 - 5}, {434 + 5, 1496}, {434 - 199.66, 1496}};
82  std::vector<zr_t> contour4(vc4, vc4 + sizeof(vc4) / sizeof(zr_t));
83  G4VSolid* part4solid = new BelleLathe("part4solid", phi0, dphi, contour4);
84  G4LogicalVolume* part4logical = new G4LogicalVolume(part4solid, Materials::get("SUS304"), "part4logical", 0, 0, 0);
85  part4logical->SetVisAttributes(att("iron"));
86  auto pv = new G4PVPlacement(gTrans * G4RotateY3D(M_PI), part4logical, "part4physical", top, false, 0, 0);
87  if (overlap) pv->CheckOverlaps(npoints);
88  }
89 
90  zr_t cont_array_in[] = {{3., 474.9}, {434., 702.27}, {434, 1496 - 20}, {434 - 214.8, 1496 - 20}, {3, 1190.2}};
91  std::vector<zr_t> contour_in(cont_array_in, cont_array_in + sizeof(cont_array_in) / sizeof(zr_t));
92  G4VSolid* innervolume_solid = new BelleLathe("innervolume_solid", 0, 2 * M_PI, contour_in);
93  G4LogicalVolume* innervolume_logical = new G4LogicalVolume(innervolume_solid, Materials::get("G4_AIR"),
94  "innervolume_logical", 0, 0, 0);
95  innervolume_logical->SetVisAttributes(att("air"));
96 
97  // Set up region for production cuts
98  G4Region* aRegion = new G4Region("ECLBackwardEnvelope");
99  innervolume_logical->SetRegion(aRegion);
100  aRegion->AddRootLogicalVolume(innervolume_logical);
101 
102  auto gpvbp = new G4PVPlacement(gTrans * G4RotateY3D(M_PI),
103  innervolume_logical, "ECLBackwardPhysical",
104  top, false, 0, 0);
105  if (overlap) gpvbp->CheckOverlaps(npoints);
106 
107  G4VSolid* innervolumesector_solid = new BelleLathe("innervolumesector_solid", -M_PI / 8, M_PI / 4, contour_in);
108  G4LogicalVolume* innervolumesector_logical = new G4LogicalVolume(innervolumesector_solid, Materials::get("G4_AIR"),
109  "innervolumesector_logical", 0, 0, 0);
110  innervolumesector_logical->SetVisAttributes(att("air"));
111  new G4PVReplica("ECLBackwardSectorPhysical", innervolumesector_logical, innervolume_logical, kPhi, 8, M_PI / 4, 0);
112 
113  // cppcheck-suppress knownConditionTrueFalse
114  if (b_ribs) {
115  double H = 60, W = 20;
116  double X0 = 702.27 + 0.001, X1 = 1496 - 20;
117  G4TwoVector r0o(X1, 0), r1o(X1 * sqrt(1 - pow(W / X1, 2)), W);
118  double beta = asin(W / X0);
119  G4TwoVector r0i(X0 / cos(beta / 2), 0), r1i(X0 * cos(beta / 2) - tan(beta / 2) * (W - X0 * sin(beta / 2)), W);
120  double dxymzp = (r0o - r0i).x(), dxypzp = (r1o - r1i).x();
121  double theta = atan(tand(27.81) / 2);
122  double dxymzm = dxymzp + tand(27.81) * H, dxypzm = dxypzp + tand(27.81) * H;
123 
124  G4TwoVector m0 = (r0i + r0o) * 0.5, m1 = (r1i + r1o) * 0.5, dm = m1 - m0;
125  double alpha = atan(dm.x() / dm.y());
126 
127  G4VSolid* solid6_p1 = new G4Trap("solid6_p1", H / 2, theta, 0, W / 2, dxymzm / 2, dxypzm / 2, alpha, W / 2, dxymzp / 2, dxypzp / 2,
128  alpha);
129  G4LogicalVolume* lsolid6_p1 = new G4LogicalVolume(solid6_p1, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
130  G4VisAttributes* asolid6 = new G4VisAttributes(G4Colour(1., 0.3, 0.2));
131  lsolid6_p1->SetVisAttributes(asolid6);
132  G4Transform3D tsolid6_p1(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, W / 2, 434 - H / 2));
133  auto pv61 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*tsolid6_p1, lsolid6_p1, "psolid6_p1", innervolumesector_logical, false, 0, 0);
134  if (overlap) pv61->CheckOverlaps(npoints);
135  auto pv62 = new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p1, lsolid6_p1, "psolid6_p2", innervolumesector_logical, false, 0, 0);
136  if (overlap) pv62->CheckOverlaps(npoints);
137 
138  H = 40;
139  dxymzm = dxymzp + tand(27.81) * H, dxypzm = dxypzp + tand(27.81) * H;
140  G4VSolid* solid6_p2 = new G4Trap("solid6_p2", H / 2, theta, 0, W / 2, dxypzm / 2, dxymzm / 2, -alpha, W / 2, dxypzp / 2, dxymzp / 2,
141  -alpha);
142  G4LogicalVolume* lsolid6_p2 = new G4LogicalVolume(solid6_p2, Materials::get("SUS304"), "lsolid6", 0, 0, 0);
143  lsolid6_p2->SetVisAttributes(asolid6);
144  G4Transform3D tsolid6_p2(G4Translate3D(X0 * cos(beta / 2) + (dxymzp / 2 + dxypzp / 2) / 2 - tan(theta)*H / 2, -W / 2, 434 - H / 2));
145  auto pv63 = new G4PVPlacement(G4RotateZ3D(0)*tsolid6_p2, lsolid6_p2, "psolid6_p3", innervolumesector_logical, false, 0, 0);
146  if (overlap) pv63->CheckOverlaps(npoints);
147  auto pv64 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid6_p2, lsolid6_p2, "psolid6_p4", innervolumesector_logical, false, 0, 0);
148  if (overlap) pv64->CheckOverlaps(npoints);
149 
150  G4VSolid* solid7_p8 = new G4Box("solid7_p8", 171. / 2, (140. - 40) / 2 / 2, 40. / 2);
151  G4LogicalVolume* lsolid7 = new G4LogicalVolume(solid7_p8, Materials::get("SUS304"), "lsolid7", 0, 0, 0);
152  G4VisAttributes* asolid7 = new G4VisAttributes(G4Colour(1., 0.3, 0.2));
153  lsolid7->SetVisAttributes(asolid7);
154  double dx = sqrt(X1 * X1 - 70 * 70) - 171. / 2;
155  G4Transform3D tsolid7_p1(G4Translate3D(dx, -20 - 25, 434 - 40. / 2));
156  auto pv71 = new G4PVPlacement(tsolid7_p1, lsolid7, "psolid7_p1", innervolumesector_logical, false, 0, 0);
157  if (overlap) pv71->CheckOverlaps(npoints);
158  G4Transform3D tsolid7_p2(G4Translate3D(dx, 20 + 25, 434 - 40. / 2));
159  auto pv72 = new G4PVPlacement(tsolid7_p2, lsolid7, "psolid7_p2", innervolumesector_logical, false, 0, 0);
160  if (overlap) pv72->CheckOverlaps(npoints);
161 
162  double L = X1 - (X0 - tand(27.81) * 40) - 10;
163  G4VSolid* solid13 = new G4Box("solid13", L / 2, 5. / 2, 18. / 2);
164  G4LogicalVolume* lsolid13 = new G4LogicalVolume(solid13, Materials::get("SUS304"), "lsolid13", 0, 0, 0);
165  G4VisAttributes* asolid13 = new G4VisAttributes(G4Colour(1., 0.5, 0.5));
166  lsolid13->SetVisAttributes(asolid13);
167  G4Transform3D tsolid13(G4TranslateZ3D(434 - 60 + 18. / 2)*G4TranslateY3D(-5. / 2 - 0.5 / 2)*G4TranslateX3D(X0 - tand(
168  27.81) * 40 + L / 2 + 5));
169  auto pv131 = new G4PVPlacement(tsolid13, lsolid13, "psolid13_p1", innervolumesector_logical, false, 0, 0);
170  if (overlap) pv131->CheckOverlaps(npoints);
171  auto pv132 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*tsolid13, lsolid13, "psolid13_p2", innervolumesector_logical, false, 0, 0);
172  if (overlap) pv132->CheckOverlaps(npoints);
173  }
174 
175 
176  double zsep = 135;
177  // cppcheck-suppress knownConditionTrueFalse
178  if (b_septum_wall) {
179  double d = 5;
180  Point_t vin[] = {{434. - zsep, 702.27 - tand(27.81)* zsep}, {434. - 60, 702.27 - tand(27.81) * 60}, {434. - 60, 1496 - 20 - d}, {434. - zsep, 1496 - 20 - d}};
181  const int n = sizeof(vin) / sizeof(Point_t);
182  Point_t c = centerofgravity(vin, vin + n);
183  G4ThreeVector contour_swall[n * 2];
184  for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2);
185  for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2);
186 
187  G4VSolid* septumwall_solid = new BelleCrystal("septumwall_solid", n, contour_swall);
188 
189  G4LogicalVolume* septumwall_logical = new G4LogicalVolume(septumwall_solid, Materials::get("A5052"),
190  "septumwall_logical", 0, 0, 0);
191  septumwall_logical->SetVisAttributes(att("alum2"));
192  auto pv = new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0), septumwall_logical,
193  "septumwall_physical", innervolumesector_logical, false, 0, 0);
194  if (overlap) pv->CheckOverlaps(npoints);
195 
196  for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
197  for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
198 
199  G4VSolid* septumwall2_solid = new BelleCrystal("septumwall2_solid", n, contour_swall);
200 
201  G4LogicalVolume* septumwall2_logical = new G4LogicalVolume(septumwall2_solid, Materials::get("A5052"),
202  "septumwall2_logical", 0, 0, 0);
203  septumwall2_logical->SetVisAttributes(att("alum2"));
204  auto pv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
205  0.5 / 2 / 2),
206  septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 0, 0);
207  if (overlap) pv0->CheckOverlaps(npoints);
208  auto pv1 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
209  -0.5 / 2 / 2),
210  septumwall2_logical, "septumwall2_physical", innervolumesector_logical, false, 1, 0);
211  if (overlap) pv1->CheckOverlaps(npoints);
212  }
213 
214  zr_t vcr[] = {{3., 474.9}, {434. - zsep, 702.27 - tand(27.81)* zsep}, {434 - zsep, 1496 - 20}, {434 - 214.8, 1496 - 20}, {3, 1190.2}};
215  std::vector<zr_t> ccr(vcr, vcr + sizeof(vcr) / sizeof(zr_t));
216  G4VSolid* crystalvolume_solid = new BelleLathe("crystalvolume_solid", 0, M_PI / 8, ccr);
217  G4LogicalVolume* crystalvolume_logical = new G4LogicalVolume(crystalvolume_solid, Materials::get("G4_AIR"),
218  "crystalvolume_logical", 0, 0, 0);
219  crystalvolume_logical->SetVisAttributes(att("air"));
220  auto gpv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 8), crystalvolume_logical, "ECLBackwardCrystalSectorPhysical_0",
221  innervolumesector_logical,
222  false, 0, 0);
223  if (overlap) gpv0->CheckOverlaps(npoints);
224  auto gpv1 = new G4PVPlacement(G4RotateZ3D(0), crystalvolume_logical, "ECLBackwardCrystalSectorPhysical_1",
225  innervolumesector_logical, false, 1,
226  0);
227  if (overlap) gpv1->CheckOverlaps(npoints);
228 
229  if (b_septum_wall) {
230  double d = 5, dr = 0.001;
231  Point_t vin[] = {{3., 474.9}, {434. - zsep, 702.27 - tand(27.81)* zsep}, {434 - zsep, 1496 - 20 - d - dr}, {434 - 214.8 - d / tand(52.90), 1496 - 20 - d - dr}, {3, 1190.2 - dr}};
232  const int n = sizeof(vin) / sizeof(Point_t);
233  Point_t c = centerofgravity(vin, vin + n);
234  G4ThreeVector contour_swall[n * 2];
235 
236  for (int i = 0; i < n; i++) contour_swall[i + 0] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, -0.5 / 2 / 2);
237  for (int i = 0; i < n; i++) contour_swall[i + n] = G4ThreeVector(vin[i].x - c.x, vin[i].y - c.y, 0.5 / 2 / 2);
238 
239  G4VSolid* septumwall3_solid = new BelleCrystal("septumwall3_solid", n, contour_swall);
240 
241  G4LogicalVolume* septumwall3_logical = new G4LogicalVolume(septumwall3_solid, Materials::get("A5052"),
242  "septumwall3_logical", 0, 0, 0);
243  septumwall3_logical->SetVisAttributes(att("alum2"));
244  auto pv0 = new G4PVPlacement(G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y, 0.5 / 2 / 2),
245  septumwall3_logical,
246  "septumwall3_physical_0", crystalvolume_logical, false, 0, overlap);
247  if (overlap) pv0->CheckOverlaps(npoints);
248  auto pv1 = new G4PVPlacement(G4RotateZ3D(M_PI / 8)*G4RotateZ3D(-M_PI / 2)*G4RotateY3D(-M_PI / 2)*G4Translate3D(c.x, c.y,
249  -0.5 / 2 / 2),
250  septumwall3_logical, "septumwall3_physical_1", crystalvolume_logical, false, 1, overlap);
251  if (overlap) pv1->CheckOverlaps(npoints);
252  }
253 
254  // cppcheck-suppress knownConditionTrueFalse
255  if (b_crystals) {
256  // vector<shape_t*> cryst = load_shapes("/ecl/data/crystal_shape_backward.dat");
257  vector<shape_t*> cryst = load_shapes(m_sap, ECLParts::backward);
258  vector<G4LogicalVolume*> wrapped_crystals;
259  for (auto it = cryst.begin(); it != cryst.end(); it++) {
260  shape_t* s = *it;
261  wrapped_crystals.push_back(wrapped_crystal(s, "backward", 0.20 - 0.02));
262  }
263  for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
264  const cplacement_t& t = *it;
265  auto s = find_if(cryst.begin(), cryst.end(), [&t](const shape_t* shape) {return shape->nshape == t.nshape;});
266  if (s == cryst.end()) continue;
267 
268  G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(t);
269  int indx = it - bp.begin();
270  auto pv = new G4PVPlacement(twc, wrapped_crystals[s - cryst.begin()], suf("ECLBackwardWrappedCrystal_Physical", indx),
271  crystalvolume_logical,
272  false, (1152 + 6624) / 16 + indx, 0);
273  if (overlap)pv->CheckOverlaps(npoints);
274  }
275  }
276 
277  // cppcheck-suppress knownConditionTrueFalse
278  if (b_preamplifier) {
279  for (vector<cplacement_t>::const_iterator it = bp.begin(); it != bp.end(); ++it) {
280  G4Transform3D twc = G4Translate3D(0, 0, 3) * get_transform(*it);
281  int indx = it - bp.begin();
282  auto pv = new G4PVPlacement(twc * G4TranslateZ3D(300 / 2 + 0.20 + get_pa_box_height() / 2)*G4RotateZ3D(-M_PI / 2), get_preamp(),
283  suf("phys_backward_preamplifier", indx), crystalvolume_logical, false, (1152 + 6624) / 16 + indx, 0);
284  if (overlap)pv->CheckOverlaps(npoints);
285  }
286  }
287 
288  // cppcheck-suppress knownConditionTrueFalse
289  if (b_support_leg) {
290  const G4VisAttributes* batt = att("iron");
291 
292  G4VSolid* s1 = new G4Box("leg_p1", 130. / 2, 185. / 2, (40. - 5) / 2);
293  G4LogicalVolume* l1 = new G4LogicalVolume(s1, Materials::get("SUS304"), "l1", 0, 0, 0);
294  G4Transform3D t1 = G4Translate3D(0, 185. / 2, (40. - 5) / 2);
295  l1->SetVisAttributes(batt);
296 
297  Point_t v3[] = {{ -212. / 2, -135. / 2}, {212. / 2 - 30, -135. / 2}, {212. / 2, -135. / 2 + 30}, {212. / 2, 135. / 2} , { -212. / 2, 135. / 2}};
298  const int n3 = sizeof(v3) / sizeof(Point_t);
299  G4ThreeVector c3[n3 * 2];
300 
301  for (int i = 0; i < n3; i++) c3[i + 0] = G4ThreeVector(v3[i].x, v3[i].y, -60. / 2);
302  for (int i = 0; i < n3; i++) c3[i + n3] = G4ThreeVector(v3[i].x, v3[i].y, 60. / 2);
303 
304  G4VSolid* s3 = new BelleCrystal("leg_p3", n3, c3);
305  G4LogicalVolume* l3 = new G4LogicalVolume(s3, Materials::get("SUS304"), "l3", 0, 0, 0);
306  G4Transform3D t3 = G4Translate3D(0, 135. / 2 + 35, 40. - 5. + 212. / 2) * G4RotateY3D(-M_PI / 2);
307  l3->SetVisAttributes(batt);
308 
309  G4VSolid* s6 = new G4Box("leg_p6", 140. / 2, 189. / 2, 160. / 2);
310  G4LogicalVolume* l6 = new G4LogicalVolume(s6, Materials::get("G4_AIR"), "l6", 0, 0, 0);
311  G4Transform3D t6 = G4Translate3D(0, 170. + 189. / 2, 57. + 35. + 160. / 2);
312  l6->SetVisAttributes(att("air"));
313 
314  G4VSolid* s6a = new G4Box("leg_p6a", 140. / 2, (189. - 45.) / 2, 160. / 2);
315  G4LogicalVolume* l6a = new G4LogicalVolume(s6a, Materials::get("SUS304"), "l6a", 0, 0, 0);
316  l6a->SetVisAttributes(batt);
317  new G4PVPlacement(G4TranslateY3D(-45. / 2), l6a, "l6a_physical", l6, false, 0, overlap);
318 
319  G4VSolid* s6b = new G4Box("leg_p6b", 60. / 2, 45. / 2, 160. / 2);
320  G4LogicalVolume* l6b = new G4LogicalVolume(s6b, Materials::get("SUS304"), "l6b", 0, 0, 0);
321  l6b->SetVisAttributes(batt);
322  double dy = 189. / 2 - 45 + 45. / 2;
323  new G4PVPlacement(G4TranslateY3D(dy), l6b, "l6b_physical", l6, false, 0, overlap);
324 
325  G4VSolid* s6c = new G4Box("leg_p6c", 40. / 2, 45. / 2, 22.5 / 2);
326  G4LogicalVolume* l6c = new G4LogicalVolume(s6c, Materials::get("SUS304"), "l6c", 0, 0, 0);
327  l6c->SetVisAttributes(batt);
328  new G4PVPlacement(G4Translate3D(30 + 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 0, overlap);
329  new G4PVPlacement(G4Translate3D(30 + 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 1, overlap);
330  new G4PVPlacement(G4Translate3D(-30 - 20, dy, 20 + 22.5 / 2), l6c, "l6c_physical", l6, false, 2, overlap);
331  new G4PVPlacement(G4Translate3D(-30 - 20, dy, -20 - 22.5 / 2), l6c, "l6c_physical", l6, false, 3, overlap);
332 
333  G4AssemblyVolume* support_leg = new G4AssemblyVolume();
334 
335  support_leg->AddPlacedVolume(l1, t1);
336  // support_leg->AddPlacedVolume(l2,t2);
337  support_leg->AddPlacedVolume(l3, t3);
338  // support_leg->AddPlacedVolume(l4,t4);
339  // support_leg->AddPlacedVolume(l5,t5);
340  support_leg->AddPlacedVolume(l6, t6);
341 
342  G4VSolid* s_all = new G4Box("leg_all", 140. / 2, 359. / 2, (257. - 5.) / 2);
343  G4LogicalVolume* l_all = new G4LogicalVolume(s_all, Materials::get("G4_AIR"), "l_all", 0, 0, 0);
344  l_all->SetVisAttributes(att("air"));
345  G4Transform3D tp = G4Translate3D(0, -359. / 2, -(257. - 5.) / 2);
346  support_leg->MakeImprint(l_all, tp, 0, overlap);
347 
348 
349  for (int i = 0; i < 8; i++)
350  new G4PVPlacement(gTrans * G4RotateX3D(M_PI)*G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0,
351  1496 - 185 + 359. / 2,
352  434 + 5 + (257. - 5.) / 2), l_all, "lall_physical", top, false, i, overlap);
353 
354 
355  G4VSolid* s1a = new G4Box("leg_p1a", 130. / 2, 178. / 2, 5. / 2);
356  G4LogicalVolume* l1a = new G4LogicalVolume(s1a, Materials::get("SUS304"), "l1a", 0, 0, 0);
357  l1a->SetVisAttributes(batt);
358  for (int i = 0; i < 8; i++)
359  new G4PVPlacement(gTrans * G4RotateX3D(M_PI)*G4RotateZ3D(-M_PI / 2 + M_PI / 8 + i * M_PI / 4)*G4Translate3D(0,
360  1496 - 185 + 178. / 2,
361  434 + 5 - 5. / 2), l1a, "l1a_physical", top, false, i, overlap);
362 
363  }
364 
365 }
Belle2::ECL::shape_t
Definition: shapes.h:39
Belle2::geometry
Common code concerning the geometry representation of the detector.
Definition: CreatorBase.h:28
Belle2::ECL::Point_t
Definition: BelleCrystal.h:40
Belle2::ECL::GeoECLCreator::backward
void backward(G4LogicalVolume &)
Place elements inside the backward endcap.
Definition: backward.cc:24
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::ECL::BelleCrystal
a Belle crystal in Geant4
Definition: BelleCrystal.h:45
Belle2::ECL::zr_t
Definition: BelleLathe.h:35
Belle2::ECL::BelleLathe
Definition: BelleLathe.h:60
Belle2::ECL::cplacement_t
Definition: shapes.h:52