Belle II Software  release-05-01-25
Ph1sustrCreator.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter, Igal Jaegle *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <beast/ph1sustr/geometry/Ph1sustrCreator.h>
12 #include <beast/ph1sustr/simulation/SensitiveDetector.h>
13 
14 #include <geometry/Materials.h>
15 #include <geometry/CreatorFactory.h>
16 #include <framework/gearbox/GearDir.h>
17 #include <framework/logging/Logger.h>
18 
19 #include <cmath>
20 #include <boost/format.hpp>
21 #include <boost/foreach.hpp>
22 #include <boost/algorithm/string.hpp>
23 
24 #include <G4LogicalVolume.hh>
25 #include <G4PVPlacement.hh>
26 
27 //Shapes
28 #include <G4Box.hh>
29 #include "G4UnionSolid.hh"
30 #include "G4SubtractionSolid.hh"
31 
32 //Visualization
33 #include "G4Colour.hh"
34 #include <G4VisAttributes.hh>
35 
36 using namespace std;
37 using namespace boost;
38 
39 namespace Belle2 {
46  namespace ph1sustr {
47 
48  // Register the creator
50  geometry::CreatorFactory<Ph1sustrCreator> Ph1sustrFactory("PH1SUSTRCreator");
51 
52  Ph1sustrCreator::Ph1sustrCreator(): m_sensitive(0)
53  {
55  }
56 
57  Ph1sustrCreator::~Ph1sustrCreator()
58  {
59  if (m_sensitive) delete m_sensitive;
60  }
61 
62  void Ph1sustrCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
63  {
64  //lets get the stepsize parameter with a default value of 5 µm
65  //double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
66 
67 
68  //no get the array. Notice that the default framework unit is cm, so the
69  //values will be automatically converted
70  vector<double> bar = content.getArray("bar");
71  B2INFO("Contents of bar: ");
72  BOOST_FOREACH(double value, bar) {
73  B2INFO("value: " << value);
74  }
75  /*double x_tpcbeamR = 0;
76  double y_tpcbeamR = 0;
77  double z_tpcbeamR = 0;
78  double x_tpcbeamL = 0;
79  double y_tpcbeamL = 0;
80  double z_tpcbeamL = 0;
81  double x_tpcbeamT = 0;
82  double y_tpcbeamT = 0;
83  double z_tpcbeamT = 0;
84  double x_tpcbeamB = 0;
85  double y_tpcbeamB = 0;
86  double z_tpcbeamB = 0;*/
87  //Lets loop over all the Active nodes
88  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
89 
90  //plate positions
91  double x_tpcbeamR = activeParams.getLength("x_tpcbeamR") * CLHEP::cm;
92  double y_tpcbeamR = activeParams.getLength("y_tpcbeamR") * CLHEP::cm;
93  double z_tpcbeamR = activeParams.getLength("z_tpcbeamR") * CLHEP::cm;
94  double x_tpcbeamL = activeParams.getLength("x_tpcbeamL") * CLHEP::cm;
95  double y_tpcbeamL = activeParams.getLength("y_tpcbeamL") * CLHEP::cm;
96  double z_tpcbeamL = activeParams.getLength("z_tpcbeamL") * CLHEP::cm;
97  double x_tpcbeamT = activeParams.getLength("x_tpcbeamT") * CLHEP::cm;
98  double y_tpcbeamT = activeParams.getLength("y_tpcbeamT") * CLHEP::cm;
99  double z_tpcbeamT = activeParams.getLength("z_tpcbeamT") * CLHEP::cm;
100  double x_tpcbeamB = activeParams.getLength("x_tpcbeamB") * CLHEP::cm;
101  double y_tpcbeamB = activeParams.getLength("y_tpcbeamB") * CLHEP::cm;
102  double z_tpcbeamB = activeParams.getLength("z_tpcbeamB") * CLHEP::cm;
103 
104  //TPC vertical: 4x @ 1614/ea
105  //TPC horizontal + BGO base: 8x @ 1583/ea
106  //TPC railroad: 8x @ 2200/ea
107  //BGO vertical: 4x @ 928/ea
108  //BGO horizontal: 4x @ 318/ea
109  //G4double dz_20V2100bgov = 843.72 / 2.*CLHEP::mm;
110  //G4double dz_20V2100bgoh = 280.00 / 2.*CLHEP::mm;
111 
112  //Beam supporting the TPC-Tube-plate
113  //define tpc beam and plate dimensions
114  double betpcbeam = 190.8 / 2. * CLHEP::mm;
115  G4double dx_tpcbeam = 2.54 * 1.63 / 2.*CLHEP::cm;
116  G4double dy_tpcbeam = 2.54 * 1.63 / 2.*CLHEP::cm;
117  G4double dz_tpcbeam = 2200. / 2.*CLHEP::mm;
118  G4double dw_tpcbeam = 2.54 * 0.25 / 2.*CLHEP::cm;
119  //G4double dx_plate = 2.54 * 0.35 / 2.*CLHEP::cm;
120 
121  //G4double dx_plate_short = 0.4765 * CLHEP::cm;
122  G4double dx_plate_short = 0.5 / 2. * CLHEP::cm;
123  G4double dy_plate_short = 27.47788 / 2.*CLHEP::cm;
124  G4double dz_plate_short = 40. / 2.*CLHEP::cm;
125 
126  G4double dx_plate = 0.5 / 2. * CLHEP::cm;
127  G4double dy_plate = 32. / 2.*CLHEP::cm;
128  G4double dz_plate = 50. / 2.*CLHEP::cm;
129 
130  //Right from e^-~--~beam
131  //create plate volume
132  G4VSolid* s_plate = new G4Box("s_plate", dx_plate, dy_plate, dz_plate);
133  G4VSolid* s_plate_short = new G4Box("s_plate_short", dx_plate_short, dy_plate_short, dz_plate_short);
134 
135  //place plate volume
136  G4LogicalVolume* l_plate = new G4LogicalVolume(s_plate, geometry::Materials::get("Al") , "l_plate", 0, 0);
137  G4LogicalVolume* l_plate_short = new G4LogicalVolume(s_plate_short, geometry::Materials::get("Al") , "l_plate_short", 0, 0);
138  G4VisAttributes* white = new G4VisAttributes(G4Colour(1, 1, 1));
139  white->SetForceAuxEdgeVisible(true);
140  l_plate->SetVisAttributes(white);
141  l_plate_short->SetVisAttributes(white);
142 
143  G4VSolid* s_tpcbeam_a = new G4Box("s_tpcbeam_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeam);
144  G4VSolid* s_tpcbeam_b = new G4Box("s_tpcbeam_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeam);
145  G4VSolid* s_tpcbeampos = new G4SubtractionSolid("s_tpcbeampos", s_tpcbeam_a, s_tpcbeam_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
146  G4VSolid* s_tpcbeamneg = new G4SubtractionSolid("s_tpcbeamneg", s_tpcbeam_a, s_tpcbeam_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
147  G4VSolid* s_tpcbeam = new G4UnionSolid("s_tpcbeam", s_tpcbeampos, s_tpcbeamneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
148 
149  //create tpc beam volumes
150  G4LogicalVolume* l_tpcbeam = new G4LogicalVolume(s_tpcbeam, geometry::Materials::get("FG_Epoxy") , "l_tpcbeam", 0, 0);
151 
152  //place plate volume
153  G4ThreeVector PH1SUSTRpos = G4ThreeVector(
154  x_tpcbeamR,
155  y_tpcbeamR,
156  z_tpcbeamR
157  );
158  new G4PVPlacement(0, PH1SUSTRpos, l_plate, "p_plateR", &topVolume, false, 1);
159 
160  //offsets hori
161  //G4double offset_h = /*110.*CLHEP::cm*/dz_tpcbeam + /*80.*CLHEP::cm*/ dz_plate - 2. * dz_tpcbeam;
162  G4double offset_h = /*110.*CLHEP::cm*/dz_tpcbeam + 80.*CLHEP::cm - 2. * dz_tpcbeam;
163 
164  //place 1st tpc beam volume
165  /*PH1SUSTRpos = G4ThreeVector(
166  x_tpcbeamR + dx_plate + dx_tpcbeam,
167  y_tpcbeamR + betpcbeam + 2. * dy_tpcbeam,
168  -offset_h / 2.
169  );
170  new G4PVPlacement(0, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
171  G4Transform3D TransForm = G4Translate3D(x_tpcbeamR + dx_plate + dx_tpcbeam,
172  y_tpcbeamR + betpcbeam + 2. * dy_tpcbeam,
173  -offset_h / 2.) * G4RotateZ3D(90.*CLHEP::deg);
174  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamR1", &topVolume, false, 1);
175 
176  //place 2nd tpc beam volume
177  /*PH1SUSTRpos = G4ThreeVector(
178  x_tpcbeamR + dx_plate + dx_tpcbeam,
179  y_tpcbeamR - betpcbeam - 2. * dy_tpcbeam,
180  -offset_h / 2.
181  );
182  new G4PVPlacement(0, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
183  TransForm = G4Translate3D(x_tpcbeamR + dx_plate + dx_tpcbeam,
184  y_tpcbeamR - betpcbeam - 2. * dy_tpcbeam,
185  -offset_h / 2.) * G4RotateZ3D(90.*CLHEP::deg);
186  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamR2", &topVolume, false, 1);
187 
188  //Left from e^-~--~beam
189  //place plate volume
190  PH1SUSTRpos = G4ThreeVector(
191  x_tpcbeamL,
192  y_tpcbeamL,
193  z_tpcbeamL
194  );
195  new G4PVPlacement(0, PH1SUSTRpos, l_plate, "p_plateL", &topVolume, false, 1);
196 
197  //place 1st tpc beam volume
198  /*PH1SUSTRpos = G4ThreeVector(
199  x_tpcbeamL - dx_plate - dx_tpcbeam,
200  y_tpcbeamL + betpcbeam + 2. * dy_tpcbeam,
201  -offset_h / 2.
202  );
203  new G4PVPlacement(0, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
204  TransForm = G4Translate3D(x_tpcbeamL - dx_plate - 3. * dx_tpcbeam,
205  y_tpcbeamL + betpcbeam + 2. * dy_tpcbeam,
206  -offset_h / 2.) * G4RotateZ3D(90.*CLHEP::deg);
207  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamL1", &topVolume, false, 1);
208 
209  //place 2nd tpc beam volume
210  /*PH1SUSTRpos = G4ThreeVector(
211  x_tpcbeamL - dx_plate - dx_tpcbeam,
212  y_tpcbeamL - betpcbeam - 2. * dy_tpcbeam,
213  -offset_h / 2.
214  );
215  new G4PVPlacement(0, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
216  TransForm = G4Translate3D(x_tpcbeamL - dx_plate - 3. * dx_tpcbeam,
217  y_tpcbeamL - betpcbeam - 2. * dy_tpcbeam,
218  -offset_h / 2.) * G4RotateZ3D(90.*CLHEP::deg);
219  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamL2", &topVolume, false, 1);
220 
221  //Bottom
222  G4RotationMatrix* rotXx = new G4RotationMatrix();
223  double Angle = 90. * CLHEP::deg;
224  rotXx->rotateZ(Angle);
225  //place bottom plate
226  PH1SUSTRpos = G4ThreeVector(
227  x_tpcbeamB,
228  y_tpcbeamB,
229  z_tpcbeamB
230  );
231  new G4PVPlacement(rotXx, PH1SUSTRpos, l_plate_short, "p_plateB", &topVolume, false, 1);
232 
233  //place 1st tpc beam volume
234  /*PH1SUSTRpos = G4ThreeVector(
235  x_tpcbeamB + betpcbeam + 2. * dx_tpcbeam,
236  y_tpcbeamB - dx_plate - dy_tpcbeam,
237  -offset_h / 2.
238  );
239  new G4PVPlacement(rotXx, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
240  TransForm = G4Translate3D(x_tpcbeamB + betpcbeam + 1.5 * dx_tpcbeam,
241  y_tpcbeamB - dx_plate - 1.5 * dy_tpcbeam,
242  -offset_h / 2.) /* G4RotateZ3D(90.*CLHEP::deg)*/;
243  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamB1", &topVolume, false, 1);
244 
245  //place 2nd tpc beam volume
246  /*PH1SUSTRpos = G4ThreeVector(
247  x_tpcbeamB - betpcbeam - 2. * dx_tpcbeam,
248  y_tpcbeamB - dx_plate - dy_tpcbeam,
249  -offset_h / 2.
250  );
251  new G4PVPlacement(rotXx, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
252  TransForm = G4Translate3D(x_tpcbeamB - betpcbeam - 1.5 * dx_tpcbeam,
253  y_tpcbeamB - dx_plate - 1.5 * dy_tpcbeam,
254  -offset_h / 2.) /* G4RotateZ3D(90.*CLHEP::deg)*/;
255  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamB2", &topVolume, false, 1);
256 
257  //Top
258  //rotXx = new G4RotationMatrix();
259  //rotXx->rotateZ( 90. );
260  //place top plate
261  PH1SUSTRpos = G4ThreeVector(
262  x_tpcbeamT,
263  y_tpcbeamT,
264  z_tpcbeamT
265  );
266  new G4PVPlacement(rotXx, PH1SUSTRpos, l_plate, "p_plateT", &topVolume, false, 1);
267 
268  //place 1st tpc beam volume
269  /*PH1SUSTRpos = G4ThreeVector(
270  x_tpcbeamT + betpcbeam + 2. * dx_tpcbeam,
271  y_tpcbeamT + dx_plate + dy_tpcbeam,
272  -offset_h / 2.
273  );
274  new G4PVPlacement(rotXx, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
275  TransForm = G4Translate3D(x_tpcbeamT + betpcbeam + 1.5 * dx_tpcbeam,
276  y_tpcbeamT + dx_plate + 3. * dy_tpcbeam,
277  -offset_h / 2.) /* G4RotateZ3D(90.*CLHEP::deg)*/;
278  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamT1", &topVolume, false, 1);
279 
280  //place 2nd tpc beam volume
281  /*PH1SUSTRpos = G4ThreeVector(
282  x_tpcbeamT - betpcbeam - 2. * dx_tpcbeam,
283  y_tpcbeamT + dx_plate + dy_tpcbeam,
284  -offset_h / 2.
285  );
286  new G4PVPlacement(rotXx, PH1SUSTRpos, l_tpcbeam, "p_tpcbeam", &topVolume, false, 1);*/
287  TransForm = G4Translate3D(x_tpcbeamT - betpcbeam - 1.5 * dx_tpcbeam,
288  y_tpcbeamT + dx_plate + 3 * dy_tpcbeam,
289  -offset_h / 2.) /* G4RotateZ3D(90.*CLHEP::deg)*/;
290  new G4PVPlacement(TransForm, l_tpcbeam, "p_tpcbeamT2", &topVolume, false, 1);
291 
292  //vertical beams
293  //G4double dz_tpcbeamv = 1537.37 / 2.*CLHEP::mm;
294  G4double dz_tpcbeamv = 1583. / 2.*CLHEP::mm;
295  G4VSolid* s_tpcbeamv_a = new G4Box("s_tpcbeamv_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeamv);
296  G4VSolid* s_tpcbeamv_b = new G4Box("s_tpcbeamv_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeamv);
297  G4VSolid* s_tpcbeamvpos = new G4SubtractionSolid("s_tpcbeamvpos", s_tpcbeamv_a, s_tpcbeamv_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
298  G4VSolid* s_tpcbeamvneg = new G4SubtractionSolid("s_tpcbeamvneg", s_tpcbeamv_a, s_tpcbeamv_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
299  G4VSolid* s_tpcbeamv = new G4UnionSolid("s_tpcbeamv", s_tpcbeamvpos, s_tpcbeamvneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
300  G4LogicalVolume* l_tpcbeamv = new G4LogicalVolume(s_tpcbeamv, geometry::Materials::get("FG_Epoxy") , "l_tpcbeamv", 0, 0);
301 
302  //offset verti
303  //G4double offset_v = fabs(76.*CLHEP::cm - 2. * dz_tpcbeamv) / 2.;
304  G4double offset_v = fabs(110.*CLHEP::cm - 2. * dz_tpcbeamv) / 2.;
305 
306  //place 1st vertical TPC beam
307  G4RotationMatrix* rotX = new G4RotationMatrix();
308  rotX->rotateX(90.*CLHEP::deg);
309  PH1SUSTRpos = G4ThreeVector(
310  x_tpcbeamL - dx_plate - 3 * dx_tpcbeam - 2.*dx_tpcbeam,
311  offset_v / 2.,
312  -800.*CLHEP::mm
313  );
314  new G4PVPlacement(rotX, PH1SUSTRpos, l_tpcbeamv, "p_tpcbeamv1", &topVolume, false, 0);
315  //place 2nd vertical TPC beam
316  PH1SUSTRpos = G4ThreeVector(
317  x_tpcbeamR + dx_plate + 3.*dx_tpcbeam + 2.*dx_tpcbeam,
318  offset_v / 2.,
319  -800.*CLHEP::mm
320  );
321  new G4PVPlacement(rotX, PH1SUSTRpos, l_tpcbeamv, "p_tpcbeamv2", &topVolume, false, 0);
322  //place 3rd vertical TPC beam
323  PH1SUSTRpos = G4ThreeVector(
324  x_tpcbeamL - dx_plate - 3.*dx_tpcbeam - 2.*dx_tpcbeam,
325  offset_v / 2.,
326  1100.*CLHEP::mm
327  );
328  new G4PVPlacement(rotX, PH1SUSTRpos, l_tpcbeamv, "p_tpcbeamv3", &topVolume, false, 0);
329  //place 4th vertical TPC beam
330  PH1SUSTRpos = G4ThreeVector(
331  x_tpcbeamR + dx_plate + 3.*dx_tpcbeam + 2.*dx_tpcbeam,
332  offset_v / 2.,
333  1100.*CLHEP::mm
334  );
335  new G4PVPlacement(rotX, PH1SUSTRpos, l_tpcbeamv, "p_tpcbeamv4", &topVolume, false, 0);
336 
337  //horizontal beams
338  //G4double dz_tpcbeamh = 1792. / 2.*CLHEP::mm;
339  G4double dz_tpcbeamh = 1614. / 2.*CLHEP::mm;
340  G4VSolid* s_tpcbeamh_a = new G4Box("s_tpcbeamh_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeamh);
341  G4VSolid* s_tpcbeamh_b = new G4Box("s_tpcbeamh_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeamh);
342  G4VSolid* s_tpcbeamhpos = new G4SubtractionSolid("s_tpcbeamhpos", s_tpcbeamh_a, s_tpcbeamh_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
343  G4VSolid* s_tpcbeamhneg = new G4SubtractionSolid("s_tpcbeanhneg", s_tpcbeamh_a, s_tpcbeamh_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
344  G4VSolid* s_tpcbeamh = new G4UnionSolid("s_tpcbeamh", s_tpcbeamhpos, s_tpcbeamhneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
345 
346  G4LogicalVolume* l_tpcbeamh = new G4LogicalVolume(s_tpcbeamh, geometry::Materials::get("FG_Epoxy") , "l_tpcbeamh", 0, 0);
347 
348  //place 1st horizontal TPC beam
349  G4RotationMatrix* rotY = new G4RotationMatrix();
350  rotY->rotateY(90.*CLHEP::deg);
351  PH1SUSTRpos = G4ThreeVector(
352  0 * CLHEP::mm,
353  y_tpcbeamB - dx_plate - 4. * dy_tpcbeam - 2. * dy_tpcbeam,
354  1100.*CLHEP::mm - 2. * dy_tpcbeam
355  );
356  new G4PVPlacement(rotY, PH1SUSTRpos, l_tpcbeamh, "p_tpcbeamh1", &topVolume, false, 0);
357 
358  //place 2nd horizontal TPC beam
359  PH1SUSTRpos = G4ThreeVector(
360  0 * CLHEP::mm,
361  y_tpcbeamT + dx_plate + 5. * dy_tpcbeam + 2. * dy_tpcbeam,
362  1100.*CLHEP::mm - 2. * dy_tpcbeam
363  );
364  new G4PVPlacement(rotY, PH1SUSTRpos, l_tpcbeamh, "p_tpcbeamh2", &topVolume, false, 0);
365 
366  //place 3rd horizontal TPC beam
367  PH1SUSTRpos = G4ThreeVector(
368  0 * CLHEP::mm,
369  y_tpcbeamB - dx_plate - 4. * dy_tpcbeam - 2. * dy_tpcbeam,
370  -800.*CLHEP::mm + 4. * dy_tpcbeam
371  );
372  new G4PVPlacement(rotY, PH1SUSTRpos, l_tpcbeamh, "p_tpcbeamh3", &topVolume, false, 0);
373 
374  //place 4th horizontal TPC beam
375  PH1SUSTRpos = G4ThreeVector(
376  0 * CLHEP::mm,
377  y_tpcbeamT + dx_plate + 5. * dy_tpcbeam + 2. * dy_tpcbeam,
378  -800.*CLHEP::mm + 4. * dy_tpcbeam
379  );
380  new G4PVPlacement(rotY, PH1SUSTRpos, l_tpcbeamh, "p_tpcbeamh4", &topVolume, false, 0);
381 
382  G4VisAttributes* brown = new G4VisAttributes(G4Colour(.5, .5, 0));
383  brown->SetForceAuxEdgeVisible(true);
384  l_tpcbeam->SetVisAttributes(brown);
385  l_tpcbeamv->SetVisAttributes(brown);
386  l_tpcbeamh->SetVisAttributes(brown);
387 
388  //CsI box beams
389  G4double dz_csibeamh = activeParams.getLength("lcsiBeamh") * CLHEP::cm / 2.;
390  G4VSolid* s_csibeamh_a = new G4Box("s_csibeamh_a", dx_tpcbeam, dy_tpcbeam, dz_csibeamh);
391  G4VSolid* s_csibeamh_b = new G4Box("s_csibeamh_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_csibeamh);
392  G4VSolid* s_csibeamhpos = new G4SubtractionSolid("s_csibeamhpos", s_csibeamh_a, s_csibeamh_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
393  G4VSolid* s_csibeamhneg = new G4SubtractionSolid("s_csibeamhneg", s_csibeamh_a, s_csibeamh_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
394  //G4VSolid* s_csibeamhpos = new G4SubtractionSolid("s_csibeamhpos", s_csibeamh_a, s_csibeamh_b, 0, G4ThreeVector(dw_tpcbeam, 0, 0));
395  //G4VSolid* s_csibeamhneg = new G4SubtractionSolid("s_csibeamhneg", s_csibeamh_a, s_csibeamh_b, 0, G4ThreeVector(-dw_tpcbeam, 0, 0));
396  G4VSolid* s_csibeamh = new G4UnionSolid("s_csibeamh", s_csibeamhpos, s_csibeamhneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
397  G4LogicalVolume* l_csibeamh = new G4LogicalVolume(s_csibeamh, geometry::Materials::get("FG_Epoxy") , "l_csibeamh", 0, 0);
398  int xdimbase = 0;
399  int ydimbase = 0;
400  int zdimbase = 0;
401  int xdimcsiBeamh = 0;
402  int ydimcsiBeamh = 0;
403  int zdimcsiBeamh = 0;
404  int xdimcsiBeamvF = 0;
405  int ydimcsiBeamvF = 0;
406  int zdimcsiBeamvF = 0;
407  int xdimcsiBeamvB = 0;
408  int ydimcsiBeamvB = 0;
409  int zdimcsiBeamvB = 0;
410  int xdimbgobeamv = 0;
411  int ydimbgobeamv = 0;
412  int zdimbgobeamv = 0;
413  int xdimbgobeamh = 0;
414  int ydimbgobeamh = 0;
415  int zdimbgobeamh = 0;
416  int xdimbgobeamt = 0;
417  int ydimbgobeamt = 0;
418  int zdimbgobeamt = 0;
419  int xdimbgobeamb = 0;
420  int ydimbgobeamb = 0;
421  int zdimbgobeamb = 0;
422  int xdimtpcbeamb = 0;
423  int ydimtpcbeamb = 0;
424  int zdimtpcbeamb = 0;
425  /*int xdimtpcbeamhh = 0;
426  int ydimtpcbeamhh = 0;
427  int zdimtpcbeamhh = 0;
428  int xdimtpcbeamvv = 0;
429  int ydimtpcbeamvv = 0;
430  int zdimtpcbeamvv = 0;*/
431  double xbase[100];
432  double ybase[100];
433  double zbase[100];
434  double xcsibeamh[100];
435  double ycsibeamh[100];
436  double zcsibeamh[100];
437  double xcsibeamvF[100];
438  double ycsibeamvF[100];
439  double zcsibeamvF[100];
440  double xcsibeamvB[100];
441  double ycsibeamvB[100];
442  double zcsibeamvB[100];
443  double xbgobeamv[100];
444  double ybgobeamv[100];
445  double zbgobeamv[100];
446  double xbgobeamh[100];
447  double ybgobeamh[100];
448  double zbgobeamh[100];
449  double xbgobeamt[100];
450  double ybgobeamt[100];
451  double zbgobeamt[100];
452  double xbgobeamb[100];
453  double ybgobeamb[100];
454  double zbgobeamb[100];
455  double xtpcbeamb[100];
456  double ytpcbeamb[100];
457  double ztpcbeamb[100];
458  /*double xtpcbeamhh[100];
459  double ytpcbeamhh[100];
460  double ztpcbeamhh[100];
461  double xtpcbeamvv[100];
462  double ytpcbeamvv[100];
463  double ztpcbeamvv[100];*/
464 
465  double x_offset = activeParams.getLength("x_offset") * CLHEP::cm;
466  double y_offset = activeParams.getLength("y_offset") * CLHEP::cm;
467 
468  for (double xcsiBeamh : activeParams.getArray("xcsiBeamh", {0})) {
469  xcsiBeamh *= CLHEP::cm;
470  xcsibeamh[xdimcsiBeamh] = xcsiBeamh - x_offset;
471  xdimcsiBeamh++;
472  }
473  for (double ycsiBeamh : activeParams.getArray("ycsiBeamh", {0})) {
474  ycsiBeamh *= CLHEP::cm;
475  ycsibeamh[ydimcsiBeamh] = ycsiBeamh - y_offset;
476  ydimcsiBeamh++;
477  }
478  for (double zcsiBeamh : activeParams.getArray("zcsiBeamh", {0})) {
479  zcsiBeamh *= CLHEP::cm;
480  zcsibeamh[zdimcsiBeamh] = zcsiBeamh;
481  zdimcsiBeamh++;
482  }
483  for (int i = 0; i < xdimcsiBeamh; i++) {
484  G4Transform3D transform = G4Translate3D(xcsibeamh[i], ycsibeamh[i] - dy_tpcbeam,
485  zcsibeamh[i] + 2. * dy_tpcbeam) * G4RotateY3D(90.*CLHEP::deg) * G4RotateZ3D(90.*CLHEP::deg);
486  new G4PVPlacement(transform, l_csibeamh, TString::Format("p_csibeamh1_%d", i).Data(), &topVolume, false, 1);
487  transform = G4Translate3D(xcsibeamh[i], ycsibeamh[i] + dy_tpcbeam,
488  zcsibeamh[i] + 2. * dy_tpcbeam) * G4RotateY3D(90.*CLHEP::deg) * G4RotateZ3D(90.*CLHEP::deg);
489  new G4PVPlacement(transform, l_csibeamh, TString::Format("p_csibeamh2_%d", i).Data(), &topVolume, false, 1);
490  //PH1SUSTRpos = G4ThreeVector(xcsibeamh[i],ycsibeamh[i],zcsibeamh[i] - dy_tpcbeam);
491  //new G4PVPlacement(rotY, PH1SUSTRpos, l_csibeamh, "p_csibeamh", &topVolume, false, 0);
492  //PH1SUSTRpos = G4ThreeVector(xcsibeamh[i],ycsibeamh[i],zcsibeamh[i] + dy_tpcbeam);
493  //new G4PVPlacement(rotY+rotXx, PH1SUSTRpos, l_csibeamh, "p_csibeamh", &topVolume, false, 0);
494  }
495 
496  G4double dz_csibeamvF = activeParams.getLength("lcsiBeamvF") * CLHEP::cm / 2.;
497  G4VSolid* s_csibeamvF_a = new G4Box("s_csibeamvF_a", dx_tpcbeam, dy_tpcbeam, dz_csibeamvF);
498  G4VSolid* s_csibeamvF_b = new G4Box("s_csibeamvF_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_csibeamvF);
499  G4VSolid* s_csibeamvFpos = new G4SubtractionSolid("s_csibeamvFpos", s_csibeamvF_a, s_csibeamvF_b, 0, G4ThreeVector(0, dw_tpcbeam,
500  0));
501  G4VSolid* s_csibeamvFneg = new G4SubtractionSolid("s_csibeamvFneg", s_csibeamvF_a, s_csibeamvF_b, 0, G4ThreeVector(0, -dw_tpcbeam,
502  0));
503  G4VSolid* s_csibeamvF = new G4UnionSolid("s_csibeamvF", s_csibeamvFpos, s_csibeamvFneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
504  G4LogicalVolume* l_csibeamvF = new G4LogicalVolume(s_csibeamvF, geometry::Materials::get("FG_Epoxy") , "l_csibeamvF", 0, 0);
505 
506  for (double xcsiBeamvF : activeParams.getArray("xcsiBeamvF", {0})) {
507  xcsiBeamvF *= CLHEP::cm;
508  xcsibeamvF[xdimcsiBeamvF] = xcsiBeamvF - x_offset;
509  xdimcsiBeamvF++;
510  }
511  for (double ycsiBeamvF : activeParams.getArray("ycsiBeamvF", {0})) {
512  ycsiBeamvF *= CLHEP::cm;
513  ycsibeamvF[ydimcsiBeamvF] = ycsiBeamvF - y_offset ;
514  ydimcsiBeamvF++;
515  }
516  for (double zcsiBeamvF : activeParams.getArray("zcsiBeamvF", {0})) {
517  zcsiBeamvF *= CLHEP::cm;
518  zcsibeamvF[zdimcsiBeamvF] = zcsiBeamvF ;
519  zdimcsiBeamvF++;
520  }
521  for (int i = 0; i < xdimcsiBeamvF; i++) {
522  //PH1SUSTRpos = G4ThreeVector(xcsibeamvF[i],0,zcsibeamh[0] - 2. * dx_tpcbeam);
523  PH1SUSTRpos = G4ThreeVector(xcsibeamvF[i], ycsibeamvF[i], zcsibeamvF[i]);
524  new G4PVPlacement(rotX, PH1SUSTRpos, l_csibeamvF, TString::Format("p_csibeamvF_%d", i).Data(), &topVolume, false, 0);
525  }
526 
527  G4double dz_csibeamvB = activeParams.getLength("lcsiBeamvB") * CLHEP::cm / 2.;
528  G4VSolid* s_csibeamvB_a = new G4Box("s_csibeamvB_a", dx_tpcbeam, dy_tpcbeam, dz_csibeamvB);
529  G4VSolid* s_csibeamvB_b = new G4Box("s_csibeamvB_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_csibeamvB);
530  G4VSolid* s_csibeamvBpos = new G4SubtractionSolid("s_csibeamvBpos", s_csibeamvB_a, s_csibeamvB_b, 0, G4ThreeVector(0, dw_tpcbeam,
531  0));
532  G4VSolid* s_csibeamvBneg = new G4SubtractionSolid("s_csibeamvBneg", s_csibeamvB_a, s_csibeamvB_b, 0, G4ThreeVector(0, -dw_tpcbeam,
533  0));
534  G4VSolid* s_csibeamvB = new G4UnionSolid("s_csibeamvB", s_csibeamvBpos, s_csibeamvBneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
535  G4LogicalVolume* l_csibeamvB = new G4LogicalVolume(s_csibeamvB, geometry::Materials::get("FG_Epoxy") , "l_csibeamvB", 0, 0);
536 
537  for (double xcsiBeamvB : activeParams.getArray("xcsiBeamvB", {0})) {
538  xcsiBeamvB *= CLHEP::cm;
539  xcsibeamvB[xdimcsiBeamvB] = xcsiBeamvB - x_offset;
540  xdimcsiBeamvB++;
541  }
542  for (double ycsiBeamvB : activeParams.getArray("ycsiBeamvB", {0})) {
543  ycsiBeamvB *= CLHEP::cm;
544  ycsibeamvB[ydimcsiBeamvB] = ycsiBeamvB - y_offset ;
545  ydimcsiBeamvB++;
546  }
547  for (double zcsiBeamvB : activeParams.getArray("zcsiBeamvB", {0})) {
548  zcsiBeamvB *= CLHEP::cm;
549  zcsibeamvB[zdimcsiBeamvB] = zcsiBeamvB ;
550  zdimcsiBeamvB++;
551  }
552  for (int i = 0; i < xdimcsiBeamvB; i++) {
553  PH1SUSTRpos = G4ThreeVector(xcsibeamvB[i], ycsibeamvB[i], zcsibeamvB[i]);
554  new G4PVPlacement(rotX, PH1SUSTRpos, l_csibeamvB, TString::Format("p_csibeamvB_%d", i).Data(), &topVolume, false, 0);
555  }
556 
557  G4double dz_base = activeParams.getLength("lBase") * CLHEP::cm / 2.;
558  G4VSolid* s_base_a = new G4Box("s_base_a", dx_tpcbeam, dy_tpcbeam, dz_base);
559  G4VSolid* s_base_b = new G4Box("s_base_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_base);
560  G4VSolid* s_base = new G4SubtractionSolid("s_base", s_base_a, s_base_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
561  G4LogicalVolume* l_base = new G4LogicalVolume(s_base, geometry::Materials::get("FG_Epoxy") , "l_base", 0, 0);
562 
563  for (double xBase : activeParams.getArray("xBase", {0})) {
564  xBase *= CLHEP::cm;
565  xbase[xdimbase] = xBase - x_offset;
566  xdimbase++;
567  }
568  for (double yBase : activeParams.getArray("yBase", {0})) {
569  yBase *= CLHEP::cm;
570  ybase[ydimbase] = yBase - y_offset;
571  ydimbase++;
572  }
573  for (double zBase : activeParams.getArray("zBase", {0})) {
574  zBase *= CLHEP::cm;
575  zbase[zdimbase] = zBase ;
576  zdimbase++;
577  }
578  for (int i = 0; i < xdimbase; i++) {
579  PH1SUSTRpos = G4ThreeVector(xbase[i], ybase[i], zbase[i]);
580  new G4PVPlacement(0, PH1SUSTRpos, l_base, TString::Format("p_base_%d", i).Data(), &topVolume, false, 0);
581  }
582 
583  G4double dz_bgobeamv = activeParams.getLength("lbgoBeamv") * CLHEP::cm / 2.;
584  G4VSolid* s_bgobeamv_a = new G4Box("s_bgobeamv_a", dx_tpcbeam, dy_tpcbeam, dz_bgobeamv);
585  G4VSolid* s_bgobeamv_b = new G4Box("s_bgobeamv_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_bgobeamv);
586  G4VSolid* s_bgobeamvpos = new G4SubtractionSolid("s_bgobeamvpos", s_bgobeamv_a, s_bgobeamv_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
587  G4VSolid* s_bgobeamvneg = new G4SubtractionSolid("s_bgobeamvneg", s_bgobeamv_a, s_bgobeamv_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
588  G4VSolid* s_bgobeamv = new G4UnionSolid("s_bgobeamv", s_bgobeamvpos, s_bgobeamvneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
589  G4LogicalVolume* l_bgobeamv = new G4LogicalVolume(s_bgobeamv, geometry::Materials::get("FG_Epoxy") , "l_bgobeamv", 0, 0);
590 
591  for (double xbgoBeamv : activeParams.getArray("xbgoBeamv", {0})) {
592  xbgoBeamv *= CLHEP::cm;
593  xbgobeamv[xdimbgobeamv] = xbgoBeamv;
594  xdimbgobeamv++;
595  }
596  for (double ybgoBeamv : activeParams.getArray("ybgoBeamv", {0})) {
597  ybgoBeamv *= CLHEP::cm;
598  ybgobeamv[ydimbgobeamv] = ybgoBeamv - y_offset;
599  ydimbgobeamv++;
600  }
601  for (double zbgoBeamv : activeParams.getArray("zbgoBeamv", {0})) {
602  zbgoBeamv *= CLHEP::cm;
603  zbgobeamv[zdimbgobeamv] = zbgoBeamv ;
604  zdimbgobeamv++;
605  }
606  for (int i = 0; i < xdimbgobeamv; i++) {
607  PH1SUSTRpos = G4ThreeVector(xbgobeamv[i], ybgobeamv[i], zbgobeamv[i]);
608  new G4PVPlacement(rotX, PH1SUSTRpos, l_bgobeamv, TString::Format("p_bgobeamv_%d", i).Data(), &topVolume, false, 0);
609  }
610 
611  G4double dz_bgobeamh = activeParams.getLength("lbgoBeamh") * CLHEP::cm / 2.;
612  G4VSolid* s_bgobeamh_a = new G4Box("s_bgobeamh_a", dx_tpcbeam, dy_tpcbeam, dz_bgobeamh);
613  G4VSolid* s_bgobeamh_b = new G4Box("s_bgobeamh_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_bgobeamh);
614  G4VSolid* s_bgobeamhpos = new G4SubtractionSolid("s_bgobeamhpos", s_bgobeamh_a, s_bgobeamh_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
615  G4VSolid* s_bgobeamhneg = new G4SubtractionSolid("s_bgobeamhneg", s_bgobeamh_a, s_bgobeamh_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
616  G4VSolid* s_bgobeamh = new G4UnionSolid("s_bgobeamh", s_bgobeamhpos, s_bgobeamhneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
617  G4LogicalVolume* l_bgobeamh = new G4LogicalVolume(s_bgobeamh, geometry::Materials::get("FG_Epoxy") , "l_bgobeamh", 0, 0);
618 
619  for (double xbgoBeamh : activeParams.getArray("xbgoBeamh", {0})) {
620  xbgoBeamh *= CLHEP::cm;
621  xbgobeamh[xdimbgobeamh] = xbgoBeamh;
622  xdimbgobeamh++;
623  }
624  for (double ybgoBeamh : activeParams.getArray("ybgoBeamh", {0})) {
625  ybgoBeamh *= CLHEP::cm;
626  ybgobeamh[ydimbgobeamh] = ybgoBeamh - y_offset;
627  ydimbgobeamh++;
628  }
629  for (double zbgoBeamh : activeParams.getArray("zbgoBeamh", {0})) {
630  zbgoBeamh *= CLHEP::cm;
631  zbgobeamh[zdimbgobeamh] = zbgoBeamh ;
632  zdimbgobeamh++;
633  }
634  for (int i = 0; i < xdimbgobeamh; i++) {
635  G4Transform3D transform = G4Translate3D(xbgobeamh[i], ybgobeamh[i],
636  zbgobeamh[i] + 2. * dy_tpcbeam) * G4RotateY3D(90.*CLHEP::deg) * G4RotateZ3D(90.*CLHEP::deg);
637  new G4PVPlacement(transform, l_bgobeamh, TString::Format("p_bgobeamh_%d", i).Data(), &topVolume, false, 0);
638  }
639 
640  G4double dz_bgobeamt = activeParams.getLength("lbgoBeamt") * CLHEP::cm / 2.;
641  G4VSolid* s_bgobeamt_a = new G4Box("s_bgobeamt_a", dx_tpcbeam, dy_tpcbeam, dz_bgobeamt);
642  G4VSolid* s_bgobeamt_b = new G4Box("s_bgobeamt_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_bgobeamt);
643  G4VSolid* s_bgobeamt = new G4SubtractionSolid("s_bgobeamt", s_bgobeamt_a, s_bgobeamt_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
644  G4LogicalVolume* l_bgobeamt = new G4LogicalVolume(s_bgobeamt, geometry::Materials::get("FG_Epoxy") , "l_bgobeamt", 0, 0);
645 
646  for (double xbgoBeamt : activeParams.getArray("xbgoBeamt", {0})) {
647  xbgoBeamt *= CLHEP::cm;
648  xbgobeamt[xdimbgobeamt] = xbgoBeamt;
649  xdimbgobeamt++;
650  }
651  for (double ybgoBeamt : activeParams.getArray("ybgoBeamt", {0})) {
652  ybgoBeamt *= CLHEP::cm;
653  ybgobeamt[ydimbgobeamt] = ybgoBeamt - y_offset;
654  ydimbgobeamt++;
655  }
656  for (double zbgoBeamt : activeParams.getArray("zbgoBeamt", {0})) {
657  zbgoBeamt *= CLHEP::cm;
658  zbgobeamt[zdimbgobeamt] = zbgoBeamt ;
659  zdimbgobeamt++;
660  }
661  for (int i = 0; i < xdimbgobeamt; i++) {
662  PH1SUSTRpos = G4ThreeVector(xbgobeamt[i], ybgobeamt[i], zbgobeamt[i]);
663  new G4PVPlacement(0, PH1SUSTRpos, l_bgobeamt, TString::Format("p_bgobeamt_%d", i).Data(), &topVolume, false, 0);
664  }
665 
666  G4double dz_bgobeamb = activeParams.getLength("lbgoBeamb") * CLHEP::cm / 2.;
667  G4VSolid* s_bgobeamb_a = new G4Box("s_bgobeamb_a", dx_tpcbeam, dy_tpcbeam, dz_bgobeamb);
668  G4VSolid* s_bgobeamb_b = new G4Box("s_bgobeamb_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_bgobeamb);
669  G4VSolid* s_bgobeambpos = new G4SubtractionSolid("s_bgobeambpos", s_bgobeamb_a, s_bgobeamb_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
670  G4VSolid* s_bgobeambneg = new G4SubtractionSolid("s_bgobeambneg", s_bgobeamb_a, s_bgobeamb_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
671  G4VSolid* s_bgobeamb = new G4UnionSolid("s_bgobeamb", s_bgobeambpos, s_bgobeambneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
672  G4LogicalVolume* l_bgobeamb = new G4LogicalVolume(s_bgobeamb, geometry::Materials::get("FG_Epoxy") , "l_bgobeamb", 0, 0);
673  for (double xbgoBeamb : activeParams.getArray("xbgoBeamb", {0})) {
674  xbgoBeamb *= CLHEP::cm;
675  xbgobeamb[xdimbgobeamb] = xbgoBeamb - x_offset;
676  xdimbgobeamb++;
677  }
678  for (double ybgoBeamb : activeParams.getArray("ybgoBeamb", {0})) {
679  ybgoBeamb *= CLHEP::cm;
680  ybgobeamb[ydimbgobeamb] = ybgoBeamb - y_offset;
681  ydimbgobeamb++;
682  }
683  for (double zbgoBeamb : activeParams.getArray("zbgoBeamb", {0})) {
684  zbgoBeamb *= CLHEP::cm;
685  zbgobeamb[zdimbgobeamb] = zbgoBeamb ;
686  zdimbgobeamb++;
687  }
688  for (int i = 0; i < xdimbgobeamb; i++) {
689  G4Transform3D transform = G4Translate3D(xbgobeamb[i], ybgobeamb[i],
690  zbgobeamb[i] + 2. * dy_tpcbeam) * G4RotateY3D(90.*CLHEP::deg) * G4RotateZ3D(90.*CLHEP::deg);
691  new G4PVPlacement(transform, l_bgobeamb, TString::Format("p_bgobeamb_%d", i).Data(), &topVolume, false, 1);
692  }
693 
694  G4double dz_tpcbeamb = activeParams.getLength("ltpcBeamb") * CLHEP::cm / 2.;
695  G4VSolid* s_tpcbeamb_a = new G4Box("s_tpcbeamb_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeamb);
696  G4VSolid* s_tpcbeamb_b = new G4Box("s_tpcbeamb_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeamb);
697  G4VSolid* s_tpcbeambpos = new G4SubtractionSolid("s_tpcbeambpos", s_tpcbeamb_a, s_tpcbeamb_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
698  G4VSolid* s_tpcbeambneg = new G4SubtractionSolid("s_tpcbeambneg", s_tpcbeamb_a, s_tpcbeamb_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
699  G4VSolid* s_tpcbeamb = new G4UnionSolid("s_tpcbeamb", s_tpcbeambpos, s_tpcbeambneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
700  G4LogicalVolume* l_tpcbeamb = new G4LogicalVolume(s_tpcbeamb, geometry::Materials::get("FG_Epoxy") , "l_tpcbeamb", 0, 0);
701  for (double xtpcBeamb : activeParams.getArray("xtpcBeamb", {0})) {
702  xtpcBeamb *= CLHEP::cm;
703  xtpcbeamb[xdimtpcbeamb] = xtpcBeamb - x_offset;
704  xdimtpcbeamb++;
705  }
706  for (double ytpcBeamb : activeParams.getArray("ytpcBeamb", {0})) {
707  ytpcBeamb *= CLHEP::cm;
708  ytpcbeamb[ydimtpcbeamb] = ytpcBeamb - y_offset;
709  ydimtpcbeamb++;
710  }
711  for (double ztpcBeamb : activeParams.getArray("ztpcBeamb", {0})) {
712  ztpcBeamb *= CLHEP::cm;
713  ztpcbeamb[zdimtpcbeamb] = ztpcBeamb ;
714  zdimtpcbeamb++;
715  }
716  for (int i = 0; i < xdimtpcbeamb; i++) {
717  G4Transform3D transform = G4Translate3D(xtpcbeamb[i], ytpcbeamb[i],
718  ztpcbeamb[i] + 2. * dy_tpcbeam) * G4RotateY3D(90.*CLHEP::deg) * G4RotateZ3D(90.*CLHEP::deg);
719  new G4PVPlacement(transform, l_tpcbeamb, TString::Format("p_tpcbeamb_%d", i).Data(), &topVolume, false, 1);
720  }
721 
722  /*
723  G4double dz_tpcbeamhh = activeParams.getLength("ltpcBeamh") * CLHEP::cm / 2.;
724  G4VSolid* s_tpcbeamhh_a = new G4Box("s_tpcbeamhh_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeamhh);
725  G4VSolid* s_tpcbeamhh_b = new G4Box("s_tpcbeamhh_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeamhh);
726  G4VSolid* s_tpcbeamhhpos = new G4SubtractionSolid("s_tpcbeamhhpos", s_tpcbeamhh_a, s_tpcbeamhh_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
727  G4VSolid* s_tpcbeamhhneg = new G4SubtractionSolid("s_tpcbeamhhneg", s_tpcbeamhh_a, s_tpcbeamhh_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
728  G4VSolid* s_tpcbeamhh = new G4UnionSolid("s_tpcbeamhh", s_tpcbeamhhpos, s_tpcbeamhhneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
729  G4LogicalVolume* l_tpcbeamhh = new G4LogicalVolume(s_tpcbeamhh, geometry::Materials::get("FG_Epoxy") , "l_tpcbeamhh", 0, 0);
730 
731  for (double xtpcBeamh : activeParams.getArray("xtpcBeamh", {0})) {
732  xtpcBeamh *= CLHEP::cm;
733  xtpcbeamhh[xdimtpcbeamhh] = xtpcBeamh - x_offset;
734  xdimtpcbeamhh++;
735  }
736  for (double ytpcBeamh : activeParams.getArray("ytpcBeamh", {0})) {
737  ytpcBeamh *= CLHEP::cm;
738  ytpcbeamhh[ydimtpcbeamhh] = ytpcBeamh - y_offset;
739  ydimtpcbeamhh++;
740  }
741  for (double ztpcBeamh : activeParams.getArray("ztpcBeamh", {0})) {
742  ztpcBeamh *= CLHEP::cm;
743  ztpcbeamhh[zdimtpcbeamhh] = ztpcBeamh ;
744  zdimtpcbeamhh++;
745  }
746  for(int i = 0; i < xdimtpcbeamhh; i++) {
747  double z_tmp = 0;
748  if(i<1)z_tmp = ztpcbeamhh[i] + dy_tpcbeam/2.;
749  else z_tmp = ztpcbeamhh[i] - dy_tpcbeam *1.1;
750  G4Transform3D transform = G4Translate3D(xtpcbeamhh[i],ytpcbeamhh[i],z_tmp) * G4RotateY3D(90.*CLHEP::deg) ;
751  new G4PVPlacement(transform, l_tpcbeamhh, TString::Format("p_tpcbeamhh_%d",i).Data(), &topVolume, false, 1);
752  }
753 
754  G4double dz_tpcbeamvv = activeParams.getLength("ltpcBeamv") * CLHEP::cm / 2.;
755  G4VSolid* s_tpcbeamvv_a = new G4Box("s_tpcbeamvv_a", dx_tpcbeam, dy_tpcbeam, dz_tpcbeamvv);
756  G4VSolid* s_tpcbeamvv_b = new G4Box("s_tpcbeamvv_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_tpcbeamvv);
757  G4VSolid* s_tpcbeamvvpos = new G4SubtractionSolid("s_tpcbeamvvpos", s_tpcbeamvv_a, s_tpcbeamvv_b, 0, G4ThreeVector(0, dw_tpcbeam, 0));
758  G4VSolid* s_tpcbeamvvneg = new G4SubtractionSolid("s_tpcbeamvvneg", s_tpcbeamvv_a, s_tpcbeamvv_b, 0, G4ThreeVector(0, -dw_tpcbeam, 0));
759  G4VSolid* s_tpcbeamvv = new G4UnionSolid("s_tpcbeamvv", s_tpcbeamvvpos, s_tpcbeamvvneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam, 0));
760  G4LogicalVolume* l_tpcbeamvv = new G4LogicalVolume(s_tpcbeamvv, geometry::Materials::get("FG_Epoxy") , "l_tpcbeamvv", 0, 0);
761 
762  for (double xtpcBeamv : activeParams.getArray("xtpcBeamv", {0})) {
763  xtpcBeamv *= CLHEP::cm;
764  xtpcbeamvv[xdimtpcbeamvv] = xtpcBeamv - x_offset;
765  xdimtpcbeamvv++;
766  }
767  for (double ytpcBeamv : activeParams.getArray("ytpcBeamv", {0})) {
768  ytpcBeamv *= CLHEP::cm;
769  ytpcbeamvv[ydimtpcbeamvv] = ytpcBeamv - y_offset;
770  ydimtpcbeamvv++;
771  }
772  for (double ztpcBeamv : activeParams.getArray("ztpcBeamv", {0})) {
773  ztpcBeamv *= CLHEP::cm;
774  ztpcbeamvv[zdimtpcbeamvv] = ztpcBeamv ;
775  zdimtpcbeamvv++;
776  }
777  for(int i = 0; i < xdimtpcbeamvv; i++) {
778  //G4Transform3D transform = G4Translate3D(xtpcbeamvv[i],ytpcbeamvv[i],ztpcbeamvv[i]) * G4RotateY3D(90.*CLHEP::deg) ;
779  //new G4PVPlacement(transform, l_tpcbeamvv, TString::Format("p_tpcbeamvv_%d",i).Data(), &topVolume, false, 1);
780  PH1SUSTRpos = G4ThreeVector(xtpcbeamvv[i],ytpcbeamvv[i],ztpcbeamvv[i]);
781  new G4PVPlacement(rotX, PH1SUSTRpos, l_tpcbeamvv, TString::Format("p_tpcbeamvv_%d",i).Data(), &topVolume, false, 0);
782  }
783  */
784  G4double dz_fangsbeamhf = activeParams.getLength("lfangsBeamhf") * CLHEP::cm / 2.;
785  G4VSolid* s_fangsbeamhf_a = new G4Box("s_fangsbeamhf_a", dx_tpcbeam, dy_tpcbeam, dz_fangsbeamhf);
786  G4VSolid* s_fangsbeamhf_b = new G4Box("s_fangsbeamhf_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_fangsbeamhf);
787  G4VSolid* s_fangsbeamhf = new G4SubtractionSolid("s_fangsbeamhf", s_fangsbeamhf_a, s_fangsbeamhf_b, 0, G4ThreeVector(0, dw_tpcbeam,
788  0));
789  G4LogicalVolume* l_fangsbeamhf = new G4LogicalVolume(s_fangsbeamhf, geometry::Materials::get("FG_Epoxy") , "l_fangsbeamhf", 0, 0);
790  double xfangsbeamhf = activeParams.getLength("xfangsBeamhf") * CLHEP::cm;
791  double yfangsbeamhf = activeParams.getLength("yfangsBeamhf") * CLHEP::cm;
792  double zfangsbeamhf = activeParams.getLength("zfangsBeamhf") * CLHEP::cm;
793  PH1SUSTRpos = G4ThreeVector(xfangsbeamhf, yfangsbeamhf, zfangsbeamhf);
794  new G4PVPlacement(rotY, PH1SUSTRpos, l_fangsbeamhf, "p_fangsbeamhf", &topVolume, false, 0);
795 
796  G4double dz_fangsbeamhb = activeParams.getLength("lfangsBeamhb") * CLHEP::cm / 2.;
797  G4VSolid* s_fangsbeamhb_a = new G4Box("s_fangsbeamhb_a", dx_tpcbeam, dy_tpcbeam, dz_fangsbeamhb);
798  G4VSolid* s_fangsbeamhb_b = new G4Box("s_fangsbeamhb_b", dx_tpcbeam - 2.*dw_tpcbeam, dy_tpcbeam - dw_tpcbeam, dz_fangsbeamhb);
799  G4VSolid* s_fangsbeamhbpos = new G4SubtractionSolid("s_fangsbeamhbpos", s_fangsbeamhb_a, s_fangsbeamhb_b, 0, G4ThreeVector(0,
800  dw_tpcbeam, 0));
801  G4VSolid* s_fangsbeamhbneg = new G4SubtractionSolid("s_fangsbeamhbneg", s_fangsbeamhb_a, s_fangsbeamhb_b, 0, G4ThreeVector(0,
802  -dw_tpcbeam, 0));
803  G4VSolid* s_fangsbeamhb = new G4UnionSolid("s_fangsbeamhb", s_fangsbeamhbpos, s_fangsbeamhbneg, 0, G4ThreeVector(0, -2.*dy_tpcbeam,
804  0));
805  G4LogicalVolume* l_fangsbeamhb = new G4LogicalVolume(s_fangsbeamhb, geometry::Materials::get("FG_Epoxy") , "l_fangsbeamhb", 0, 0);
806  double xfangsbeamhb = activeParams.getLength("xfangsBeamhb") * CLHEP::cm;
807  double yfangsbeamhb = activeParams.getLength("yfangsBeamhb") * CLHEP::cm;
808  double zfangsbeamhb = activeParams.getLength("zfangsBeamhb") * CLHEP::cm;
809  TransForm = G4Translate3D(xfangsbeamhb, yfangsbeamhb, zfangsbeamhb) /* G4RotateX3D(90.*CLHEP::deg)*/ * G4RotateY3D(90.*CLHEP::deg);
810  new G4PVPlacement(TransForm, l_fangsbeamhb, "p_fangsbeamhb", &topVolume, false, 1);
811  }
812  }
813  } // ph1sustr namespace
815 } // Belle2 namespace
Belle2::geometry::Materials::get
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:65
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::GearDir
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:41
Belle2::ph1sustr::Ph1sustrFactory
geometry::CreatorFactory< Ph1sustrCreator > Ph1sustrFactory("PH1SUSTRCreator")
Creator creates the phase 1 support structure geometry.
Belle2::gearbox::Interface::getLength
double getLength(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard length unit.
Definition: Interface.h:261
Belle2::gearbox::Interface::getArray
std::vector< double > getArray(const std::string &path) const noexcept(false)
Get the parameter path as a list of double values converted to the standard unit.
Definition: Interface.cc:133
Belle2::ph1sustr::Ph1sustrCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector phase 1 support structure.
Definition: Ph1sustrCreator.h:36
Belle2::PXD::SensitiveDetector
VXD::SensitiveDetector< PXDSimHit, PXDTrueHit > SensitiveDetector
The PXD Sensitive Detector class.
Definition: SensitiveDetector.h:36
Belle2::ph1sustr::Ph1sustrCreator::create
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Function to actually create the geometry, has to be overridden by derived classes.
Definition: Ph1sustrCreator.cc:62
Belle2::geometry::GeometryTypes
GeometryTypes
Flag indiciating the type of geometry to be used.
Definition: GeometryManager.h:39