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