Analysis Software
Documentation for sPHENIX simulation software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4BbcDetector.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4BbcDetector.cc
1 #include "PHG4BbcDetector.h"
2 
3 #include "PHG4BbcDisplayAction.h"
4 
5 #include <phparameter/PHParameters.h>
6 
7 #include <g4main/PHG4Detector.h>
9 #include <g4main/PHG4Subsystem.h>
10 
11 #include <phool/recoConsts.h>
12 
13 #include <Geant4/G4LogicalVolume.hh>
14 #include <Geant4/G4Material.hh>
15 #include <Geant4/G4NistManager.hh>
16 #include <Geant4/G4PVPlacement.hh>
17 #include <Geant4/G4Polyhedra.hh>
18 #include <Geant4/G4RotationMatrix.hh> // for G4RotationMatrix
19 #include <Geant4/G4String.hh> // for G4String
20 #include <Geant4/G4SystemOfUnits.hh>
21 #include <Geant4/G4ThreeVector.hh>
22 #include <Geant4/G4Tubs.hh>
23 #include <Geant4/G4Box.hh>
24 #include <Geant4/G4SubtractionSolid.hh>
25 #include <Geant4/G4Trap.hh>
26 
27 #include <Geant4/G4Types.hh> // for G4double, G4int
28 #include <Geant4/G4VPhysicalVolume.hh>
29 
30 #include <cmath>
31 #include <vector>
32 #include <iostream> // for operator<<, endl, bas...
33 
34 
35 class PHCompositeNode;
36 
38  : PHG4Detector(subsys, Node, dnam)
39  , m_DisplayAction(dynamic_cast<PHG4BbcDisplayAction *>(subsys->GetDisplayAction()))
40  , m_Params(params)
41  , m_ActiveFlag(m_Params->get_int_param("active"))
42  , m_SupportActiveFlag(m_Params->get_int_param("supportactive"))
43  , front_bbcz(248*cm)
44 {
45 }
46 
47 //_______________________________________________________________
48 //_______________________________________________________________
50 {
51  G4LogicalVolume *mylogvol = volume->GetLogicalVolume();
52 
53  if (m_ActiveFlag)
54  {
55  if (m_PhysLogicalVolSet.find(mylogvol) != m_PhysLogicalVolSet.end())
56  {
57  return 1;
58  }
59  }
61  {
62  if (m_SupportLogicalVolSet.find(mylogvol) != m_SupportLogicalVolSet.end())
63  {
64  return -2;
65  }
66  }
67  return 0;
68 }
69 
70 void PHG4BbcDetector::ConstructMe(G4LogicalVolume *logicWorld)
71 {
72  //std::cout << "PHG4BbcDetector::ConstructMe()" << std::endl;
73 
75 
76  //========================
77  // Build a single BBC PMT
78  //========================
79 
80  // BBC Detector Tube (BBCD) Logical Volume (contains all parts of PMT+Radiator)
81  G4Material *WorldMaterial = GetDetectorMaterial(rc->get_StringFlag("WorldMaterial"));
82  G4double z_bbcd[] = {-6.15 * cm, 6.15 * cm};
83  G4double len_bbcd = z_bbcd[1] - z_bbcd[0];
84  G4double rin_bbcd[] = {0 * cm, 0 * cm};
85  G4double rout_bbcd[] = {1.4 * cm, 1.4 * cm};
86  G4Polyhedra *bbcd = new G4Polyhedra("bbcd", 0., 2 * M_PI, 6, 2, z_bbcd, rin_bbcd, rout_bbcd);
87  G4LogicalVolume *bbcd_lv = new G4LogicalVolume(bbcd, WorldMaterial, G4String("Bbc_tube"));
88 
89  //
90  // Place the BBCA, BBCQ, BBCP, BBCR and BBCH in to BBCD.
91  //
92 
93  // BBC Attachment Plate (BBCA)
94  const G4double z_bbca[] = {-0.5 * cm, -0.201 * cm, -0.2 * cm, 0.5 * cm};
95  const G4double rInner_bbca[] = {0.2 * cm, 0.2 * cm, 0.2 * cm, 0.2 * cm};
96  const G4double rOuter_bbca[] = {1.4 * cm, 1.4 * cm, 1.375 * cm, 1.375 * cm};
97  G4Polyhedra *bbca = new G4Polyhedra("bbca", 0., 2 * M_PI, 6, 4, z_bbca, rInner_bbca, rOuter_bbca);
98 
99  G4Material *Aluminum = GetDetectorMaterial("G4_Al");
100 
101  G4LogicalVolume *bbca_lv = new G4LogicalVolume(bbca, Aluminum, G4String("Bbc_attach_plate"));
102  m_SupportLogicalVolSet.insert(bbca_lv);
103  GetDisplayAction()->AddVolume(bbca_lv, "Bbc_attach_plate");
104 
105  G4double xpos = 0. * cm;
106  G4double ypos = 0. * cm;
107  G4double len_bbca = z_bbca[3] - z_bbca[0];
108  G4double zpos = z_bbcd[0] + len_bbca * 0.5;
109 
110  G4VPhysicalVolume *bbca_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbca_lv, "BBCA", bbcd_lv, false, 0);
111  if (!bbca_phys) // more to prevent compiler warnings
112  {
113  std::cout << "placement of BBCA failed" << std::endl;
114  }
115 
116  G4Material *Quartz = GetDetectorMaterial("G4_SILICON_DIOXIDE");
117 
118  // BBC Quartz Radiator
119  const G4double z_bbcq[] = {-1.5 * cm, 1.5 * cm};
120  const G4double rInner_bbcq[] = {0., 0.};
121  const G4double rOuter_bbcq[] = {1.27 * cm, 1.27 * cm};
122  G4Polyhedra *bbcq = new G4Polyhedra("bbcq", 0., 2 * M_PI, 6, 2, z_bbcq, rInner_bbcq, rOuter_bbcq);
123 
124  G4LogicalVolume *bbcq_lv = new G4LogicalVolume(bbcq, Quartz, G4String("Bbc_quartz"));
125  GetDisplayAction()->AddVolume(bbcq_lv, "Bbc_quartz");
126 
127  xpos = 0. * cm;
128  ypos = 0. * cm;
129  G4double len_bbcq = z_bbcq[1] - z_bbcq[0];
130  zpos += len_bbca * 0.5 + len_bbcq * 0.5;
131 
132  G4VPhysicalVolume *bbcq_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcq_lv, "BBCQ", bbcd_lv, false, 0);
133  if (!bbcq_phys) // more to prevent compiler warnings
134  {
135  std::cout << "placement of BBCQ failed" << std::endl;
136  }
137 
138  // BBC PMT
139  const G4double len_bbcp = 4.4 * cm;
140  const G4double rInner_bbcp = 1.09 * cm;
141  const G4double rOuter_bbcp = 1.29 * cm;
142  G4Tubs *bbcp = new G4Tubs("bbcp", rInner_bbcp, rOuter_bbcp, len_bbcp * 0.5, 0 * deg, 360 * deg);
143 
144  G4LogicalVolume *bbcp_lv = new G4LogicalVolume(bbcp, Quartz, G4String("Bbc_PMT"));
145  GetDisplayAction()->AddVolume(bbcp_lv, "Bbc_PMT");
146 
147  xpos = 0. * cm;
148  ypos = 0. * cm;
149  zpos += len_bbcq * 0.5 + len_bbcp * 0.5;
150 
151  G4VPhysicalVolume *bbcp_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcp_lv, "BBCP", bbcd_lv, false, 0);
152  if (!bbcp_phys) // more to prevent compiler warnings
153  {
154  std::cout << "placement of BBCP failed" << std::endl;
155  }
156 
157  // BBC Breeder Module
158  const G4double len_bbcr = 3.9 * cm;
159  const G4double rInner_bbcr = 0.0 * cm;
160  const G4double rOuter_bbcr = 1.2 * cm;
161  G4Tubs *bbcr = new G4Tubs("bbcr", rInner_bbcr, rOuter_bbcr, len_bbcr * 0.5, 0 * deg, 360 * deg);
162 
163  G4int natoms;
164  G4int ncomponents;
165  G4double density;
166  G4Material *G10 = new G4Material("BBC_G10", density = 1.700 * g / cm3, ncomponents = 4);
167  G4NistManager *manager = G4NistManager::Instance();
168  G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("Si"), natoms = 1);
169  G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("O"), natoms = 2);
170  G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("C"), natoms = 3);
171  G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("H"), natoms = 3);
172 
173  G4LogicalVolume *bbcr_lv = new G4LogicalVolume(bbcr, G10, "Bbc_Breeder_Module");
174  GetDisplayAction()->AddVolume(bbcr_lv, "Bbc_Breeder_Module");
175 
176  xpos = 0. * cm;
177  ypos = 0. * cm;
178  zpos += len_bbcp * 0.5 + len_bbcr * 0.5;
179 
180  G4PVPlacement *plcmt = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcr_lv, "BBCR", bbcd_lv, false, 0);
181  if (!plcmt) // more to prevent compiler warnings
182  {
183  std::cout << "placement of BBCR failed" << std::endl;
184  }
185 
186  // BBC Mu Metal Shield
187  const G4double z_bbch[] = {-5.65 * cm, 5.65 * cm};
188  const G4double rInner_bbch[] = {1.375 * cm, 1.375 * cm};
189  const G4double rOuter_bbch[] = {1.4 * cm, 1.4 * cm};
190  G4Polyhedra *bbch = new G4Polyhedra("bbch", 0., 2 * M_PI, 6, 2, z_bbch, rInner_bbch, rOuter_bbch);
191 
192  G4Material *MuMetal = manager->FindOrBuildMaterial("G4_STAINLESS-STEEL");
193 
194  G4LogicalVolume *bbch_lv = new G4LogicalVolume(bbch, MuMetal, G4String("Bbc_Shield"));
195 
196  xpos = 0. * cm;
197  ypos = 0. * cm;
198  G4double len_bbch = z_bbch[1] - z_bbch[0];
199  zpos = z_bbcd[0] + 0.3 * cm + len_bbch * 0.5;
200 
201  G4VPhysicalVolume *bbch_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbch_lv, "BBCH", bbcd_lv, false, 0);
202  if (!bbch_phys) // more to prevent compiler warnings
203  {
204  std::cout << "placement of BBCH failed" << std::endl;
205  }
206 
207  // Locations of the 64 PMT tubes in an arm.
208  // These are the x,y for the south BBC.
209  // The north inverts the x coordinate (x -> -x)
210  // (NB: Should probably move this to a geometry object...)
211  float TubeLoc[64][2] = {
212  {-12.2976, 4.26},
213  {-12.2976, 1.42},
214  {-9.83805, 8.52},
215  {-9.83805, 5.68},
216  {-9.83805, 2.84},
217  {-7.37854, 9.94},
218  {-7.37854, 7.1},
219  {-7.37854, 4.26},
220  {-7.37854, 1.42},
221  {-4.91902, 11.36},
222  {-4.91902, 8.52},
223  {-4.91902, 5.68},
224  {-2.45951, 12.78},
225  {-2.45951, 9.94},
226  {-2.45951, 7.1},
227  {0, 11.36},
228  {0, 8.52},
229  {2.45951, 12.78},
230  {2.45951, 9.94},
231  {2.45951, 7.1},
232  {4.91902, 11.36},
233  {4.91902, 8.52},
234  {4.91902, 5.68},
235  {7.37854, 9.94},
236  {7.37854, 7.1},
237  {7.37854, 4.26},
238  {7.37854, 1.42},
239  {9.83805, 8.52},
240  {9.83805, 5.68},
241  {9.83805, 2.84},
242  {12.2976, 4.26},
243  {12.2976, 1.42},
244  {12.2976, -4.26},
245  {12.2976, -1.42},
246  {9.83805, -8.52},
247  {9.83805, -5.68},
248  {9.83805, -2.84},
249  {7.37854, -9.94},
250  {7.37854, -7.1},
251  {7.37854, -4.26},
252  {7.37854, -1.42},
253  {4.91902, -11.36},
254  {4.91902, -8.52},
255  {4.91902, -5.68},
256  {2.45951, -12.78},
257  {2.45951, -9.94},
258  {2.45951, -7.1},
259  {0, -11.36},
260  {0, -8.52},
261  {-2.45951, -12.78},
262  {-2.45951, -9.94},
263  {-2.45951, -7.1},
264  {-4.91902, -11.36},
265  {-4.91902, -8.52},
266  {-4.91902, -5.68},
267  {-7.37854, -9.94},
268  {-7.37854, -7.1},
269  {-7.37854, -4.26},
270  {-7.37854, -1.42},
271  {-9.83805, -8.52},
272  {-9.83805, -5.68},
273  {-9.83805, -2.84},
274  {-12.2976, -4.26},
275  {-12.2976, -1.42}};
276 
277  //m_bbcz = m_Params->get_double_param("z") * cm;
278  m_bbcz = 250.0 * cm; // The front face of the quartz is at 250 cm
279 
280  const float tube_zpos = m_bbcz + len_bbcd / 2.0 - len_bbca;
281 
282  const int NPMT = 64; // No. PMTs per arm
283  G4RotationMatrix *arm_rot[2];
284  for (int iarm = 0; iarm < 2; iarm++)
285  {
286  arm_rot[iarm] = new G4RotationMatrix;
287  float xside = 1.0;
288  float zside = 1.0;
289  if (iarm == 0) // South Arm = 0
290  {
291  xside = 1.0;
292  zside = -1.0;
293  arm_rot[iarm]->rotateY(180 * deg);
294  }
295  else // North Arm = 1
296  {
297  xside = -1.0;
298  zside = 1.0;
299  }
300 
301  // Add BBC PMT's
302  for (int itube = 0; itube < NPMT; itube++)
303  {
304  // Full PMT Housing with Active Quartz Cerenkov Radiators
305  float tube_xpos = xside * TubeLoc[itube][0] * cm;
306  float tube_ypos = TubeLoc[itube][1] * cm;
307  new G4PVPlacement(arm_rot[iarm], G4ThreeVector(tube_xpos, tube_ypos, zside * tube_zpos),
308  bbcd_lv, "BBCD", logicWorld, false, iarm * NPMT + itube, OverlapCheck());
309  }
310  }
311 
312  // Now Build the BBC Housing
313 
314  // BBC Outer and Inner Cylindrical Shells
315  G4Tubs *bbc_outer_shell = new G4Tubs("bbc_outer_shell", 14.9 * cm, 15 * cm, 11.5 * cm, 0, 2 * M_PI);
316  G4LogicalVolume *bbc_outer_shell_lv = new G4LogicalVolume(bbc_outer_shell, Aluminum, G4String("Bbc_Outer_Shell"));
317  GetDisplayAction()->AddVolume(bbc_outer_shell_lv, "Bbc_Outer_Shell");
318 
319  G4Tubs *bbc_inner_shell = new G4Tubs("bbc_inner_shell", 5.0 * cm, 5.5 * cm, 11.5 * cm, 0, 2 * M_PI);
320  G4LogicalVolume *bbc_inner_shell_lv = new G4LogicalVolume(bbc_inner_shell, Aluminum, G4String("Bbc_Inner_Shell"));
321  GetDisplayAction()->AddVolume(bbc_inner_shell_lv, "Bbc_Inner_Shell");
322 
323  G4VPhysicalVolume *outer_shell_vol[2] = {nullptr};
324  G4VPhysicalVolume *inner_shell_vol[2] = {nullptr};
325 
326  // Place South Shells
327  outer_shell_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 - 11.5) * cm),
328  bbc_outer_shell_lv, "BBC_OUTER_SHELL", logicWorld, false, 0, OverlapCheck());
329  inner_shell_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 - 11.5) * cm),
330  bbc_inner_shell_lv, "BBC_INNER_SHELL", logicWorld, false, 0, OverlapCheck());
331 
332  // Place North Shells
333  outer_shell_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 + 11.5) * cm),
334  bbc_outer_shell_lv, "BBC_OUTER_SHELL", logicWorld, false, 1, OverlapCheck());
335  inner_shell_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 + 11.5) * cm),
336  bbc_inner_shell_lv, "BBC_INNER_SHELL", logicWorld, false, 0, OverlapCheck());
337  // this is more to prevent compiler warnings about unused variables
338  if (!outer_shell_vol[0] || !outer_shell_vol[1] || !inner_shell_vol[0] || !inner_shell_vol[1])
339  {
340  std::cout << "problem placing BBC Sheels" << std::endl;
341  }
342 
343  // BBC Front and Back Plates
344  G4Tubs *bbc_plate = new G4Tubs("bbc_fplate", 5 * cm, 15 * cm, 0.5 * cm, 0, 2 * M_PI);
345  G4LogicalVolume *bbc_plate_lv = new G4LogicalVolume(bbc_plate, Aluminum, G4String("Bbc_Cover_Plates"));
346  GetDisplayAction()->AddVolume(bbc_plate_lv, "Bbc_Cover_Plates");
347 
348  G4VPhysicalVolume *fplate_vol[2] = {nullptr}; // Front Plates
349  G4VPhysicalVolume *bplate_vol[2] = {nullptr}; // Back Plates
350 
351  // Place South Plates
352  fplate_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 + 0.5) * cm),
353  bbc_plate_lv, "BBC_FPLATE", logicWorld, false, 0, OverlapCheck());
354  bplate_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 + 0.5 - 24.0) * cm),
355  bbc_plate_lv, "BBC_BPLATE", logicWorld, false, 0, OverlapCheck());
356 
357  // Place North Plates
358  fplate_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 - 0.5) * cm),
359  bbc_plate_lv, "BBC_FPLATE", logicWorld, false, 1, OverlapCheck());
360  bplate_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 - 0.5 + 24.0) * cm),
361  bbc_plate_lv, "BBC_BPLATE", logicWorld, false, 0, OverlapCheck());
362 
363  // Place BBC Cables
364  G4Material *Cu = manager->FindOrBuildMaterial("G4_Cu");
365  const G4double len_cable = 120 * cm;
366  const G4double r_CableConductor = 0.09525 * cm;
367  G4Tubs *bbc_cablecond = new G4Tubs("bbc_cablecond", 0., r_CableConductor, len_cable * 0.5, 0 * deg, 360 * deg);
368 
369  G4LogicalVolume *bbc_cablecond_lv = new G4LogicalVolume(bbc_cablecond, Cu, G4String("Bbc_CableCond"));
370  GetDisplayAction()->AddVolume(bbc_cablecond_lv, "Bbc_CableCond");
371 
372  const G4double rIn_CableShield = 0.302876 * cm;
373  const G4double rOut_CableShield = 0.3175 * cm;
374  G4Tubs *bbc_cableshield = new G4Tubs("bbc_cableshield", rIn_CableShield, rOut_CableShield, len_cable * 0.5, 0 * deg, 360 * deg);
375 
376  G4LogicalVolume *bbc_cableshield_lv = new G4LogicalVolume(bbc_cableshield, Cu, G4String("Bbc_CableShield"));
377  GetDisplayAction()->AddVolume(bbc_cableshield_lv, "Bbc_CableShield");
378 
379  ypos = len_cable / 2 + 5 * cm;
380 
381  // For now we make this vertical, but they should be slanted toward the endcap
382  G4RotationMatrix *rot_cable = new G4RotationMatrix();
383  rot_cable->rotateX(90 * deg);
384 
385  int icable = 0;
386 
387  for (int iarm = 0; iarm < 2; iarm++)
388  {
389  float zsign = -1.0;
390  if (iarm == 1)
391  {
392  zsign = 1.0;
393  }
394 
395  for (int iring = 1; iring < 5; iring++)
396  {
397  float ring_radius = iring * 0.67 * cm;
398  int ncables = 2 * M_PI * ring_radius / (0.635 * cm);
399  double dphi = 2 * M_PI / ncables;
400 
401  //G4cout << "BBC_CABLE " << iring << "\t" << ring_radius << "\t" << ncables << "\t" << dphi*180/3.14 << G4endl;
402 
403  // place cables in ring
404  for (int ic = 0; ic < ncables; ic++)
405  {
406  xpos = cos(dphi * ic) * ring_radius;
407  zpos = sin(dphi * ic) * ring_radius + zsign * (m_bbcz + 30 * cm);
408  // Place Inner Conductor
409  new G4PVPlacement(rot_cable, G4ThreeVector(xpos, ypos, zpos), bbc_cablecond_lv, "BBC_Cable_Cond", logicWorld, false, icable, OverlapCheck());
410  // Place Shield
411  new G4PVPlacement(rot_cable, G4ThreeVector(xpos, ypos, zpos), bbc_cableshield_lv, "BBC_Cable_Shield", logicWorld, false, icable++, OverlapCheck());
412  }
413  }
414  }
415 
416  // Now make the Supports
417  ConstructSupport(logicWorld);
418 
419  // this is more to prevent compiler warnings about unused variables
420  if (!fplate_vol[0] || !fplate_vol[1] || !bplate_vol[0] || !bplate_vol[1])
421  {
422  std::cout << "problem placing BBC Sheets" << std::endl;
423  }
424 
425  // bbcq is the active detector element
426  m_PhysLogicalVolSet.insert(bbcq_lv);
427 
428  // GetDisplayAction()->AddVolume(bbc_plate_lv, "Bbc_plate");
429  // GetDisplayAction()->AddVolume(bbc_outer_shell_lv, "Bbc_oshell");
430  // GetDisplayAction()->AddVolume(bbc_inner_shell_lv, "Bbc_ishell");
431 
432  std::cout << "INSIDE BBC" << std::endl;
433 
434  return;
435 }
436 
437 void PHG4BbcDetector::ConstructSupport(G4LogicalVolume *logicWorld)
438 {
439  //std::cout << "PHG4BbcDetector::ConstructSupport()" << std::endl;
440 
441  G4double fractionmass;
442  G4int ncomponents;
443  G4double density;
444 
446  G4Material *WorldMaterial = GetDetectorMaterial(rc->get_StringFlag("WorldMaterial"));
447  G4Material *Delrin = GetDetectorMaterial("G4_POLYOXYMETHYLENE");
448  //G4Material *Aluminum = GetDetectorMaterial("G4_Al");
449 
450  // McMaster-Carr #8548K36, implemented as PDG G10, https://pdg.lbl.gov/2022/AtomicNuclearProperties/HTML/G10.html
451  G4Material *Fiberglass = new G4Material("BBC_Fiberglass", density = 1.800 * g / cm3, ncomponents = 9);
452  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("B"), fractionmass = 0.018640);
453  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Mg"), fractionmass = 0.010842);
454  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("O"), fractionmass = 0.044453);
455  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Al"), fractionmass = 0.151423);
456  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("H"), fractionmass = 0.081496);
457  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Ca"), fractionmass = 0.385764);
458  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("C"), fractionmass = 0.275853);
459  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("N"), fractionmass = 0.003583);
460  Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Si"), fractionmass = 0.027945);
461 
462  // BBC Base Plates
463  G4double basep_width = 35.56 * cm;
464  G4double basep_height = 1.91 * cm;
465  G4double basep_len = 46.99 * cm;
466  G4double basep_zpos = 228.6 * cm; // z-pos of front edge of base plate
467  G4Box *bbc_base_plate = new G4Box("bbc_base_plate", basep_width/2, basep_height/2, basep_len/2);
468  G4LogicalVolume *bbc_base_plate_lv = new G4LogicalVolume(bbc_base_plate, Delrin, G4String("Bbc_Base_Plates"));
469  GetDisplayAction()->AddVolume(bbc_base_plate_lv, "Bbc_Base_Plates");
470 
471  // Place South Base Plates
472  new G4PVPlacement(nullptr, G4ThreeVector( 0 *cm, -15.*cm - basep_height/2, (-basep_zpos - basep_len/2 ) ),
473  bbc_base_plate_lv, "BBC_BASE_PLATE", logicWorld, false, 0, OverlapCheck());
474 
475  // Place North Base Plates
476  new G4PVPlacement(nullptr, G4ThreeVector( 0 *cm, -15.*cm - basep_height/2, (basep_zpos + basep_len/2 ) ),
477  bbc_base_plate_lv, "BBC_BASE_PLATE", logicWorld, false, 1, OverlapCheck());
478 
479 
480  // BBC Side Support Plates
481  G4double sidesupportp_width = 1.27 * cm;
482  G4double sidesupportp_height = 14.57 * cm;
483  G4double sidesupportp_len = 25.00 * cm;
484 
485  G4Box *bbc_sidesupport_plate = new G4Box("bbc_sidesupport_plate", sidesupportp_width/2, sidesupportp_height/2, sidesupportp_len/2);
486  G4LogicalVolume *bbc_sidesupport_plate_lv = new G4LogicalVolume(bbc_sidesupport_plate, Delrin, G4String("Bbc_Sidesupport_Plates"));
487  GetDisplayAction()->AddVolume(bbc_sidesupport_plate_lv, "Bbc_Sidesupport_Plates");
488 
489  // Make and Place Holes in Side Support Plates
490  G4Tubs *bbc_sidesupport_hole = new G4Tubs("bbc_sidesupport_hole", 0., (7.62/2)*cm, sidesupportp_width/2, 0, 2.0*M_PI);
491  G4LogicalVolume *bbc_sidesupport_hole_lv = new G4LogicalVolume(bbc_sidesupport_hole, WorldMaterial, G4String("Bbc_Sidesupport_Holes"));
492 
493  G4RotationMatrix *rot_sideholes = new G4RotationMatrix;
494  rot_sideholes->rotateY(90.*deg);
495  new G4PVPlacement(rot_sideholes, G4ThreeVector(0,-0.935*cm,11.43*cm/2), bbc_sidesupport_hole_lv, "BBC_SIDESUPPORT_HOLE", bbc_sidesupport_plate_lv, false, 0, OverlapCheck());
496  new G4PVPlacement(rot_sideholes, G4ThreeVector(0,-0.935*cm,-11.43*cm/2), bbc_sidesupport_hole_lv, "BBC_SIDESUPPORT_HOLE", bbc_sidesupport_plate_lv, false, 1, OverlapCheck());
497 
498  // Place South Side Support Plates
499  new G4PVPlacement(nullptr, G4ThreeVector(-15*cm-sidesupportp_width/2, -15*cm+sidesupportp_height/2, -front_bbcz-sidesupportp_len/2),
500  bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 0, OverlapCheck());
501  new G4PVPlacement(nullptr, G4ThreeVector(15*cm+sidesupportp_width/2, -15*cm+sidesupportp_height/2, -front_bbcz-sidesupportp_len/2),
502  bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 1, OverlapCheck());
503 
504  // Place North Side Support Plate
505  new G4PVPlacement(nullptr, G4ThreeVector(-15*cm-sidesupportp_width/2, -15*cm+sidesupportp_height/2, front_bbcz+sidesupportp_len/2),
506  bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 2, OverlapCheck());
507  new G4PVPlacement(nullptr, G4ThreeVector(15*cm+sidesupportp_width/2, -15*cm+sidesupportp_height/2, front_bbcz+sidesupportp_len/2),
508  bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 3, OverlapCheck());
509 
510  // BBC Support Post
511  G4double supportp_width = 10.16 * cm;
512  G4double supportp_height = 97.16 * cm;
513  G4double supportp_len = 10.16 * cm;
514  G4double supportp_thick = 0.64 * cm;
515 
516  G4Box *bbc_support_post_outside = new G4Box("bbc_support_post_outside", supportp_width/2, supportp_height/2, supportp_len/2);
517  G4Box *bbc_support_post_inside = new G4Box("bbc_support_post_inside", supportp_width/2-supportp_thick, supportp_height/2-supportp_thick, supportp_len/2-supportp_thick);
518  G4SubtractionSolid* bbc_support_post = new G4SubtractionSolid("support_post", bbc_support_post_outside, bbc_support_post_inside);
519 
520  G4LogicalVolume *bbc_support_post_lv = new G4LogicalVolume(bbc_support_post, Fiberglass, G4String("Bbc_Support_Post"));
521 
522  GetDisplayAction()->AddVolume(bbc_support_post_lv, "Bbc_Support_Post");
523 
524  // Place South Support Post
525  new G4PVPlacement(nullptr, G4ThreeVector(0, -15.*cm-basep_height-supportp_height/2, -228.*cm - supportp_len/2),
526  bbc_support_post_lv, "BBC_SUPPORT_POST", logicWorld, false, 0, OverlapCheck());
527 
528  // Place North Support Post
529  new G4PVPlacement(nullptr, G4ThreeVector(0, -15.*cm-basep_height-supportp_height/2, 228.*cm + supportp_len/2),
530  bbc_support_post_lv, "BBC_SUPPORT_POST", logicWorld, false, 1, OverlapCheck());
531 
532 
533  // BBC Support Arm
534  G4double supporta_width = 148.59 * cm;
535  G4double supporta_height = 10.16 * cm;
536  G4double supporta_len = 10.16 * cm;
537  G4double supporta_thick = 0.64 * cm;
538 
539  G4Box *bbc_support_arm_outside = new G4Box("bbc_support_arm_outside", supporta_width/2, supporta_height/2, supporta_len/2);
540  G4Box *bbc_support_arm_inside = new G4Box("bbc_support_arm_inside", supporta_width/2-supporta_thick, supporta_height/2-supporta_thick, supporta_len/2-supporta_thick);
541  G4SubtractionSolid* bbc_support_arm = new G4SubtractionSolid("support_arm", bbc_support_arm_outside, bbc_support_arm_inside);
542 
543  G4LogicalVolume *bbc_support_arm_lv = new G4LogicalVolume(bbc_support_arm, Fiberglass, G4String("Bbc_Support_Arm"));
544 
545  GetDisplayAction()->AddVolume(bbc_support_arm_lv, "Bbc_Support_Arm");
546 
547  // Place South Support Arms
548  new G4PVPlacement(nullptr, G4ThreeVector( -supporta_width/2 - supportp_width/2, -15.*cm-basep_height-20.48*cm - supporta_height/2, -228.6*cm - supporta_len/2),
549  bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 0, OverlapCheck());
550  new G4PVPlacement(nullptr, G4ThreeVector( supporta_width/2 + supportp_width/2, -15.*cm-basep_height-20.48*cm - supporta_height/2, -228.6*cm - supporta_len/2),
551  bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 1, OverlapCheck());
552 
553  // Place North Support Arms
554  new G4PVPlacement(nullptr, G4ThreeVector( -supporta_width/2 - supportp_width/2, -15.*cm-basep_height-20.48*cm - supporta_height/2, 228.6*cm + supporta_len/2),
555  bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 2, OverlapCheck());
556  new G4PVPlacement(nullptr, G4ThreeVector( supporta_width/2 + supportp_width/2, -15.*cm-basep_height-20.48*cm - supporta_height/2, 228.6*cm + supporta_len/2),
557  bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 3, OverlapCheck());
558 
559 
560  // BBC Gusset Plates (implementation is broken up into 3 parts, in Trap and 2 Boxes)
561  G4double gussetp0_pz = 1.27 * cm;
562  G4double gussetp0_py = (45.72 - 11.11) * cm;
563  G4double gussetp0_px = 44.62 * cm; // measured off drawing
564  G4double gussetp0_pltx = 5.08 * cm;
565  G4Trap *bbc_gusset_plate0 = new G4Trap("bbc_gusset_plate0", gussetp0_pz, gussetp0_py, gussetp0_px, gussetp0_pltx);
566  G4LogicalVolume *bbc_gusset0_plate_lv = new G4LogicalVolume(bbc_gusset_plate0, Delrin, G4String("Bbc_Gusset0_Plates"));
567  GetDisplayAction()->AddVolume(bbc_gusset0_plate_lv, "Bbc_Gusset0_Plates");
568 
569  // Place South Gusset Plates (Trapezoid part)
570  G4RotationMatrix *rot_sgusset = new G4RotationMatrix;
571  rot_sgusset->rotateY(90.*deg);
572  rot_sgusset->rotateZ(180.*deg);
573 
574  G4double xpos = supportp_width/2 + gussetp0_pz/2;
575  G4double ypos = -gussetp0_py/2 - 15*cm - basep_height;
576  G4double zpos = 0.25*(gussetp0_px + gussetp0_pltx) + 228.6*cm + 11.11*cm;
577  new G4PVPlacement(rot_sgusset, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 0, OverlapCheck());
578  new G4PVPlacement(rot_sgusset, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 1, OverlapCheck());
579 
580  // Place North Gusset Plates (Trapezoid part)
581  G4RotationMatrix *rot_ngusset = new G4RotationMatrix;
582  rot_ngusset->rotateY(-90.*deg);
583  rot_ngusset->rotateZ(180.*deg);
584 
585  new G4PVPlacement(rot_ngusset, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 2, OverlapCheck());
586  new G4PVPlacement(rot_ngusset, G4ThreeVector(xpos, ypos, zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 3, OverlapCheck());
587 
588  // Gusset Plate Box 1
589  G4double gussetp1_x = 1.27 * cm; // measured off drawing
590  G4double gussetp1_y = 20.48 * cm;
591  G4double gussetp1_z = 11.11 * cm;
592 
593  G4Box *bbc_gusset_plate1 = new G4Box("bbc_gusset_plate1", gussetp1_x/2, gussetp1_y/2, gussetp1_z/2);
594  G4LogicalVolume *bbc_gusset1_plate_lv = new G4LogicalVolume(bbc_gusset_plate1, Delrin, G4String("Bbc_Gusset1_Plates"));
595  GetDisplayAction()->AddVolume(bbc_gusset1_plate_lv, "Bbc_Gusset1_Plates");
596 
597  // Place South Gusset Plates (Box 1)
598  ypos = -15*cm - basep_height - gussetp1_y/2;
599  zpos = 228.6*cm + gussetp1_z/2;
600  new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 0, OverlapCheck());
601  new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 1, OverlapCheck());
602 
603  // Place North Gusset Plates (Box 1)
604  new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 2, OverlapCheck());
605  new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 3, OverlapCheck());
606 
607  // Gusset Plate Box 2
608  G4double gussetp2_x = 1.27 * cm; // measured off drawing
609  G4double gussetp2_y = (45.72 - 10.80 - 20.48) * cm;
610  G4double gussetp2_z = 11.11 * cm;
611 
612  G4Box *bbc_gusset_plate2 = new G4Box("bbc_gusset_plate2", gussetp2_x/2, gussetp2_y/2, gussetp2_z/2);
613  G4LogicalVolume *bbc_gusset2_plate_lv = new G4LogicalVolume(bbc_gusset_plate2, Delrin, G4String("Bbc_Gusset2_Plates"));
614  GetDisplayAction()->AddVolume(bbc_gusset2_plate_lv, "Bbc_Gusset2_Plates");
615 
616  // Place South Gusset Plates (Box 2)
617  ypos = -15*cm - basep_height - 20.48*cm - 10.80*cm - gussetp2_y/2;
618  zpos = 228.6*cm + gussetp2_z/2;
619  new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 0, OverlapCheck());
620  new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 1, OverlapCheck());
621 
622  // Place North Gusset Plates (Box 2)
623  new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 2, OverlapCheck());
624  new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 3, OverlapCheck());
625 
626  // BBC Splice Plates
627  G4double splicep_x = 35.56 * cm;
628  G4double splicep_y = 10.16 * cm;
629  G4double splicep_z = 0.64 * cm;
630 
631  G4Box *bbc_splice_plate = new G4Box("bbc_splice_plate", splicep_x/2, splicep_y/2, splicep_z/2);
632  G4LogicalVolume *bbc_splice_plate_lv = new G4LogicalVolume(bbc_splice_plate, Delrin, G4String("Bbc_Splice_Plates"));
633  GetDisplayAction()->AddVolume(bbc_splice_plate_lv, "Bbc_Splice_Plates");
634 
635  // Make and Place Holes in Splice Plates
636  G4Tubs *bbc_splice_hole = new G4Tubs("bbc_splice_hole", 0., (6.35/2)*cm, splicep_z/2, 0, 2.0*M_PI);
637  G4LogicalVolume *bbc_splice_hole_lv = new G4LogicalVolume(bbc_splice_hole, WorldMaterial, G4String("Bbc_Splice_Holes"));
638 
639  new G4PVPlacement(nullptr, G4ThreeVector(-12.7*cm,0,0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 0, OverlapCheck());
640  new G4PVPlacement(nullptr, G4ThreeVector(0,0,0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 1, OverlapCheck());
641  new G4PVPlacement(nullptr, G4ThreeVector(12.7*cm,0,0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 2, OverlapCheck());
642 
643  // Place South Splice Plate
644  ypos = -15*cm - basep_height - 20.48*cm - splicep_y/2;
645  zpos = 228.6*cm + supportp_len + splicep_z/2;
646  new G4PVPlacement(nullptr, G4ThreeVector(0,ypos,-zpos), bbc_splice_plate_lv, "BBC_SPLICE_PLATE", logicWorld, false, 0, OverlapCheck());
647 
648  // Place North Splice Plate
649  new G4PVPlacement(nullptr, G4ThreeVector(0,ypos,zpos), bbc_splice_plate_lv, "BBC_SPLICE_PLATE", logicWorld, false, 1, OverlapCheck());
650 
651 
652 }
653 
654 void PHG4BbcDetector::Print(const std::string &what) const
655 {
656  std::cout << "Bbc Detector:" << std::endl;
657  if (what == "ALL" || what == "VOLUME")
658  {
659  std::cout << "Version 0.1" << std::endl;
660  }
661  return;
662 }