Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:18:53

0001 #include "PHG4BbcDetector.h"
0002 
0003 #include "PHG4BbcDisplayAction.h"
0004 
0005 #include <phparameter/PHParameters.h>
0006 
0007 #include <g4main/PHG4Detector.h>
0008 #include <g4main/PHG4DisplayAction.h>
0009 #include <g4main/PHG4Subsystem.h>
0010 
0011 #include <phool/recoConsts.h>
0012 
0013 #include <Geant4/G4Box.hh>
0014 #include <Geant4/G4LogicalVolume.hh>
0015 #include <Geant4/G4Material.hh>
0016 #include <Geant4/G4NistManager.hh>
0017 #include <Geant4/G4PVPlacement.hh>
0018 #include <Geant4/G4Polyhedra.hh>
0019 #include <Geant4/G4RotationMatrix.hh>  // for G4RotationMatrix
0020 #include <Geant4/G4String.hh>          // for G4String
0021 #include <Geant4/G4SubtractionSolid.hh>
0022 #include <Geant4/G4SystemOfUnits.hh>
0023 #include <Geant4/G4ThreeVector.hh>
0024 #include <Geant4/G4Trap.hh>
0025 #include <Geant4/G4Tubs.hh>
0026 
0027 #include <Geant4/G4Types.hh>  // for G4double, G4int
0028 #include <Geant4/G4VPhysicalVolume.hh>
0029 
0030 #include <cmath>
0031 #include <iostream>  // for operator<<, endl, bas...
0032 #include <vector>
0033 
0034 class PHCompositeNode;
0035 
0036 PHG4BbcDetector::PHG4BbcDetector(PHG4Subsystem *subsys, PHCompositeNode *Node, PHParameters *params, const std::string &dnam)
0037   : PHG4Detector(subsys, Node, dnam)
0038   , m_DisplayAction(dynamic_cast<PHG4BbcDisplayAction *>(subsys->GetDisplayAction()))
0039   , m_Params(params)
0040   , m_ActiveFlag(m_Params->get_int_param("active"))
0041   , m_SupportActiveFlag(m_Params->get_int_param("supportactive"))
0042   , front_bbcz(248 * cm)
0043 {
0044 }
0045 
0046 //_______________________________________________________________
0047 //_______________________________________________________________
0048 int PHG4BbcDetector::IsInBbc(G4VPhysicalVolume *volume) const
0049 {
0050   G4LogicalVolume *mylogvol = volume->GetLogicalVolume();
0051 
0052   if (m_ActiveFlag)
0053   {
0054     if (m_PhysLogicalVolSet.find(mylogvol) != m_PhysLogicalVolSet.end())
0055     {
0056       return 1;
0057     }
0058   }
0059   if (m_SupportActiveFlag)
0060   {
0061     if (m_SupportLogicalVolSet.find(mylogvol) != m_SupportLogicalVolSet.end())
0062     {
0063       return -2;
0064     }
0065   }
0066   return 0;
0067 }
0068 
0069 void PHG4BbcDetector::ConstructMe(G4LogicalVolume *logicWorld)
0070 {
0071   // std::cout << "PHG4BbcDetector::ConstructMe()" << std::endl;
0072 
0073   recoConsts *rc = recoConsts::instance();
0074 
0075   //========================
0076   // Build a single BBC PMT
0077   //========================
0078 
0079   // BBC Detector Tube (BBCD) Logical Volume (contains all parts of PMT+Radiator)
0080   G4Material *WorldMaterial = GetDetectorMaterial(rc->get_StringFlag("WorldMaterial"));
0081   G4double z_bbcd[] = {-6.15 * cm, 6.15 * cm};
0082   G4double len_bbcd = z_bbcd[1] - z_bbcd[0];
0083   G4double rin_bbcd[] = {0 * cm, 0 * cm};
0084   G4double rout_bbcd[] = {1.4 * cm, 1.4 * cm};
0085   G4Polyhedra *bbcd = new G4Polyhedra("bbcd", 0., 2 * M_PI, 6, 2, z_bbcd, rin_bbcd, rout_bbcd);
0086   G4LogicalVolume *bbcd_lv = new G4LogicalVolume(bbcd, WorldMaterial, G4String("Bbc_tube"));
0087 
0088   //
0089   //  Place the BBCA, BBCQ, BBCP, BBCR and BBCH in to BBCD.
0090   //
0091 
0092   // BBC Attachment Plate (BBCA)
0093   const G4double z_bbca[] = {-0.5 * cm, -0.201 * cm, -0.2 * cm, 0.5 * cm};
0094   const G4double rInner_bbca[] = {0.2 * cm, 0.2 * cm, 0.2 * cm, 0.2 * cm};
0095   const G4double rOuter_bbca[] = {1.4 * cm, 1.4 * cm, 1.375 * cm, 1.375 * cm};
0096   G4Polyhedra *bbca = new G4Polyhedra("bbca", 0., 2 * M_PI, 6, 4, z_bbca, rInner_bbca, rOuter_bbca);
0097 
0098   G4Material *Aluminum = GetDetectorMaterial("G4_Al");
0099 
0100   G4LogicalVolume *bbca_lv = new G4LogicalVolume(bbca, Aluminum, G4String("Bbc_attach_plate"));
0101   m_SupportLogicalVolSet.insert(bbca_lv);
0102   GetDisplayAction()->AddVolume(bbca_lv, "Bbc_attach_plate");
0103 
0104   G4double xpos = 0. * cm;
0105   G4double ypos = 0. * cm;
0106   G4double len_bbca = z_bbca[3] - z_bbca[0];
0107   G4double zpos = z_bbcd[0] + len_bbca * 0.5;
0108 
0109   G4VPhysicalVolume *bbca_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbca_lv, "BBCA", bbcd_lv, false, 0);
0110   if (!bbca_phys)  // more to prevent compiler warnings
0111   {
0112     std::cout << "placement of BBCA failed" << std::endl;
0113   }
0114 
0115   G4Material *Quartz = GetDetectorMaterial("G4_SILICON_DIOXIDE");
0116 
0117   // BBC Quartz Radiator
0118   const G4double z_bbcq[] = {-1.5 * cm, 1.5 * cm};
0119   const G4double rInner_bbcq[] = {0., 0.};
0120   const G4double rOuter_bbcq[] = {1.27 * cm, 1.27 * cm};
0121   G4Polyhedra *bbcq = new G4Polyhedra("bbcq", 0., 2 * M_PI, 6, 2, z_bbcq, rInner_bbcq, rOuter_bbcq);
0122 
0123   G4LogicalVolume *bbcq_lv = new G4LogicalVolume(bbcq, Quartz, G4String("Bbc_quartz"));
0124   GetDisplayAction()->AddVolume(bbcq_lv, "Bbc_quartz");
0125 
0126   xpos = 0. * cm;
0127   ypos = 0. * cm;
0128   G4double len_bbcq = z_bbcq[1] - z_bbcq[0];
0129   zpos += len_bbca * 0.5 + len_bbcq * 0.5;
0130 
0131   G4VPhysicalVolume *bbcq_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcq_lv, "BBCQ", bbcd_lv, false, 0);
0132   if (!bbcq_phys)  // more to prevent compiler warnings
0133   {
0134     std::cout << "placement of BBCQ failed" << std::endl;
0135   }
0136 
0137   // BBC PMT
0138   const G4double len_bbcp = 4.4 * cm;
0139   const G4double rInner_bbcp = 1.09 * cm;
0140   const G4double rOuter_bbcp = 1.29 * cm;
0141   G4Tubs *bbcp = new G4Tubs("bbcp", rInner_bbcp, rOuter_bbcp, len_bbcp * 0.5, 0 * deg, 360 * deg);
0142 
0143   G4LogicalVolume *bbcp_lv = new G4LogicalVolume(bbcp, Quartz, G4String("Bbc_PMT"));
0144   GetDisplayAction()->AddVolume(bbcp_lv, "Bbc_PMT");
0145 
0146   xpos = 0. * cm;
0147   ypos = 0. * cm;
0148   zpos += len_bbcq * 0.5 + len_bbcp * 0.5;
0149 
0150   G4VPhysicalVolume *bbcp_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcp_lv, "BBCP", bbcd_lv, false, 0);
0151   if (!bbcp_phys)  // more to prevent compiler warnings
0152   {
0153     std::cout << "placement of BBCP failed" << std::endl;
0154   }
0155 
0156   // BBC Breeder Module
0157   const G4double len_bbcr = 3.9 * cm;
0158   const G4double rInner_bbcr = 0.0 * cm;
0159   const G4double rOuter_bbcr = 1.2 * cm;
0160   G4Tubs *bbcr = new G4Tubs("bbcr", rInner_bbcr, rOuter_bbcr, len_bbcr * 0.5, 0 * deg, 360 * deg);
0161 
0162   G4int natoms;
0163   G4int ncomponents;
0164   G4double density;
0165   G4Material *G10 = new G4Material("BBC_G10", density = 1.700 * g / cm3, ncomponents = 4);
0166   G4NistManager *manager = G4NistManager::Instance();
0167   G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("Si"), natoms = 1);
0168   G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("O"), natoms = 2);
0169   G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("C"), natoms = 3);
0170   G10->AddElement(G4NistManager::Instance()->FindOrBuildElement("H"), natoms = 3);
0171 
0172   G4LogicalVolume *bbcr_lv = new G4LogicalVolume(bbcr, G10, "Bbc_Breeder_Module");
0173   GetDisplayAction()->AddVolume(bbcr_lv, "Bbc_Breeder_Module");
0174 
0175   xpos = 0. * cm;
0176   ypos = 0. * cm;
0177   zpos += len_bbcp * 0.5 + len_bbcr * 0.5;
0178 
0179   G4PVPlacement *plcmt = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbcr_lv, "BBCR", bbcd_lv, false, 0);
0180   if (!plcmt)  // more to prevent compiler warnings
0181   {
0182     std::cout << "placement of BBCR failed" << std::endl;
0183   }
0184 
0185   // BBC Mu Metal Shield
0186   const G4double z_bbch[] = {-5.65 * cm, 5.65 * cm};
0187   const G4double rInner_bbch[] = {1.375 * cm, 1.375 * cm};
0188   const G4double rOuter_bbch[] = {1.4 * cm, 1.4 * cm};
0189   G4Polyhedra *bbch = new G4Polyhedra("bbch", 0., 2 * M_PI, 6, 2, z_bbch, rInner_bbch, rOuter_bbch);
0190 
0191   G4Material *MuMetal = manager->FindOrBuildMaterial("G4_STAINLESS-STEEL");
0192 
0193   G4LogicalVolume *bbch_lv = new G4LogicalVolume(bbch, MuMetal, G4String("Bbc_Shield"));
0194 
0195   xpos = 0. * cm;
0196   ypos = 0. * cm;
0197   G4double len_bbch = z_bbch[1] - z_bbch[0];
0198   zpos = z_bbcd[0] + 0.3 * cm + len_bbch * 0.5;
0199 
0200   G4VPhysicalVolume *bbch_phys = new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbch_lv, "BBCH", bbcd_lv, false, 0);
0201   if (!bbch_phys)  // more to prevent compiler warnings
0202   {
0203     std::cout << "placement of BBCH failed" << std::endl;
0204   }
0205 
0206   // Locations of the 64 PMT tubes in an arm.
0207   // These are the x,y for the south BBC.
0208   // The north inverts the x coordinate (x -> -x)
0209   // (NB: Should probably move this to a geometry object...)
0210   float TubeLoc[64][2] = {
0211       {-12.2976, 4.26},
0212       {-12.2976, 1.42},
0213       {-9.83805, 8.52},
0214       {-9.83805, 5.68},
0215       {-9.83805, 2.84},
0216       {-7.37854, 9.94},
0217       {-7.37854, 7.1},
0218       {-7.37854, 4.26},
0219       {-7.37854, 1.42},
0220       {-4.91902, 11.36},
0221       {-4.91902, 8.52},
0222       {-4.91902, 5.68},
0223       {-2.45951, 12.78},
0224       {-2.45951, 9.94},
0225       {-2.45951, 7.1},
0226       {0, 11.36},
0227       {0, 8.52},
0228       {2.45951, 12.78},
0229       {2.45951, 9.94},
0230       {2.45951, 7.1},
0231       {4.91902, 11.36},
0232       {4.91902, 8.52},
0233       {4.91902, 5.68},
0234       {7.37854, 9.94},
0235       {7.37854, 7.1},
0236       {7.37854, 4.26},
0237       {7.37854, 1.42},
0238       {9.83805, 8.52},
0239       {9.83805, 5.68},
0240       {9.83805, 2.84},
0241       {12.2976, 4.26},
0242       {12.2976, 1.42},
0243       {12.2976, -4.26},
0244       {12.2976, -1.42},
0245       {9.83805, -8.52},
0246       {9.83805, -5.68},
0247       {9.83805, -2.84},
0248       {7.37854, -9.94},
0249       {7.37854, -7.1},
0250       {7.37854, -4.26},
0251       {7.37854, -1.42},
0252       {4.91902, -11.36},
0253       {4.91902, -8.52},
0254       {4.91902, -5.68},
0255       {2.45951, -12.78},
0256       {2.45951, -9.94},
0257       {2.45951, -7.1},
0258       {0, -11.36},
0259       {0, -8.52},
0260       {-2.45951, -12.78},
0261       {-2.45951, -9.94},
0262       {-2.45951, -7.1},
0263       {-4.91902, -11.36},
0264       {-4.91902, -8.52},
0265       {-4.91902, -5.68},
0266       {-7.37854, -9.94},
0267       {-7.37854, -7.1},
0268       {-7.37854, -4.26},
0269       {-7.37854, -1.42},
0270       {-9.83805, -8.52},
0271       {-9.83805, -5.68},
0272       {-9.83805, -2.84},
0273       {-12.2976, -4.26},
0274       {-12.2976, -1.42}};
0275 
0276   // m_bbcz = m_Params->get_double_param("z") * cm;
0277   m_bbcz = 250.0 * cm;  // The front face of the quartz is at 250 cm
0278 
0279   const float tube_zpos = m_bbcz + len_bbcd / 2.0 - len_bbca;
0280 
0281   const int NPMT = 64;  // No. PMTs per arm
0282   G4RotationMatrix *arm_rot[2];
0283   for (int iarm = 0; iarm < 2; iarm++)
0284   {
0285     arm_rot[iarm] = new G4RotationMatrix;
0286     float xside = 1.0;
0287     float zside = 1.0;
0288     if (iarm == 0)  // South Arm = 0
0289     {
0290       xside = 1.0;
0291       zside = -1.0;
0292       arm_rot[iarm]->rotateY(180 * deg);
0293     }
0294     else  // North Arm = 1
0295     {
0296       xside = -1.0;
0297       zside = 1.0;
0298     }
0299 
0300     // Add BBC PMT's
0301     for (int itube = 0; itube < NPMT; itube++)
0302     {
0303       // Full PMT Housing with Active Quartz Cerenkov Radiators
0304       float tube_xpos = xside * TubeLoc[itube][0] * cm;
0305       float tube_ypos = TubeLoc[itube][1] * cm;
0306       new G4PVPlacement(arm_rot[iarm], G4ThreeVector(tube_xpos, tube_ypos, zside * tube_zpos),
0307                         bbcd_lv, "BBCD", logicWorld, false, iarm * NPMT + itube, OverlapCheck());
0308     }
0309   }
0310 
0311   // Now Build the BBC Housing
0312 
0313   // BBC Outer and Inner Cylindrical Shells
0314   G4Tubs *bbc_outer_shell = new G4Tubs("bbc_outer_shell", 14.9 * cm, 15 * cm, 11.5 * cm, 0, 2 * M_PI);
0315   G4LogicalVolume *bbc_outer_shell_lv = new G4LogicalVolume(bbc_outer_shell, Aluminum, G4String("Bbc_Outer_Shell"));
0316   GetDisplayAction()->AddVolume(bbc_outer_shell_lv, "Bbc_Outer_Shell");
0317 
0318   G4Tubs *bbc_inner_shell = new G4Tubs("bbc_inner_shell", 5.0 * cm, 5.5 * cm, 11.5 * cm, 0, 2 * M_PI);
0319   G4LogicalVolume *bbc_inner_shell_lv = new G4LogicalVolume(bbc_inner_shell, Aluminum, G4String("Bbc_Inner_Shell"));
0320   GetDisplayAction()->AddVolume(bbc_inner_shell_lv, "Bbc_Inner_Shell");
0321 
0322   G4VPhysicalVolume *outer_shell_vol[2] = {nullptr};
0323   G4VPhysicalVolume *inner_shell_vol[2] = {nullptr};
0324 
0325   // Place South Shells
0326   outer_shell_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 - 11.5) * cm),
0327                                          bbc_outer_shell_lv, "BBC_OUTER_SHELL", logicWorld, false, 0, OverlapCheck());
0328   inner_shell_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 - 11.5) * cm),
0329                                          bbc_inner_shell_lv, "BBC_INNER_SHELL", logicWorld, false, 0, OverlapCheck());
0330 
0331   // Place North Shells
0332   outer_shell_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 + 11.5) * cm),
0333                                          bbc_outer_shell_lv, "BBC_OUTER_SHELL", logicWorld, false, 1, OverlapCheck());
0334   inner_shell_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 + 11.5) * cm),
0335                                          bbc_inner_shell_lv, "BBC_INNER_SHELL", logicWorld, false, 0, OverlapCheck());
0336   // this is more to prevent compiler warnings about unused variables
0337   if (!outer_shell_vol[0] || !outer_shell_vol[1] || !inner_shell_vol[0] || !inner_shell_vol[1])
0338   {
0339     std::cout << "problem placing BBC Sheels" << std::endl;
0340   }
0341 
0342   // BBC Front and Back Plates
0343   G4Tubs *bbc_plate = new G4Tubs("bbc_fplate", 5 * cm, 15 * cm, 0.5 * cm, 0, 2 * M_PI);
0344   G4LogicalVolume *bbc_plate_lv = new G4LogicalVolume(bbc_plate, Aluminum, G4String("Bbc_Cover_Plates"));
0345   GetDisplayAction()->AddVolume(bbc_plate_lv, "Bbc_Cover_Plates");
0346 
0347   G4VPhysicalVolume *fplate_vol[2] = {nullptr};  // Front Plates
0348   G4VPhysicalVolume *bplate_vol[2] = {nullptr};  // Back Plates
0349 
0350   // Place South Plates
0351   fplate_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 + 0.5) * cm),
0352                                     bbc_plate_lv, "BBC_FPLATE", logicWorld, false, 0, OverlapCheck());
0353   bplate_vol[0] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (-250 + 1.0 + 0.5 - 24.0) * cm),
0354                                     bbc_plate_lv, "BBC_BPLATE", logicWorld, false, 0, OverlapCheck());
0355 
0356   // Place North Plates
0357   fplate_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 - 0.5) * cm),
0358                                     bbc_plate_lv, "BBC_FPLATE", logicWorld, false, 1, OverlapCheck());
0359   bplate_vol[1] = new G4PVPlacement(nullptr, G4ThreeVector(0, 0, (250 - 1.0 - 0.5 + 24.0) * cm),
0360                                     bbc_plate_lv, "BBC_BPLATE", logicWorld, false, 0, OverlapCheck());
0361 
0362   // Place BBC Cables
0363   G4Material *Cu = manager->FindOrBuildMaterial("G4_Cu");
0364   const G4double len_cable = 120 * cm;
0365   const G4double r_CableConductor = 0.09525 * cm;
0366   G4Tubs *bbc_cablecond = new G4Tubs("bbc_cablecond", 0., r_CableConductor, len_cable * 0.5, 0 * deg, 360 * deg);
0367 
0368   G4LogicalVolume *bbc_cablecond_lv = new G4LogicalVolume(bbc_cablecond, Cu, G4String("Bbc_CableCond"));
0369   GetDisplayAction()->AddVolume(bbc_cablecond_lv, "Bbc_CableCond");
0370 
0371   const G4double rIn_CableShield = 0.302876 * cm;
0372   const G4double rOut_CableShield = 0.3175 * cm;
0373   G4Tubs *bbc_cableshield = new G4Tubs("bbc_cableshield", rIn_CableShield, rOut_CableShield, len_cable * 0.5, 0 * deg, 360 * deg);
0374 
0375   G4LogicalVolume *bbc_cableshield_lv = new G4LogicalVolume(bbc_cableshield, Cu, G4String("Bbc_CableShield"));
0376   GetDisplayAction()->AddVolume(bbc_cableshield_lv, "Bbc_CableShield");
0377 
0378   ypos = len_cable / 2 + 5 * cm;
0379 
0380   // For now we make this vertical, but they should be slanted toward the endcap
0381   G4RotationMatrix *rot_cable = new G4RotationMatrix();
0382   rot_cable->rotateX(90 * deg);
0383 
0384   int icable = 0;
0385 
0386   for (int iarm = 0; iarm < 2; iarm++)
0387   {
0388     float zsign = -1.0;
0389     if (iarm == 1)
0390     {
0391       zsign = 1.0;
0392     }
0393 
0394     for (int iring = 1; iring < 5; iring++)
0395     {
0396       float ring_radius = iring * 0.67 * cm;
0397       int ncables = 2 * M_PI * ring_radius / (0.635 * cm);
0398       double dphi = 2 * M_PI / ncables;
0399 
0400       // G4cout << "BBC_CABLE " << iring << "\t" << ring_radius << "\t" << ncables << "\t" << dphi*180/3.14 << G4endl;
0401 
0402       // place cables in ring
0403       for (int ic = 0; ic < ncables; ic++)
0404       {
0405         xpos = cos(dphi * ic) * ring_radius;
0406         zpos = sin(dphi * ic) * ring_radius + zsign * (m_bbcz + 30 * cm);
0407         // Place Inner Conductor
0408         new G4PVPlacement(rot_cable, G4ThreeVector(xpos, ypos, zpos), bbc_cablecond_lv, "BBC_Cable_Cond", logicWorld, false, icable, OverlapCheck());
0409         // Place Shield
0410         new G4PVPlacement(rot_cable, G4ThreeVector(xpos, ypos, zpos), bbc_cableshield_lv, "BBC_Cable_Shield", logicWorld, false, icable++, OverlapCheck());
0411       }
0412     }
0413   }
0414 
0415   // Now make the Supports
0416   ConstructSupport(logicWorld);
0417 
0418   // this is more to prevent compiler warnings about unused variables
0419   if (!fplate_vol[0] || !fplate_vol[1] || !bplate_vol[0] || !bplate_vol[1])
0420   {
0421     std::cout << "problem placing BBC Sheets" << std::endl;
0422   }
0423 
0424   // bbcq is the active detector element
0425   m_PhysLogicalVolSet.insert(bbcq_lv);
0426 
0427   //  GetDisplayAction()->AddVolume(bbc_plate_lv, "Bbc_plate");
0428   //  GetDisplayAction()->AddVolume(bbc_outer_shell_lv, "Bbc_oshell");
0429   //  GetDisplayAction()->AddVolume(bbc_inner_shell_lv, "Bbc_ishell");
0430 
0431   std::cout << "INSIDE BBC" << std::endl;
0432 
0433   return;
0434 }
0435 
0436 void PHG4BbcDetector::ConstructSupport(G4LogicalVolume *logicWorld)
0437 {
0438   // std::cout << "PHG4BbcDetector::ConstructSupport()" << std::endl;
0439 
0440   G4double fractionmass;
0441   G4int ncomponents;
0442   G4double density;
0443 
0444   recoConsts *rc = recoConsts::instance();
0445   G4Material *WorldMaterial = GetDetectorMaterial(rc->get_StringFlag("WorldMaterial"));
0446   G4Material *Delrin = GetDetectorMaterial("G4_POLYOXYMETHYLENE");
0447   // G4Material *Aluminum = GetDetectorMaterial("G4_Al");
0448 
0449   // McMaster-Carr #8548K36, implemented as PDG G10, https://pdg.lbl.gov/2022/AtomicNuclearProperties/HTML/G10.html
0450   G4Material *Fiberglass = new G4Material("BBC_Fiberglass", density = 1.800 * g / cm3, ncomponents = 9);
0451   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("B"), fractionmass = 0.018640);
0452   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Mg"), fractionmass = 0.010842);
0453   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("O"), fractionmass = 0.044453);
0454   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Al"), fractionmass = 0.151423);
0455   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("H"), fractionmass = 0.081496);
0456   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Ca"), fractionmass = 0.385764);
0457   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("C"), fractionmass = 0.275853);
0458   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("N"), fractionmass = 0.003583);
0459   Fiberglass->AddElement(G4NistManager::Instance()->FindOrBuildElement("Si"), fractionmass = 0.027945);
0460 
0461   // BBC Base Plates
0462   G4double basep_width = 35.56 * cm;
0463   G4double basep_height = 1.91 * cm;
0464   G4double basep_len = 46.99 * cm;
0465   G4double basep_zpos = 228.6 * cm;  // z-pos of front edge of base plate
0466   G4Box *bbc_base_plate = new G4Box("bbc_base_plate", basep_width / 2, basep_height / 2, basep_len / 2);
0467   G4LogicalVolume *bbc_base_plate_lv = new G4LogicalVolume(bbc_base_plate, Delrin, G4String("Bbc_Base_Plates"));
0468   GetDisplayAction()->AddVolume(bbc_base_plate_lv, "Bbc_Base_Plates");
0469 
0470   // Place South Base Plates
0471   new G4PVPlacement(nullptr, G4ThreeVector(0 * cm, -15. * cm - basep_height / 2, (-basep_zpos - basep_len / 2)),
0472                     bbc_base_plate_lv, "BBC_BASE_PLATE", logicWorld, false, 0, OverlapCheck());
0473 
0474   // Place North Base Plates
0475   new G4PVPlacement(nullptr, G4ThreeVector(0 * cm, -15. * cm - basep_height / 2, (basep_zpos + basep_len / 2)),
0476                     bbc_base_plate_lv, "BBC_BASE_PLATE", logicWorld, false, 1, OverlapCheck());
0477 
0478   // BBC Side Support Plates
0479   G4double sidesupportp_width = 1.27 * cm;
0480   G4double sidesupportp_height = 14.57 * cm;
0481   G4double sidesupportp_len = 25.00 * cm;
0482 
0483   G4Box *bbc_sidesupport_plate = new G4Box("bbc_sidesupport_plate", sidesupportp_width / 2, sidesupportp_height / 2, sidesupportp_len / 2);
0484   G4LogicalVolume *bbc_sidesupport_plate_lv = new G4LogicalVolume(bbc_sidesupport_plate, Delrin, G4String("Bbc_Sidesupport_Plates"));
0485   GetDisplayAction()->AddVolume(bbc_sidesupport_plate_lv, "Bbc_Sidesupport_Plates");
0486 
0487   // Make and Place Holes in Side Support Plates
0488   G4Tubs *bbc_sidesupport_hole = new G4Tubs("bbc_sidesupport_hole", 0., (7.62 / 2) * cm, sidesupportp_width / 2, 0, 2.0 * M_PI);
0489   G4LogicalVolume *bbc_sidesupport_hole_lv = new G4LogicalVolume(bbc_sidesupport_hole, WorldMaterial, G4String("Bbc_Sidesupport_Holes"));
0490 
0491   G4RotationMatrix *rot_sideholes = new G4RotationMatrix;
0492   rot_sideholes->rotateY(90. * deg);
0493   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());
0494   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());
0495 
0496   // Place South Side Support Plates
0497   new G4PVPlacement(nullptr, G4ThreeVector(-15 * cm - sidesupportp_width / 2, -15 * cm + sidesupportp_height / 2, -front_bbcz - sidesupportp_len / 2),
0498                     bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 0, OverlapCheck());
0499   new G4PVPlacement(nullptr, G4ThreeVector(15 * cm + sidesupportp_width / 2, -15 * cm + sidesupportp_height / 2, -front_bbcz - sidesupportp_len / 2),
0500                     bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 1, OverlapCheck());
0501 
0502   // Place North Side Support Plate
0503   new G4PVPlacement(nullptr, G4ThreeVector(-15 * cm - sidesupportp_width / 2, -15 * cm + sidesupportp_height / 2, front_bbcz + sidesupportp_len / 2),
0504                     bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 2, OverlapCheck());
0505   new G4PVPlacement(nullptr, G4ThreeVector(15 * cm + sidesupportp_width / 2, -15 * cm + sidesupportp_height / 2, front_bbcz + sidesupportp_len / 2),
0506                     bbc_sidesupport_plate_lv, "BBC_SIDESUPPORT_PLATE", logicWorld, false, 3, OverlapCheck());
0507 
0508   // BBC Support Post
0509   G4double supportp_width = 10.16 * cm;
0510   G4double supportp_height = 97.16 * cm;
0511   G4double supportp_len = 10.16 * cm;
0512   G4double supportp_thick = 0.64 * cm;
0513 
0514   G4Box *bbc_support_post_outside = new G4Box("bbc_support_post_outside", supportp_width / 2, supportp_height / 2, supportp_len / 2);
0515   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);
0516   G4SubtractionSolid *bbc_support_post = new G4SubtractionSolid("support_post", bbc_support_post_outside, bbc_support_post_inside);
0517 
0518   G4LogicalVolume *bbc_support_post_lv = new G4LogicalVolume(bbc_support_post, Fiberglass, G4String("Bbc_Support_Post"));
0519 
0520   GetDisplayAction()->AddVolume(bbc_support_post_lv, "Bbc_Support_Post");
0521 
0522   // Place South Support Post
0523   new G4PVPlacement(nullptr, G4ThreeVector(0, -15. * cm - basep_height - supportp_height / 2, -228. * cm - supportp_len / 2),
0524                     bbc_support_post_lv, "BBC_SUPPORT_POST", logicWorld, false, 0, OverlapCheck());
0525 
0526   // Place North Support Post
0527   new G4PVPlacement(nullptr, G4ThreeVector(0, -15. * cm - basep_height - supportp_height / 2, 228. * cm + supportp_len / 2),
0528                     bbc_support_post_lv, "BBC_SUPPORT_POST", logicWorld, false, 1, OverlapCheck());
0529 
0530   // BBC Support Arm
0531   G4double supporta_width = 148.59 * cm;
0532   G4double supporta_height = 10.16 * cm;
0533   G4double supporta_len = 10.16 * cm;
0534   G4double supporta_thick = 0.64 * cm;
0535 
0536   G4Box *bbc_support_arm_outside = new G4Box("bbc_support_arm_outside", supporta_width / 2, supporta_height / 2, supporta_len / 2);
0537   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);
0538   G4SubtractionSolid *bbc_support_arm = new G4SubtractionSolid("support_arm", bbc_support_arm_outside, bbc_support_arm_inside);
0539 
0540   G4LogicalVolume *bbc_support_arm_lv = new G4LogicalVolume(bbc_support_arm, Fiberglass, G4String("Bbc_Support_Arm"));
0541 
0542   GetDisplayAction()->AddVolume(bbc_support_arm_lv, "Bbc_Support_Arm");
0543 
0544   // Place South Support Arms
0545   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),
0546                     bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 0, OverlapCheck());
0547   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),
0548                     bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 1, OverlapCheck());
0549 
0550   // Place North Support Arms
0551   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),
0552                     bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 2, OverlapCheck());
0553   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),
0554                     bbc_support_arm_lv, "BBC_SUPPORT_ARM", logicWorld, false, 3, OverlapCheck());
0555 
0556   // BBC Gusset Plates (implementation is broken up into 3 parts, in Trap and 2 Boxes)
0557   G4double gussetp0_pz = 1.27 * cm;
0558   G4double gussetp0_py = (45.72 - 11.11) * cm;
0559   G4double gussetp0_px = 44.62 * cm;  // measured off drawing
0560   G4double gussetp0_pltx = 5.08 * cm;
0561   G4Trap *bbc_gusset_plate0 = new G4Trap("bbc_gusset_plate0", gussetp0_pz, gussetp0_py, gussetp0_px, gussetp0_pltx);
0562   G4LogicalVolume *bbc_gusset0_plate_lv = new G4LogicalVolume(bbc_gusset_plate0, Delrin, G4String("Bbc_Gusset0_Plates"));
0563   GetDisplayAction()->AddVolume(bbc_gusset0_plate_lv, "Bbc_Gusset0_Plates");
0564 
0565   // Place South Gusset Plates (Trapezoid part)
0566   G4RotationMatrix *rot_sgusset = new G4RotationMatrix;
0567   rot_sgusset->rotateY(90. * deg);
0568   rot_sgusset->rotateZ(180. * deg);
0569 
0570   G4double xpos = supportp_width / 2 + gussetp0_pz / 2;
0571   G4double ypos = -gussetp0_py / 2 - 15 * cm - basep_height;
0572   G4double zpos = 0.25 * (gussetp0_px + gussetp0_pltx) + 228.6 * cm + 11.11 * cm;
0573   new G4PVPlacement(rot_sgusset, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 0, OverlapCheck());
0574   new G4PVPlacement(rot_sgusset, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 1, OverlapCheck());
0575 
0576   // Place North Gusset Plates (Trapezoid part)
0577   G4RotationMatrix *rot_ngusset = new G4RotationMatrix;
0578   rot_ngusset->rotateY(-90. * deg);
0579   rot_ngusset->rotateZ(180. * deg);
0580 
0581   new G4PVPlacement(rot_ngusset, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 2, OverlapCheck());
0582   new G4PVPlacement(rot_ngusset, G4ThreeVector(xpos, ypos, zpos), bbc_gusset0_plate_lv, "BBC_GUSSET_PLATE0", logicWorld, false, 3, OverlapCheck());
0583 
0584   // Gusset Plate Box 1
0585   G4double gussetp1_x = 1.27 * cm;  // measured off drawing
0586   G4double gussetp1_y = 20.48 * cm;
0587   G4double gussetp1_z = 11.11 * cm;
0588 
0589   G4Box *bbc_gusset_plate1 = new G4Box("bbc_gusset_plate1", gussetp1_x / 2, gussetp1_y / 2, gussetp1_z / 2);
0590   G4LogicalVolume *bbc_gusset1_plate_lv = new G4LogicalVolume(bbc_gusset_plate1, Delrin, G4String("Bbc_Gusset1_Plates"));
0591   GetDisplayAction()->AddVolume(bbc_gusset1_plate_lv, "Bbc_Gusset1_Plates");
0592 
0593   // Place South Gusset Plates (Box 1)
0594   ypos = -15 * cm - basep_height - gussetp1_y / 2;
0595   zpos = 228.6 * cm + gussetp1_z / 2;
0596   new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 0, OverlapCheck());
0597   new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 1, OverlapCheck());
0598 
0599   // Place North Gusset Plates (Box 1)
0600   new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 2, OverlapCheck());
0601   new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbc_gusset1_plate_lv, "BBC_GUSSET_PLATE1", logicWorld, false, 3, OverlapCheck());
0602 
0603   // Gusset Plate Box 2
0604   G4double gussetp2_x = 1.27 * cm;  // measured off drawing
0605   G4double gussetp2_y = (45.72 - 10.80 - 20.48) * cm;
0606   G4double gussetp2_z = 11.11 * cm;
0607 
0608   G4Box *bbc_gusset_plate2 = new G4Box("bbc_gusset_plate2", gussetp2_x / 2, gussetp2_y / 2, gussetp2_z / 2);
0609   G4LogicalVolume *bbc_gusset2_plate_lv = new G4LogicalVolume(bbc_gusset_plate2, Delrin, G4String("Bbc_Gusset2_Plates"));
0610   GetDisplayAction()->AddVolume(bbc_gusset2_plate_lv, "Bbc_Gusset2_Plates");
0611 
0612   // Place South Gusset Plates (Box 2)
0613   ypos = -15 * cm - basep_height - 20.48 * cm - 10.80 * cm - gussetp2_y / 2;
0614   zpos = 228.6 * cm + gussetp2_z / 2;
0615   new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, -zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 0, OverlapCheck());
0616   new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, -zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 1, OverlapCheck());
0617 
0618   // Place North Gusset Plates (Box 2)
0619   new G4PVPlacement(nullptr, G4ThreeVector(-xpos, ypos, zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 2, OverlapCheck());
0620   new G4PVPlacement(nullptr, G4ThreeVector(xpos, ypos, zpos), bbc_gusset2_plate_lv, "BBC_GUSSET_PLATE2", logicWorld, false, 3, OverlapCheck());
0621 
0622   // BBC Splice Plates
0623   G4double splicep_x = 35.56 * cm;
0624   G4double splicep_y = 10.16 * cm;
0625   G4double splicep_z = 0.64 * cm;
0626 
0627   G4Box *bbc_splice_plate = new G4Box("bbc_splice_plate", splicep_x / 2, splicep_y / 2, splicep_z / 2);
0628   G4LogicalVolume *bbc_splice_plate_lv = new G4LogicalVolume(bbc_splice_plate, Delrin, G4String("Bbc_Splice_Plates"));
0629   GetDisplayAction()->AddVolume(bbc_splice_plate_lv, "Bbc_Splice_Plates");
0630 
0631   // Make and Place Holes in Splice Plates
0632   G4Tubs *bbc_splice_hole = new G4Tubs("bbc_splice_hole", 0., (6.35 / 2) * cm, splicep_z / 2, 0, 2.0 * M_PI);
0633   G4LogicalVolume *bbc_splice_hole_lv = new G4LogicalVolume(bbc_splice_hole, WorldMaterial, G4String("Bbc_Splice_Holes"));
0634 
0635   new G4PVPlacement(nullptr, G4ThreeVector(-12.7 * cm, 0, 0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 0, OverlapCheck());
0636   new G4PVPlacement(nullptr, G4ThreeVector(0, 0, 0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 1, OverlapCheck());
0637   new G4PVPlacement(nullptr, G4ThreeVector(12.7 * cm, 0, 0), bbc_splice_hole_lv, "BBC_SPLICE_HOLE", bbc_splice_plate_lv, false, 2, OverlapCheck());
0638 
0639   // Place South Splice Plate
0640   ypos = -15 * cm - basep_height - 20.48 * cm - splicep_y / 2;
0641   zpos = 228.6 * cm + supportp_len + splicep_z / 2;
0642   new G4PVPlacement(nullptr, G4ThreeVector(0, ypos, -zpos), bbc_splice_plate_lv, "BBC_SPLICE_PLATE", logicWorld, false, 0, OverlapCheck());
0643 
0644   // Place North Splice Plate
0645   new G4PVPlacement(nullptr, G4ThreeVector(0, ypos, zpos), bbc_splice_plate_lv, "BBC_SPLICE_PLATE", logicWorld, false, 1, OverlapCheck());
0646 }
0647 
0648 void PHG4BbcDetector::Print(const std::string &what) const
0649 {
0650   std::cout << "Bbc Detector:" << std::endl;
0651   if (what == "ALL" || what == "VOLUME")
0652   {
0653     std::cout << "Version 0.1" << std::endl;
0654   }
0655   return;
0656 }