Belle II Software  release-08-01-10
Ph1bpipeCreator.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/ph1bpipe/geometry/Ph1bpipeCreator.h>
10 #include <beast/ph1bpipe/simulation/SensitiveDetector.h>
11 
12 #include <geometry/Materials.h>
13 #include <geometry/CreatorFactory.h>
14 #include <framework/gearbox/GearDir.h>
15 
16 #include <cmath>
17 #include <boost/format.hpp>
18 #include <boost/foreach.hpp>
19 #include <boost/algorithm/string.hpp>
20 
21 #include <G4LogicalVolume.hh>
22 #include <G4PVPlacement.hh>
23 #include <G4VisAttributes.hh>
24 
25 //Shapes
26 #include <G4Box.hh>
27 #include <G4EllipticalTube.hh>
28 #include "G4UnionSolid.hh"
29 #include "G4SubtractionSolid.hh"
30 #include "G4IntersectionSolid.hh"
31 #include <G4UserLimits.hh>
32 #include "G4Tubs.hh"
33 #include "G4Trd.hh"
34 
35 using namespace std;
36 using namespace boost;
37 
38 namespace Belle2 {
45  namespace ph1bpipe {
46 
47  // Register the creator
50 
51  Ph1bpipeCreator::Ph1bpipeCreator(): m_sensitive(0)
52  {
53  //m_sensitive = new SensitiveDetector();
54  }
55 
56  Ph1bpipeCreator::~Ph1bpipeCreator()
57  {
58  if (m_sensitive) delete m_sensitive;
59  }
60 
61  void Ph1bpipeCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
62  {
63 
64  m_sensitive = new SensitiveDetector();
65 
66  //Beam pipes 6 pillars
67  double pillar_height = 61.9 * CLHEP::cm / 2.;
68  double pillar_length = 3. * CLHEP::cm;
69  double pillar_width = pillar_length;
70  G4VSolid* s_bppillar = new G4Box("s_bppillar", pillar_length, pillar_height, pillar_width);
71  G4LogicalVolume* l_bppillar = new G4LogicalVolume(s_bppillar, geometry::Materials::get("Al"), "l_bppillar", 0, 0);
72  G4VisAttributes* white = new G4VisAttributes(G4Colour(1, 1, 1));
73  white->SetForceAuxEdgeVisible(true);
74  l_bppillar->SetVisAttributes(white);
75  double y_offset = -76. * CLHEP::cm + pillar_height;
76  G4ThreeVector Pillarpos = G4ThreeVector(0, y_offset, -154.0 * CLHEP::cm);
77  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar1", &topVolume, false, 1);
78  Pillarpos = G4ThreeVector(0, y_offset, 154.0 * CLHEP::cm);
79  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar2", &topVolume, false, 1);
80  Pillarpos = G4ThreeVector(11.4175236299 * CLHEP::cm, y_offset, -227.758203051 * CLHEP::cm);
81  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar3", &topVolume, false, 1);
82  Pillarpos = G4ThreeVector(-10.4976636985 * CLHEP::cm, y_offset, 227.758203051 * CLHEP::cm);
83  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar4", &topVolume, false, 1);
84  Pillarpos = G4ThreeVector(-8.71250756316 * CLHEP::cm, y_offset, -209.819190072 * CLHEP::cm);
85  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar5", &topVolume, false, 1);
86  Pillarpos = G4ThreeVector(8.71248973173 * CLHEP::cm, y_offset, 209.819190072 * CLHEP::cm);
87  new G4PVPlacement(0, Pillarpos, l_bppillar, "p_bppilar6", &topVolume, false, 1);
88 
89  //Central beam pipe reinforcement
90  double x_reih = 2. * CLHEP::cm / 2.;
91  double y_reih = 2.3 * CLHEP::cm / 2.;
92  double z_reih = 48. * CLHEP::cm / 2.;
93  G4VSolid* s_reih = new G4Box("s_reih", x_reih, y_reih, z_reih);
94  G4LogicalVolume* l_reih = new G4LogicalVolume(s_reih, geometry::Materials::get("Al"), "l_reih", 0, 0);
95  l_reih->SetVisAttributes(white);
96  G4ThreeVector Reihpos = G4ThreeVector(72.8780869619 * CLHEP::mm, 0, 1.35841468498 * CLHEP::mm);
97  new G4PVPlacement(0, Reihpos, l_reih, "p_Reih1", &topVolume, false, 1);
98  Reihpos = G4ThreeVector(-72.8780869619 * CLHEP::mm, 0, 1.35841468498 * CLHEP::mm);
99  new G4PVPlacement(0, Reihpos, l_reih, "p_Reih2", &topVolume, false, 1);
100 
101  double x_reiv = 2. * CLHEP::cm / 2.;
102  double y_reiv = 5.2 * CLHEP::cm / 2.;
103  double z_reiv = 140. * CLHEP::cm / 2.;
104  G4VSolid* s_reiv = new G4Box("s_reiv", x_reiv, y_reiv, z_reiv);
105  G4LogicalVolume* l_reiv = new G4LogicalVolume(s_reiv, geometry::Materials::get("Al"), "l_reiv", 0, 0);
106  l_reiv->SetVisAttributes(white);
107  //G4ThreeVector Reivpos = G4ThreeVector(0, -77.5018052955 * CLHEP::mm, 0);
108  G4ThreeVector Reivpos = G4ThreeVector(0, -83.0 * CLHEP::mm, 0);
109  new G4PVPlacement(0, Reivpos, l_reiv, "p_Reiv1", &topVolume, false, 1);
110  //Reivpos = G4ThreeVector(0, 77.4981947045 * CLHEP::mm,0 );
111  Reivpos = G4ThreeVector(0, 83.0 * CLHEP::mm, 0);
112  new G4PVPlacement(0, Reivpos, l_reiv, "p_Reiv2", &topVolume, false, 1);
113 
114  /*
115  Central beampipe +- pipe_hz = 20cm
116  Flanges of centra BP endcap_hz = 2.2cm thick
117  Distance of the end of xshape from IP xpipe_hz1 = 200cm
118  Length of far beampipe backward xpipePOS_hz = 168.1cm
119  Length of far beampipe backward xpipeMIN_hz = 141.3cm
120  */
121 
122  //lets get the stepsize parameter with a default value of 5 µm
123  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
124 // flag_limitStep = true for SynRad simulation
125  //double stepMax = 50. *CLHEP::mm;
126  //cout << " stepMax = " << stepMax << endl;
127  //bool flag_limitStep = false;
128 // bool flag_limitStep = true;
129 
130 // double stepSize = content.getLength("stepSize", 5*CLHEP::um);
131  double SafetyLength = 0.1 * CLHEP ::cm;
132 
133  double xpipe_hz1 = content.getLength("xpipe_hz1") * CLHEP::cm;
134  double pipe_hz = content.getLength("pipe_hz") * CLHEP::cm;
135  double pipe_outerRadius_x = content.getLength("pipe_outerRadius_x") * CLHEP::cm;
136  double pipe_outerRadius_y = content.getLength("pipe_outerRadius_y") * CLHEP::cm;
137  double pipe_innerRadius_x = content.getLength("pipe_innerRadius_x") * CLHEP::cm;
138  double pipe_innerRadius_y = content.getLength("pipe_innerRadius_y") * CLHEP::cm;
139  double xpipe_innerRadius = content.getLength("xpipe_innerRadius") * CLHEP::cm;
140  double xpipe_outerRadius = content.getLength("xpipe_outerRadius") * CLHEP::cm;
141  double pipe_innerRadiusTiN_x = content.getLength("pipe_innerRadiusTiN_x") * CLHEP::cm;
142  double pipe_innerRadiusTiN_y = content.getLength("pipe_innerRadiusTiN_y") * CLHEP::cm;
143  double xpipe_innerRadiusTiN = content.getLength("xpipe_innerRadiusTiN") * CLHEP::cm;
144  double endcap_hz = content.getLength("endcap_hz") * CLHEP::cm;
145  double endcap_outerRadius = content.getLength("endcap_outerRadius") * CLHEP::cm;
146  double xpipePOS_hz = content.getLength("xpipePOS_hz") * CLHEP::cm;
147  double xpipeMIN_hz = content.getLength("xpipeMIN_hz") * CLHEP::cm;
148 
149  //create ph1bpipe volume
150  double startAngle = 0.*CLHEP::deg;
151  double spanningAngle = 360.*CLHEP::deg;
152  double A1 = 0.0415 * CLHEP::rad;
153 // //cout << " A1 = " << A1 << endl;
154  double OBtemp = sqrt(pow(xpipe_hz1, 2) + pow(xpipe_innerRadius / 2., 2));
155  // cout << " OBtemp = " << OBtemp << endl;
156  double BCtemp = OBtemp * sin(A1);
157  // cout << " BCtemp = " << BCtemp << endl;
158  double ABtemp = sqrt(pow(pipe_hz, 2) + pow(OBtemp, 2) - 2 * pipe_hz * OBtemp * cos(A1));
159  // cout << " ABtemp = " << ABtemp << endl;
160  double A2 = asin(BCtemp / ABtemp) * CLHEP::rad;
161  // cout << " A2 = " << A2 << endl;
162  double xpipe_hz = ABtemp - 2 * endcap_hz / cos(A2); // Length of xshape pipe
163  // cout << " xpipe_hz = " << xpipe_hz << endl;
164  double xcont2 = BCtemp + xpipe_outerRadius; // +Z size of AluCont_f trapezoid
165  double xcont1 = xpipe_outerRadius + 2 * endcap_hz * tan(A2); //-Z size of AluCont_f trapezoid
166  double xcont3 = xcont2 + xpipePOS_hz * sin(A1); // +Z size of AluCont_fp trapezoid
167  double xcont4 = xcont3 - 2 * xpipe_outerRadius; // +Z size of AluCont_fps trapezoid
168  double xcont5 = xcont2 - 2 * xpipe_outerRadius; // -Z size of AluCont_fps trapezoid
169  double xcont3M = xcont2 + xpipeMIN_hz * sin(A1); // +Z size of AluCont_bp trapezoid
170  double xcont4M = xcont3M - 2 * xpipe_outerRadius; // +Z size of AluCont_bps trapezoid
171  // cout << " xcont2 = " << xcont2 << " xcont1 = " << xcont1 << endl;
172  // cout << " xcont3 = " << xcont3 << " xcont4 = " << xcont4 << " xcont5 = " << xcont5 << endl;
173  // cout << " xcont3M = " << xcont3M << " xcont4M = " << xcont4M << endl;
174  double dxtr = 0.5 * (ABtemp + 2 * endcap_hz / cos(A2)) * sin(2 * A2);
175  double dztr = (ABtemp + 2 * endcap_hz / cos(A2)) * sin(A2) * sin(A2);
176  // cout << "dxtr = " << dxtr << " dztr = " << dztr << endl;
177 
178  string matPipe = content.getString("MaterialPipe");
179  string matTiN = content.getString("MaterialTiN");
180  string vacPipe = content.getString("MatVacuum");
181  G4double tubinR = 0.0 * CLHEP::cm;
182 
183  G4VSolid* s_PH1BPIPE = new G4EllipticalTube("s_PH1BPIPE",
184  pipe_outerRadius_x,
185  pipe_outerRadius_y,
186  pipe_hz);
187  G4LogicalVolume* l_PH1BPIPE = new G4LogicalVolume(s_PH1BPIPE, geometry::Materials::get(matPipe), "l_PH1BPIPE", 0, 0);
188  //Lets limit the Geant4 stepsize inside the volume
189  l_PH1BPIPE->SetUserLimits(new G4UserLimits(stepSize));
190  //position central ph1bpipe volume
191  G4ThreeVector PH1BPIPEpos = G4ThreeVector(
192  content.getLength("x_ph1bpipe") * CLHEP::cm,
193  content.getLength("y_ph1bpipe") * CLHEP::cm,
194  content.getLength("z_ph1bpipe") * CLHEP::cm
195  );
196  new G4PVPlacement(0, PH1BPIPEpos, l_PH1BPIPE, "p_PH1BPIPE", &topVolume, false, 0);
197 
198  G4VSolid* s_PH1BPIPETiN = new G4EllipticalTube("s_PH1BPIPETiN",
199  pipe_innerRadius_x,
200  pipe_innerRadius_y,
201  pipe_hz);
202  G4LogicalVolume* l_PH1BPIPETiN = new G4LogicalVolume(s_PH1BPIPETiN, geometry::Materials::get(matTiN), "l_PH1BPIPETiN", 0, 0);
203  new G4PVPlacement(0, PH1BPIPEpos, l_PH1BPIPETiN, "p_PH1BPIPETiN", l_PH1BPIPE, false, 0);
204 
205  G4VSolid* s_PH1BPIPEV = new G4EllipticalTube("s_PH1BPIPEV",
206  pipe_innerRadiusTiN_x,
207  pipe_innerRadiusTiN_y,
208  pipe_hz);
209  G4VSolid* s_PH1BPIPEVac = new G4IntersectionSolid("s_PH1BPIPEVac", s_PH1BPIPE, s_PH1BPIPEV);
210  G4LogicalVolume* l_PH1BPIPEVac = new G4LogicalVolume(s_PH1BPIPEVac, geometry::Materials::get(vacPipe), "l_PH1BPIPEVac", 0, 0);
211  //if (flag_limitStep) l_PH1BPIPEVac->SetUserLimits(new G4UserLimits(stepMax));
212  new G4PVPlacement(0, PH1BPIPEpos, l_PH1BPIPEVac, "p_PH1BPIPEVac", l_PH1BPIPETiN, false, 0);
213 
214  //create central endcaps
215  G4VSolid* s_PH1BPIPEendcap = new G4Tubs("s_PH1BPIPEendcap", 0.,
216  endcap_outerRadius,
217  endcap_hz,
218  startAngle, spanningAngle);
219  G4LogicalVolume* l_PH1BPIPEendcapTop = new G4LogicalVolume(s_PH1BPIPEendcap, geometry::Materials::get(matPipe),
220  "l_PH1BPIPEendcapTop", 0, 0);
221  G4LogicalVolume* l_PH1BPIPEendcapBot = new G4LogicalVolume(s_PH1BPIPEendcap, geometry::Materials::get(matPipe),
222  "l_PH1BPIPEendcapBot", 0, 0);
223  //position central endcaps
224  G4ThreeVector PH1BPIPEendcapposTop = G4ThreeVector(
225  content.getLength("x_ph1bpipe") * CLHEP::cm,
226  content.getLength("y_ph1bpipe") * CLHEP::cm,
227  content.getLength("z_ph1bpipe") * CLHEP::cm + pipe_hz + endcap_hz
228  );
229 
230  G4ThreeVector PH1BPIPEendcapposBot = G4ThreeVector(
231  content.getLength("x_ph1bpipe") * CLHEP::cm,
232  content.getLength("y_ph1bpipe") * CLHEP::cm,
233  content.getLength("z_ph1bpipe") * CLHEP::cm - pipe_hz - endcap_hz
234  );
235  new G4PVPlacement(0, PH1BPIPEendcapposTop, l_PH1BPIPEendcapTop, "p_PH1BPIPEendcapTop", &topVolume, false, 0);
236  new G4PVPlacement(0, PH1BPIPEendcapposBot, l_PH1BPIPEendcapBot, "p_PH1BPIPEendcapBot", &topVolume, false, 0);
237 
238  G4VSolid* s_PH1BPIPEendcapTiN = new G4EllipticalTube("s_PH1BPIPEendcapTiN",
239  pipe_innerRadius_x,
240  pipe_innerRadius_y,
241  endcap_hz);
242  G4LogicalVolume* l_PH1BPIPEendcapTiNTop = new G4LogicalVolume(s_PH1BPIPEendcapTiN, geometry::Materials::get(matTiN),
243  "l_PH1BPIPEendcapTiNTop", 0, 0);
244  G4LogicalVolume* l_PH1BPIPEendcapTiNBot = new G4LogicalVolume(s_PH1BPIPEendcapTiN, geometry::Materials::get(matTiN),
245  "l_PH1BPIPEendcapTiNBot", 0, 0);
246  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_PH1BPIPEendcapTiNTop, "p_PH1BPIPEendcapTiNTop", l_PH1BPIPEendcapTop, false, 0);
247  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_PH1BPIPEendcapTiNBot, "p_PH1BPIPEendcapTiNBot", l_PH1BPIPEendcapBot, false, 0);
248 
249  G4VSolid* s_PH1BPIPEendcapV = new G4EllipticalTube("s_PH1BPIPEendcapV",
250  pipe_innerRadiusTiN_x,
251  pipe_innerRadiusTiN_y,
252  endcap_hz);
253  G4VSolid* s_PH1BPIPEendcapVac = new G4IntersectionSolid("s_PH1BPIPEendcapVac", s_PH1BPIPEendcap, s_PH1BPIPEendcapV);
254  G4LogicalVolume* l_PH1BPIPEendcapVac = new G4LogicalVolume(s_PH1BPIPEendcapVac, geometry::Materials::get(vacPipe),
255  "l_PH1BPIPEendcapVac", 0, 0);
256  //if (flag_limitStep) l_PH1BPIPEendcapVac->SetUserLimits(new G4UserLimits(stepMax));
257  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_PH1BPIPEendcapVac, "p_PH1BPIPEendcapVacTop", l_PH1BPIPEendcapTiNTop, false, 0);
258  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_PH1BPIPEendcapVac, "p_PH1BPIPEendcapVacBot", l_PH1BPIPEendcapTiNBot, false, 0);
259 
260  //create x shape tubes forward
261  // Alu trapezoid containing x-shape tube
262  G4VSolid* AluCont_f = new G4Trd("AluCont_f", xcont1, xcont2, xcont1, xcont1, xpipe_hz / 2.0 * cos(A2));
263  G4LogicalVolume* logi_AluCont_f = new G4LogicalVolume(AluCont_f, geometry::Materials::get(matPipe), "logi_AluCont_f", 0, 0);
264  new G4PVPlacement(0, G4ThreeVector(0, 0, pipe_hz + endcap_hz * 2.0 + 0.5 * xpipe_hz * cos(A2)),
265  logi_AluCont_f, "phys_AluCont_f", &topVolume, false, 0);
266 // //cout << " AluCont_f z pos = " << pipe_hz + endcap_hz*2.0 + xpipe_hz*cos(A2)/2. << endl;
267 // //cout << " AluCont_f zmax pos = " << pipe_hz + endcap_hz*2.0 + xpipe_hz*cos(A2) << endl;
268 
269  // create x part
270  G4VSolid* s_Xshapef_11 = new G4Tubs("s_Xshapef_11",
271  tubinR,
272  xpipe_innerRadius,
273  xpipe_hz,
274  startAngle, spanningAngle);
275  G4VSolid* s_Xshapef_12 = new G4Tubs("s_Xshapef_12",
276  tubinR,
277  xpipe_innerRadius * 2.0,
278  xpipe_hz / 2.0,
279  startAngle, spanningAngle);
280  G4VSolid* s_Xshapef_21 = new G4Tubs("s_Xshapef_21",
281  tubinR,
282  xpipe_innerRadiusTiN,
283  xpipe_hz,
284  startAngle, spanningAngle);
285  G4VSolid* s_Xshapef_22 = new G4Tubs("s_Xshapef_22",
286  tubinR,
287  xpipe_innerRadiusTiN * 2.0,
288  xpipe_hz / 2.0,
289  startAngle, spanningAngle);
290  //create TiN layer on
291 
292  G4VSolid* s_XshapefTiN_11 = new G4Tubs("s_XshapefTiN_11",
293  tubinR,
294  xpipe_innerRadius,
295  xpipe_hz,
296  startAngle, spanningAngle);
297  G4VSolid* s_XshapefTiN_12 = new G4Tubs("s_XshapefTiN_12",
298  tubinR,
299  xpipe_innerRadius * 2.0,
300  xpipe_hz / 2.0,
301  startAngle, spanningAngle);
302  G4VSolid* s_XshapefTiN_21 = new G4Tubs("s_XshapefTiN_21",
303  tubinR,
304  xpipe_innerRadius,
305  xpipe_hz,
306  startAngle, spanningAngle);
307  G4VSolid* s_XshapefTiN_22 = new G4Tubs("s_XshapefTiN_22",
308  tubinR,
309  xpipe_innerRadius * 2.0,
310  xpipe_hz / 2.0,
311  startAngle, spanningAngle);
312  //Slanted tube1
313  G4Transform3D transform_sXshapef_12 = G4Translate3D(0., 0., 0.);
314  transform_sXshapef_12 = transform_sXshapef_12 * G4RotateY3D(-A2);
315  G4IntersectionSolid* Xshapef1 = new G4IntersectionSolid("Xshapef1", s_Xshapef_11, s_Xshapef_12, transform_sXshapef_12);
316  //Slanted tube2
317  G4Transform3D transform_sXshapef_22 = G4Translate3D(0., 0., 0.);
318  transform_sXshapef_22 = transform_sXshapef_22 * G4RotateY3D(A2);
319  G4IntersectionSolid* Xshapef2 = new G4IntersectionSolid("Xshapef2", s_Xshapef_21, s_Xshapef_22, transform_sXshapef_22);
320  //Slanted tube1TiN
321  G4Transform3D transform_sXshapefTiN_12 = G4Translate3D(0., 0., 0.);
322  transform_sXshapefTiN_12 = transform_sXshapefTiN_12 * G4RotateY3D(-A2);
323  G4IntersectionSolid* XshapefTiN1 = new G4IntersectionSolid("XshapefTiN1", s_XshapefTiN_11, s_XshapefTiN_12,
324  transform_sXshapefTiN_12);
325  //Slanted tube2TiN
326  G4Transform3D transform_sXshapefTiN_22 = G4Translate3D(0., 0., 0.);
327  transform_sXshapefTiN_22 = transform_sXshapefTiN_22 * G4RotateY3D(A2);
328  G4IntersectionSolid* XshapefTiN2 = new G4IntersectionSolid("XshapefTiN2", s_XshapefTiN_21, s_XshapefTiN_22,
329  transform_sXshapefTiN_22);
330 
331  // (sXshapeTiN2 + sXshapeTiN1)
332  G4Transform3D transform_XshapefTiN1 = G4Translate3D(dxtr, 0., -dztr);
333  transform_XshapefTiN1 = transform_XshapefTiN1 * G4RotateY3D(2.*A2);
334  G4UnionSolid* XshapeTiNForwx = new G4UnionSolid("XshapeTiNForwx", XshapefTiN2, XshapefTiN1, transform_XshapefTiN1);
335  // Place XshapeTiNForwx into logi_AluCont_f
336  G4RotationMatrix* yRot = new G4RotationMatrix;
337  yRot->rotateY(A2);
338  G4ThreeVector transform_XshapeTiNForwx(-dxtr / 2., 0., 0.);
339  G4IntersectionSolid* XshapeTiNForw = new G4IntersectionSolid("XshapeTiNForw", AluCont_f, XshapeTiNForwx, yRot,
340  transform_XshapeTiNForwx);
341  G4LogicalVolume* l_XshapeTiNForw = new G4LogicalVolume(XshapeTiNForw, geometry::Materials::get(matTiN), "l_XshapeTiNForw", 0, 0);
342  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapeTiNForw, "p_XshapeTiNForw", logi_AluCont_f, false, 0);
343 
344  // sXshape2 + sXshape1
345  G4Transform3D transform_Xshapef1 = G4Translate3D(dxtr, 0., -dztr);
346  transform_Xshapef1 = transform_Xshapef1 * G4RotateY3D(2.*A2);
347  G4UnionSolid* XshapeForwx = new G4UnionSolid("XshapeForwx", Xshapef2, Xshapef1, transform_Xshapef1);
348  // Place XshapeForwx into l_XshapeTiNForw
349  G4Transform3D transform_XshapeForwx = G4Translate3D(-dxtr / 2., 0., 0.);
350  transform_XshapeForwx = transform_XshapeForwx * G4RotateY3D(-A2);
351  G4IntersectionSolid* XshapeForw = new G4IntersectionSolid("XshapeForw", AluCont_f, XshapeForwx, transform_XshapeForwx);
352  G4LogicalVolume* l_XshapeForw = new G4LogicalVolume(XshapeForw, geometry::Materials::get(vacPipe), "l_XshapeForw", 0, 0);
353  //if (flag_limitStep) l_XshapeForw ->SetUserLimits(new G4UserLimits(stepMax));
354  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapeForw, "p_XshapeForw", l_XshapeTiNForw, false, 0);
355 
356  //create far forward parts
357  // Alu trapezoid containing far parts
358  G4VSolid* AluCont_fpO = new G4Trd("AluCont_fpO", xcont2, xcont3, xcont1, xcont1, xpipePOS_hz / 2.0 * cos(A1));
359  G4VSolid* AluCont_fps = new G4Trd("AluCont_fps", xcont5, xcont4, xcont1, xcont1, xpipePOS_hz / 2.0 * cos(A1));
360  G4Transform3D transform_AluCont_fps = G4Translate3D(0., 0., 0.);
361  G4SubtractionSolid* AluCont_fp = new G4SubtractionSolid("AluCont_fp", AluCont_fpO, AluCont_fps, transform_AluCont_fps);
362  G4LogicalVolume* logi_AluCont_fp = new G4LogicalVolume(AluCont_fp, geometry::Materials::get(matPipe), "logi_AluCont_fp", 0, 0);
363  new G4PVPlacement(0, G4ThreeVector(0, 0, pipe_hz + endcap_hz * 2.0 + xpipe_hz * cos(A2) + xpipePOS_hz / 2.0 * cos(A1)),
364  logi_AluCont_fp, "phys_AluCont_fp", &topVolume, false, 0);
365  //HER far forward part
366  //slanted tubes Vacuum
367  G4VSolid* s_XshapefPOS_11 = new G4Tubs("s_XshapefPOS11",
368  tubinR,
369  xpipe_innerRadius,
370  xpipePOS_hz,
371  startAngle, spanningAngle);
372  G4VSolid* s_XshapefPOS_12 = new G4Tubs("s_XshapefPOS12",
373  tubinR,
374  xpipe_innerRadius * 2,
375  xpipePOS_hz / 2. - 0.7568202457 * SafetyLength,
376  startAngle, spanningAngle);
377  G4Transform3D transform_sXshapefPOS_12 = G4Translate3D(0., 0., 0.);
378  transform_sXshapefPOS_12 = transform_sXshapefPOS_12 * G4RotateY3D(-A1);
379  G4IntersectionSolid* XshapefPOS1 = new G4IntersectionSolid("XshapefPOS1", s_XshapefPOS_11, s_XshapefPOS_12,
380  transform_sXshapefPOS_12);
381  //put HER vacuum part in AluCont_fp
382  G4LogicalVolume* l_XshapefPOS1 = new G4LogicalVolume(XshapefPOS1, geometry::Materials::get(vacPipe), "l_XshapefPOS1", 0, 0);
383  //if (flag_limitStep) l_XshapefPOS1 ->SetUserLimits(new G4UserLimits(stepMax));
384  G4RotationMatrix* yRotP = new G4RotationMatrix;
385  yRotP->rotateY(-A1);
386  new G4PVPlacement(yRotP, G4ThreeVector(2 * endcap_hz * tan(A2) + xpipe_hz * sin(A2) + xpipePOS_hz / 2.*sin(A1), 0., 0.),
387  l_XshapefPOS1, "p_XshapefPOS1", logi_AluCont_fp, false, 0);
388  //LER far forward part
389  //slanted tubes TiN
390  G4VSolid* s_XshapefPOS_21 = new G4Tubs("s_XshapefPOS21",
391  tubinR,
392  xpipe_innerRadius,
393  xpipePOS_hz,
394  startAngle, spanningAngle);
395  G4VSolid* s_XshapefPOS_22 = new G4Tubs("s_XshapefPOS22",
396  tubinR,
397  xpipe_innerRadius * 2,
398  xpipePOS_hz / 2. - 0.7568202457 * SafetyLength,
399  startAngle, spanningAngle);
400  G4Transform3D transform_sXshapefPOS_22 = G4Translate3D(0., 0., 0.);
401  transform_sXshapefPOS_22 = transform_sXshapefPOS_22 * G4RotateY3D(A1);
402  G4IntersectionSolid* XshapefPOS2 = new G4IntersectionSolid("XshapefPOS2", s_XshapefPOS_21, s_XshapefPOS_22,
403  transform_sXshapefPOS_22);
404  //put HER TiN part in AluCont_fp
405  G4LogicalVolume* l_XshapefPOS2 = new G4LogicalVolume(XshapefPOS2, geometry::Materials::get(matTiN), "l_XshapefPOS2", 0, 0);
406  G4RotationMatrix* yRotM = new G4RotationMatrix;
407  yRotM->rotateY(A1);
408  new G4PVPlacement(yRotM, G4ThreeVector(-(2 * endcap_hz * tan(A2) + xpipe_hz * sin(A2) + xpipePOS_hz / 2.*sin(A1)), 0., 0.),
409  l_XshapefPOS2, "p_XshapefPOS2", logi_AluCont_fp, false, 0);
410  //slanted tubes Vacuum
411  G4VSolid* s_XshapefPOS_31 = new G4Tubs("s_XshapefPOS31",
412  tubinR,
413  xpipe_innerRadiusTiN,
414  xpipePOS_hz,
415  startAngle, spanningAngle);
416  G4VSolid* s_XshapefPOS_32 = new G4Tubs("s_XshapefPOS32",
417  tubinR,
418  xpipe_innerRadiusTiN * 2,
419  xpipePOS_hz / 2. - 0.7568202457 * SafetyLength,
420  startAngle, spanningAngle);
421  G4Transform3D transform_sXshapefPOS_32 = G4Translate3D(0., 0., 0.);
422  transform_sXshapefPOS_32 = transform_sXshapefPOS_32 * G4RotateY3D(A1);
423  G4IntersectionSolid* XshapefPOS3 = new G4IntersectionSolid("XshapefPOS3", s_XshapefPOS_31, s_XshapefPOS_32,
424  transform_sXshapefPOS_32);
425  //put HER vacuum part into l_XshapefPOS2
426  G4LogicalVolume* l_XshapefPOS3 = new G4LogicalVolume(XshapefPOS3, geometry::Materials::get(vacPipe), "l_XshapefPOS3", 0, 0);
427  //if (flag_limitStep) l_XshapefPOS3 ->SetUserLimits(new G4UserLimits(stepMax));
428  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapefPOS3, "p_XshapefPOS3", l_XshapefPOS2, false, 0);
429 
430  // create x-shape tubes backward
431  // Alu trapezoid containing x-shape tube
432  G4VSolid* AluCont_b = new G4Trd("AluCont_b", xcont2, xcont1, xcont1, xcont1, xpipe_hz / 2.0 * cos(A2));
433  G4LogicalVolume* logi_AluCont_b = new G4LogicalVolume(AluCont_b, geometry::Materials::get(matPipe), "logi_AluCont_b", 0, 0);
434  new G4PVPlacement(0, G4ThreeVector(0, 0, -pipe_hz - endcap_hz * 2.0 - xpipe_hz * cos(A2) / 2.), logi_AluCont_b, "phys_AluCont_b",
435  &topVolume, false, 1);
436  //cout << " logi_AluCont_b z = " << -pipe_hz - endcap_hz*2.0 - xpipe_hz*cos(A2)/2. << endl;
437  //cout << " logi_AluCont_b zmax = " << -pipe_hz - endcap_hz*2.0 - xpipe_hz*cos(A2) << endl;
438 
439  // create x part
440  G4VSolid* s_Xshapeb_11 = new G4Tubs("s_Xshapeb_11",
441  tubinR,
442  xpipe_innerRadiusTiN,
443  xpipe_hz,
444  startAngle, spanningAngle);
445  G4VSolid* s_Xshapeb_12 = new G4Tubs("s_Xshapeb_12",
446  tubinR,
447  xpipe_innerRadiusTiN * 2.0,
448  xpipe_hz / 2.0,
449  startAngle, spanningAngle);
450  G4VSolid* s_Xshapeb_21 = new G4Tubs("s_Xshapeb_21",
451  tubinR,
452  xpipe_innerRadius,
453  xpipe_hz,
454  startAngle, spanningAngle);
455  G4VSolid* s_Xshapeb_22 = new G4Tubs("s_Xshapeb_22",
456  tubinR,
457  xpipe_innerRadius * 2.0,
458  xpipe_hz / 2.0,
459  startAngle, spanningAngle);
460  //create TiN layer on
461  G4VSolid* s_XshapebTiN_11 = new G4Tubs("s_XshapebTiN_11",
462  tubinR,
463  xpipe_innerRadius,
464  xpipe_hz,
465  startAngle, spanningAngle);
466  G4VSolid* s_XshapebTiN_12 = new G4Tubs("s_XshapebTiN_12",
467  tubinR,
468  xpipe_innerRadius * 2.0,
469  xpipe_hz / 2.0,
470  startAngle, spanningAngle);
471  G4VSolid* s_XshapebTiN_21 = new G4Tubs("s_XshapebTiN_21",
472  tubinR,
473  xpipe_innerRadius,
474  xpipe_hz,
475  startAngle, spanningAngle);
476  G4VSolid* s_XshapebTiN_22 = new G4Tubs("s_XshapebTiN_22",
477  tubinR,
478  xpipe_innerRadius * 2.0,
479  xpipe_hz / 2.0,
480  startAngle, spanningAngle);
481  //Slanted tube1
482  G4Transform3D transform_sXshapeb_12 = G4Translate3D(0., 0., 0.);
483  transform_sXshapeb_12 = transform_sXshapeb_12 * G4RotateY3D(A2);
484  G4IntersectionSolid* Xshapeb1 = new G4IntersectionSolid("Xshapeb1", s_Xshapeb_11, s_Xshapeb_12, transform_sXshapeb_12);
485  //Slanted tube2
486  G4Transform3D transform_sXshapeb_22 = G4Translate3D(0., 0., 0.);
487  transform_sXshapeb_22 = transform_sXshapeb_22 * G4RotateY3D(-A2);
488  G4IntersectionSolid* Xshapeb2 = new G4IntersectionSolid("Xshapeb2", s_Xshapeb_21, s_Xshapeb_22, transform_sXshapeb_22);
489  //Slanted tube1TiN
490  G4Transform3D transform_sXshapebTiN_12 = G4Translate3D(0., 0., 0.);
491  transform_sXshapebTiN_12 = transform_sXshapebTiN_12 * G4RotateY3D(A2);
492  G4IntersectionSolid* XshapebTiN1 = new G4IntersectionSolid("XshapebTiN1", s_XshapebTiN_11, s_XshapebTiN_12,
493  transform_sXshapebTiN_12);
494  //Slanted tube2TiN
495  G4Transform3D transform_sXshapebTiN_22 = G4Translate3D(0., 0., 0.);
496  transform_sXshapebTiN_22 = transform_sXshapebTiN_22 * G4RotateY3D(-A2);
497  G4IntersectionSolid* XshapebTiN2 = new G4IntersectionSolid("XshapebTiN2", s_XshapebTiN_21, s_XshapebTiN_22,
498  transform_sXshapebTiN_22);
499 
500  // (sXshapeTiN2 + sXshapeTiN1)
501  G4Transform3D transform_XshapebTiN1 = G4Translate3D(dxtr, 0., dztr);
502  //cout << "dxtr= " << dxtr << " dztr= " << dztr << endl;;
503  transform_XshapebTiN1 = transform_XshapebTiN1 * G4RotateY3D(-2.*A2);
504  G4UnionSolid* XshapeTiNBackwx = new G4UnionSolid("XshapeTiNBackwx", XshapebTiN2, XshapebTiN1, transform_XshapebTiN1);
505  // Place XshapeTiNBackwx into logi_AluCont_b
506  G4Transform3D transform_XshapeTiNBackwx = G4Translate3D(-dxtr / 2., 0., 0.);
507  transform_XshapeTiNBackwx = transform_XshapeTiNBackwx * G4RotateY3D(A2);
508  G4IntersectionSolid* XshapeTiNBackw = new G4IntersectionSolid("XshapeTiNBackw", AluCont_b, XshapeTiNBackwx,
509  transform_XshapeTiNBackwx);
510  G4LogicalVolume* l_XshapeTiNBackw = new G4LogicalVolume(XshapeTiNBackw, geometry::Materials::get(matTiN), "l_XshapeTiNBackw", 0,
511  0);
512  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapeTiNBackw, "p_XshapeTiNBackw", logi_AluCont_b, false, 0);
513 
514  // sXshape2 + sXshape1
515  G4Transform3D transform_Xshapeb1 = G4Translate3D(dxtr, 0., dztr);
516  transform_Xshapeb1 = transform_Xshapeb1 * G4RotateY3D(-2.*A2);
517  G4UnionSolid* XshapeBackwx = new G4UnionSolid("XshapebBackwx", Xshapeb2, Xshapeb1, transform_Xshapeb1);
518  G4Transform3D transform_XshapeBackwx = G4Translate3D(-dxtr / 2., 0., 0.);
519  transform_XshapeBackwx = transform_XshapeBackwx * G4RotateY3D(A2);
520  G4IntersectionSolid* XshapeBackw = new G4IntersectionSolid("XshapeBackw", AluCont_b, XshapeBackwx, transform_XshapeBackwx);
521  G4LogicalVolume* l_XshapeBackw = new G4LogicalVolume(XshapeBackw, geometry::Materials::get(vacPipe), "l_XshapeBackw", 0, 0);
522  //if (flag_limitStep) l_XshapeBackw->SetUserLimits(new G4UserLimits(stepMax));
523  // Place XshapeBackw
524  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapeBackw, "p_XshapeBackw", l_XshapeTiNBackw, false, 0);
525 
526  //create far backward parts
527  // Alu trapezoid containing far parts
528  G4VSolid* AluCont_bpO = new G4Trd("AluCont_bpO", xcont3M, xcont2, xcont1, xcont1, xpipeMIN_hz / 2.0 * cos(A1));
529  G4VSolid* AluCont_bps = new G4Trd("AluCont_bps", xcont4M, xcont5, xcont1, xcont1, xpipeMIN_hz / 2.0 * cos(A1));
530  G4Transform3D transform_AluCont_bps = G4Translate3D(0., 0., 0.);
531  G4SubtractionSolid* AluCont_bp = new G4SubtractionSolid("AluCont_bp", AluCont_bpO, AluCont_bps, transform_AluCont_bps);
532  G4LogicalVolume* logi_AluCont_bp = new G4LogicalVolume(AluCont_bp, geometry::Materials::get(matPipe), "logi_AluCont_bp", 0, 0);
533  new G4PVPlacement(0, G4ThreeVector(0, 0, -(pipe_hz + endcap_hz * 2.0 + xpipe_hz * cos(A2) + xpipeMIN_hz / 2.0 * cos(A1))),
534  logi_AluCont_bp, "phys_AluCont_bp", &topVolume, false, 0);
535  //HER far backward part
536  //slanted tubes Vacuum
537  G4VSolid* s_XshapebMIN_11 = new G4Tubs("s_XshapebMIN11",
538  tubinR,
539  xpipe_innerRadius,
540  xpipeMIN_hz,
541  startAngle, spanningAngle);
542  G4VSolid* s_XshapebMIN_12 = new G4Tubs("s_XshapebMIN12",
543  tubinR,
544  xpipe_innerRadius * 2,
545  xpipeMIN_hz / 2. - 0.525641366 * SafetyLength,
546  startAngle, spanningAngle);
547  G4Transform3D transform_sXshapebMIN_12 = G4Translate3D(0., 0., 0.);
548  transform_sXshapebMIN_12 = transform_sXshapebMIN_12 * G4RotateY3D(-A1);
549  G4IntersectionSolid* XshapebMIN1 = new G4IntersectionSolid("XshapebMIN1", s_XshapebMIN_11, s_XshapebMIN_12,
550  transform_sXshapebMIN_12);
551  //put HER vacuum part in AluCont_bp
552  G4LogicalVolume* l_XshapebMIN1 = new G4LogicalVolume(XshapebMIN1, geometry::Materials::get(vacPipe), "l_XshapebMIN1", 0, 0);
553  //if (flag_limitStep) l_XshapebMIN1 ->SetUserLimits(new G4UserLimits(stepMax));
554  new G4PVPlacement(yRotP, G4ThreeVector(-(2 * endcap_hz * tan(A2) + xpipe_hz * sin(A2) + xpipeMIN_hz / 2.*sin(A1)), 0., 0.),
555  l_XshapebMIN1, "p_XshapebMIN1", logi_AluCont_bp, false, 0);
556  //LER far backward part
557  //slanted tubes TiN
558  G4VSolid* s_XshapebMIN_21 = new G4Tubs("s_XshapebMIN21",
559  tubinR,
560  xpipe_innerRadius,
561  xpipeMIN_hz,
562  startAngle, spanningAngle);
563  G4VSolid* s_XshapebMIN_22 = new G4Tubs("s_XshapebMIN22",
564  tubinR,
565  xpipe_innerRadius * 2,
566  xpipeMIN_hz / 2. - 0.525641366 * SafetyLength,
567  startAngle, spanningAngle);
568  G4Transform3D transform_sXshapebMIN_22 = G4Translate3D(0., 0., 0.);
569  transform_sXshapebMIN_22 = transform_sXshapebMIN_22 * G4RotateY3D(A1);
570  G4IntersectionSolid* XshapebMIN2 = new G4IntersectionSolid("XshapebMIN2", s_XshapebMIN_21, s_XshapebMIN_22,
571  transform_sXshapebMIN_22);
572  //put HER TiN part in AluCont_bp
573  G4LogicalVolume* l_XshapebMIN2 = new G4LogicalVolume(XshapebMIN2, geometry::Materials::get(matTiN), "l_XshapebMIN2", 0, 0);
574  new G4PVPlacement(yRotM, G4ThreeVector(2 * endcap_hz * tan(A2) + xpipe_hz * sin(A2) + xpipeMIN_hz / 2.*sin(A1), 0., 0.),
575  l_XshapebMIN2, "p_XshapebMIN2", logi_AluCont_bp, false, 0);
576  //slanted tubes Vacuum
577  G4VSolid* s_XshapebMIN_31 = new G4Tubs("s_XshapebMIN31",
578  tubinR,
579  xpipe_innerRadiusTiN,
580  xpipeMIN_hz,
581  startAngle, spanningAngle);
582  G4VSolid* s_XshapebMIN_32 = new G4Tubs("s_XshapebMIN32",
583  tubinR,
584  xpipe_innerRadiusTiN * 2,
585  xpipeMIN_hz / 2. - 0.525641366 * SafetyLength,
586  startAngle, spanningAngle);
587  G4Transform3D transform_sXshapebMIN_32 = G4Translate3D(0., 0., 0.);
588  transform_sXshapebMIN_32 = transform_sXshapebMIN_32 * G4RotateY3D(A1);
589  G4IntersectionSolid* XshapebMIN3 = new G4IntersectionSolid("XshapebMIN3", s_XshapebMIN_31, s_XshapebMIN_32,
590  transform_sXshapebMIN_32);
591  //put HER vacuum part in AluCont_bp
592  G4LogicalVolume* l_XshapebMIN3 = new G4LogicalVolume(XshapebMIN3, geometry::Materials::get(vacPipe), "l_XshapebMIN3", 0, 0);
593  //if (flag_limitStep) l_XshapebMIN3 ->SetUserLimits(new G4UserLimits(stepMax));
594  new G4PVPlacement(0, G4ThreeVector(0., 0., 0.), l_XshapebMIN3, "p_XshapebMIN3", l_XshapebMIN2, false, 0);
595 
596  }
597  } // ph1bpipe namespace
599 } // Belle2 namespace
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
Sensitive Detector implementation of the PH1BPIPE detector.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
double tan(double a)
tan for double
Definition: beamHelpers.h:31
GeometryTypes
Flag indiciating the type of geometry to be used.
geometry::CreatorFactory< Ph1bpipeCreator > Ph1bpipeFactory("PH1BPIPECreator")
Creator creates the phase 1 beam pipe for |s| < 4 m geometry.
Abstract base class for different kinds of events.
Very simple class to provide an easy way to register creators with the CreatorManager.