Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:17:52

0001 #include "PHG4sPHENIXMagnetDetector.h"
0002 
0003 #include "PHG4sPHENIXMagnetDisplayAction.h"
0004 
0005 #include <phparameter/PHParameters.h>
0006 
0007 #include <g4gdml/PHG4GDMLUtility.hh>
0008 
0009 #include <g4main/PHG4Detector.h>       // for PHG4Detector
0010 #include <g4main/PHG4DisplayAction.h>  // for PHG4DisplayAction
0011 #include <g4main/PHG4Subsystem.h>
0012 
0013 #include <phool/recoConsts.h>
0014 
0015 #include <TSystem.h>
0016 
0017 #include <Geant4/G4Box.hh>
0018 #include <Geant4/G4LogicalVolume.hh>
0019 #include <Geant4/G4PVPlacement.hh>
0020 #include <Geant4/G4PhysicalConstants.hh>
0021 #include <Geant4/G4Polycone.hh>
0022 #include <Geant4/G4RotationMatrix.hh>  // for G4RotationMatrix
0023 #include <Geant4/G4String.hh>          // for G4String
0024 #include <Geant4/G4SubtractionSolid.hh>
0025 #include <Geant4/G4SystemOfUnits.hh>
0026 #include <Geant4/G4ThreeVector.hh>  // for G4ThreeVector
0027 #include <Geant4/G4Transform3D.hh>  // for G4Transform3D
0028 #include <Geant4/G4Tubs.hh>
0029 #include <Geant4/G4Types.hh>            // for G4double, G4int
0030 #include <Geant4/G4VPhysicalVolume.hh>  // for G4VPhysicalVolume
0031 
0032 #include <cassert>
0033 #include <cmath>
0034 #include <iostream>
0035 #include <sstream>
0036 
0037 class G4Material;
0038 class G4VSolid;
0039 class PHCompositeNode;
0040 
0041 //_______________________________________________________________________
0042 PHG4sPHENIXMagnetDetector::PHG4sPHENIXMagnetDetector(PHG4Subsystem* subsys, PHCompositeNode* Node, PHParameters* parameters, const std::string& dnam, const int detid)
0043   : PHG4Detector(subsys, Node, dnam)
0044   , m_DisplayAction(dynamic_cast<PHG4sPHENIXMagnetDisplayAction*>(subsys->GetDisplayAction()))
0045   , m_Params(parameters)
0046   , m_GdmlConfig(PHG4GDMLUtility::GetOrMakeConfigNode(Node))
0047   , m_ActiveFlag(m_Params->get_int_param("active"))
0048   , m_Layer(detid)
0049   , m_SuperDetector("NONE")
0050 {
0051   assert(m_GdmlConfig);
0052 }
0053 
0054 //_______________________________________________________________________
0055 int PHG4sPHENIXMagnetDetector::IsInsPHENIXMagnet(G4VPhysicalVolume* volume) const
0056 {
0057   G4LogicalVolume* mylogvol = volume->GetLogicalVolume();
0058 
0059   if (m_LogicalVolSet.find(mylogvol) != m_LogicalVolSet.end())
0060   {
0061     return 1;
0062   }
0063 
0064   return 0;
0065 }
0066 
0067 //_______________________________________________________________________
0068 void PHG4sPHENIXMagnetDetector::ConstructMe(G4LogicalVolume* logicWorld)
0069 {
0070   if (Verbosity() > 0)
0071   {
0072     std::cout << "PHG4sPHENIXMagnetDetector: Begin Construction" << std::endl;
0073   }
0074 
0075   recoConsts* rc = recoConsts::instance();
0076   G4Material* WorldMaterial = GetDetectorMaterial(rc->get_StringFlag("WorldMaterial"));
0077   G4Material* aluminum = GetDetectorMaterial("G4_Al");
0078   G4Material* copper = GetDetectorMaterial("G4_Cu");
0079 
0080   // Mother volume for solenoid
0081 
0082   G4Tubs* cryostat = SolenoidTubes(0);
0083   G4LogicalVolume* logCryostat =
0084       new G4LogicalVolume(cryostat, aluminum, "CRYOSTAT");
0085   GetDisplayAction()->AddVolume(logCryostat, "CRYOSTAT");
0086   m_LogicalVolSet.insert(logCryostat);
0087   new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logCryostat, "CRYOSTAT", logicWorld, false, false, OverlapCheck());
0088 
0089   // Air (or vacuum?) inside cryostat
0090 
0091   G4Polycone* cryostatInterior = SolenoidPolycones(0);
0092   G4LogicalVolume* logCryostatInterior =
0093       new G4LogicalVolume(cryostatInterior, WorldMaterial, "CRYOINT");
0094 
0095   GetDisplayAction()->AddVolume(logCryostatInterior, "CRYOINT");
0096   m_LogicalVolSet.insert(logCryostatInterior);
0097   new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logCryostatInterior, "CRYOINT", logCryostat, false, false, OverlapCheck());
0098 
0099   // Thermal shield
0100 
0101   G4Tubs* therm = SolenoidTubes(1);
0102   G4LogicalVolume* logTherm =
0103       new G4LogicalVolume(therm, aluminum, "THERM");
0104   GetDisplayAction()->AddVolume(logTherm, "THERM");
0105   m_LogicalVolSet.insert(logTherm);
0106   new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logTherm, "THERM", logCryostatInterior, false, false, OverlapCheck());
0107 
0108   G4Tubs* thermvac = SolenoidTubes(2);
0109   G4LogicalVolume* logThermVac =
0110       new G4LogicalVolume(thermvac, WorldMaterial, "THERMVAC");
0111   GetDisplayAction()->AddVolume(logThermVac, "THERMVAC");
0112   m_LogicalVolSet.insert(logThermVac);
0113   new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logThermVac, "THERMVAC", logTherm, false, false, OverlapCheck());
0114 
0115   // Coil support
0116 
0117   G4Polycone* coilSupport = SolenoidPolycones(1);
0118   G4LogicalVolume* logCoilSupport =
0119       new G4LogicalVolume(coilSupport, aluminum, "COILSUP");
0120   GetDisplayAction()->AddVolume(logCoilSupport, "COILSUP");
0121   m_LogicalVolSet.insert(logCoilSupport);
0122   new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logCoilSupport, "COILSUP", logThermVac, false, false, OverlapCheck());
0123 
0124   // Coil
0125 
0126   G4Tubs* coil = SolenoidTubes(3);
0127   G4LogicalVolume* logCoil =
0128       new G4LogicalVolume(coil, aluminum, "COIL");
0129   GetDisplayAction()->AddVolume(logCoil, "COIL");
0130   m_LogicalVolSet.insert(logCoil);
0131   G4double zpos = 1.25 * cm;  // ifr_inert.dat file line 185
0132   new G4PVPlacement(0, G4ThreeVector(0.0, 0.0, zpos),
0133                     logCoil, "COIL", logThermVac, false, false, OverlapCheck());
0134 
0135   // Tube connecting solenoid to cryotower
0136 
0137   G4Tubs* connector = CryoTubes(3);
0138   G4LogicalVolume* logConnector =
0139       new G4LogicalVolume(connector, aluminum, "CONNECTOR");
0140   GetDisplayAction()->AddVolume(logConnector, "CONNECTOR");
0141   m_LogicalVolSet.insert(logConnector);
0142 
0143   G4double ycon = 158.5 * cm;   // ifr_inert.dat file line 189
0144   G4double zcon = -197.5 * cm;  // ifr_inert.dat file line 189
0145   new G4PVPlacement(0, G4ThreeVector(0.0, ycon, zcon),
0146                     logConnector, "CONNECTOR", logicWorld, false, false, OverlapCheck());
0147 
0148   // Bus bar in connecting tube
0149 
0150   G4Box* busbarD = Block(14);
0151   G4LogicalVolume* logBusbarD =
0152       new G4LogicalVolume(busbarD, copper, "BusbarD", 0, 0, 0);
0153   GetDisplayAction()->AddVolume(logBusbarD, "BusbarD");
0154   m_LogicalVolSet.insert(logBusbarD);
0155 
0156   G4double yD = 158.5 * cm;   // ifr_inert.dat file line 75
0157   G4double zD = -197.5 * cm;  // ifr_inert.dat file line 75
0158   G4RotationMatrix rotationD;
0159   rotationD.rotateZ(90. * deg);
0160   rotationD.rotateX(90. * deg);
0161   G4Transform3D transformD(rotationD, G4ThreeVector(0.0, yD, zD));
0162   new G4PVPlacement(transformD,
0163                     logBusbarD, "BusbarD", logicWorld, false, false, OverlapCheck());
0164 
0165   return;
0166 }
0167 
0168 G4Tubs* PHG4sPHENIXMagnetDetector::SolenoidTubes(G4int itube) const
0169 {
0170   G4double tubes_r_inner[8] = {142.0, 147.6, 148.1, 151.0, 12.44, 10.15, 0.0, 10.86};    // ifr_inert.dat file lines 181-189
0171   G4double tubes_r_outer[8] = {177.0, 168.0, 167.5, 155.08, 13.35, 10.50, 0.79, 12.37};  // ifr_inert.dat file lines 181-189
0172   G4double tubes_length[8] = {385.0, 371.0, 370.0, 351.29, 121.0, 121.0, 121.0, 9.98};   // ifr_inert.dat file lines 181-189
0173 
0174   G4double r_inner = tubes_r_inner[itube] * cm;
0175   G4double r_outer = tubes_r_outer[itube] * cm;
0176   G4double half_length = 0.5 * tubes_length[itube] * cm;
0177 
0178   return new G4Tubs("SolenoidTube", r_inner, r_outer, half_length,
0179                     0.0 * deg, 360.0 * deg);
0180 }
0181 
0182 G4Polycone* PHG4sPHENIXMagnetDetector::SolenoidPolycones(G4int ipolycone) const
0183 {
0184   G4double zPlane[10];
0185   G4double rInner[10];
0186   G4double rOuter[10];
0187 
0188   G4int nzplanes = 10;  // ifr_inert.dat file lines 203-211
0189   // ifr_inert.dat file lines 203-211
0190   G4double polycones_zplane[2][10] = {{-188.5, -179.5, -172.0, -159.5, -150.2, 155.1, 164.4, 172.0, 179.5, 188.5},
0191                                       {-176.9, -174.4, -174.4, -164.0, -158.4, 163.4, 169.0, 176.9, 176.9, 181.9}};
0192   G4double polycones_r_inner[2][10] = {{145.5, 145.5, 143.0, 143.0, 143.0, 143.0, 143.0, 143.0, 145.5, 145.5},
0193                                        {151.0, 151.0, 155.1, 155.1, 155.1, 155.1, 155.1, 155.1, 151.0, 151.0}};
0194   G4double polycones_r_outer[2][10] = {{172.0, 172.0, 172.0, 172.0, 174.5, 174.5, 172.0, 172.0, 172.0, 172.0},
0195                                        {160.4, 160.4, 160.4, 160.4, 158.6, 158.6, 160.4, 160.4, 160.4, 160.4}};
0196   G4int iz;
0197   for (iz = 0; iz < nzplanes; iz++)
0198   {
0199     zPlane[iz] = polycones_zplane[ipolycone][iz] * cm;
0200     rInner[iz] = polycones_r_inner[ipolycone][iz] * cm;
0201     rOuter[iz] = polycones_r_outer[ipolycone][iz] * cm;
0202   }
0203 
0204   return new G4Polycone("SolenoidPolycone", 0.0 * deg, 360.0 * deg, nzplanes,
0205                         zPlane, rInner, rOuter);
0206 }
0207 
0208 G4Tubs* PHG4sPHENIXMagnetDetector::CryoTubes(G4int itube) const
0209 {
0210   G4double tubes_r_inner[8] = {142.0, 147.6, 148.1, 151.0, 12.44, 10.15, 0.0, 10.86};    // ifr_inert.dat file lines 181-189
0211   G4double tubes_r_outer[8] = {177.0, 168.0, 167.5, 155.08, 13.35, 10.50, 0.79, 12.37};  // ifr_inert.dat file lines 181-189
0212   G4double tubes_length[8] = {385.0, 371.0, 370.0, 351.29, 121.0, 121.0, 121.0, 9.98};   // ifr_inert.dat file lines 181-189
0213 
0214   G4double r_inner = tubes_r_inner[itube + 4] * cm;
0215   G4double r_outer = tubes_r_outer[itube + 4] * cm;
0216   G4double half_length = 0.5 * tubes_length[itube + 4] * cm;
0217 
0218   return new G4Tubs("CryoTube", r_inner, r_outer, half_length,
0219                     0.0 * deg, 360.0 * deg);
0220 }
0221 
0222 G4Box* PHG4sPHENIXMagnetDetector::Block(G4int iblock) const
0223 {
0224   G4double blocks_length[15] = {375.0, 138.9, 300.0, 315.0, 128.9, 128.9, 166.84, 160.84, 128.9, 128.9, 251.5, 121.0, 54.4, 14.6, 9.98};  // ifr_inert.dat file lines 60-75
0225   G4double blocks_width[15] = {172.52, 80.0, 290.0, 290.0, 104.9, 104.9, 104.9, 104.9, 114.9, 114.9, 15.0, 6.1, 6.1, 6.1, 6.1};           // ifr_inert.dat file lines 60-75
0226   G4double blocks_thickness[15] = {5.5, 15.0, 25.0, 25.0, 8.5, 8.5, 10.0, 10.0, 25.0, 25.0, 15.0, 3.4, 3.4, 3.4, 3.4};                    // ifr_inert.dat file lines 60-75
0227 
0228   G4double dx = 0.5 * blocks_length[iblock] * cm;
0229   G4double dy = 0.5 * blocks_width[iblock] * cm;
0230   G4double dz = 0.5 * blocks_thickness[iblock] * cm;
0231 
0232   return new G4Box("IfrBlock", dx, dy, dz);
0233 }