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