Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:58

0001 #include "PHG4InttDetector.h"
0002 
0003 #include "PHG4InttDefs.h"  // for SEGMENTATION_Z
0004 #include "PHG4InttDisplayAction.h"
0005 #include "PHG4InttFPHXParameterisation.h"
0006 
0007 #include <intt/CylinderGeomIntt.h>
0008 #include <intt/InttMapping.h>
0009 #include <intt/InttSurveyMap.h>
0010 
0011 #include <g4detectors/PHG4CylinderGeomContainer.h>
0012 
0013 #include <phparameter/PHParameters.h>
0014 #include <phparameter/PHParametersContainer.h>
0015 
0016 #include <g4main/PHG4Detector.h>       // for PHG4Detector
0017 #include <g4main/PHG4DisplayAction.h>  // for PHG4DisplayAction
0018 #include <g4main/PHG4Subsystem.h>      // for PHG4Subsystem
0019 
0020 #include <ffamodules/CDBInterface.h>
0021 
0022 #include <phool/PHCompositeNode.h>
0023 #include <phool/PHIODataNode.h>
0024 #include <phool/PHNode.h>          // for PHNode
0025 #include <phool/PHNodeIterator.h>  // for PHNodeIterator
0026 #include <phool/PHObject.h>        // for PHObject
0027 #include <phool/getClass.h>
0028 #include <phool/phool.h>  // for PHWHERE
0029 #include <phool/recoConsts.h>
0030 
0031 #include <TSystem.h>
0032 
0033 #include <Geant4/G4Box.hh>
0034 #include <Geant4/G4GenericTrap.hh>
0035 #include <Geant4/G4LogicalVolume.hh>
0036 #include <Geant4/G4PVParameterised.hh>
0037 #include <Geant4/G4PVPlacement.hh>
0038 #include <Geant4/G4RotationMatrix.hh>  // for G4RotationMatrix
0039 #include <Geant4/G4String.hh>          // for G4String
0040 #include <Geant4/G4SubtractionSolid.hh>
0041 #include <Geant4/G4SystemOfUnits.hh>
0042 #include <Geant4/G4ThreeVector.hh>  // for G4ThreeVector
0043 #include <Geant4/G4Transform3D.hh>  // for G4Transform3D
0044 #include <Geant4/G4Tubs.hh>
0045 #include <Geant4/G4TwoVector.hh>        // for G4TwoVector
0046 #include <Geant4/G4VPhysicalVolume.hh>  // for G4VPhysicalVolume
0047 #include <Geant4/geomdefs.hh>           // for kZAxis
0048 
0049 #include <algorithm>  // for fill_n
0050 #include <array>
0051 #include <cmath>
0052 #include <cstdlib>  // for exit, NULL
0053 #include <filesystem>
0054 #include <format>
0055 #include <iostream>  // for operator<<, basic...
0056 
0057 class G4VPVParameterisation;
0058 class G4VSolid;
0059 
0060 PHG4InttDetector::PHG4InttDetector(PHG4Subsystem *subsys, PHCompositeNode *Node, PHParametersContainer *parameters, const std::string &dnam,
0061                                    const std::pair<std::vector<std::pair<int, int>>::const_iterator, std::vector<std::pair<int, int>>::const_iterator> &layer_b_e)
0062   : PHG4Detector(subsys, Node, dnam)
0063   , m_DisplayAction(dynamic_cast<PHG4InttDisplayAction *>(subsys->GetDisplayAction()))
0064   , m_ParamsContainer(parameters)
0065   , m_LayerBeginEndIteratorPair(layer_b_e)
0066 {
0067   for (auto layeriter = m_LayerBeginEndIteratorPair.first; layeriter != m_LayerBeginEndIteratorPair.second; ++layeriter)
0068   {
0069     int layer = layeriter->second;
0070     const PHParameters *par = m_ParamsContainer->GetParameters(layer);
0071     m_IsActiveMap.insert(std::make_pair(layer, par->get_int_param("active")));
0072     m_IsAbsorberActiveMap.insert(std::make_pair(layer, par->get_int_param("absorberactive")));
0073   }
0074   const PHParameters *par = m_ParamsContainer->GetParameters(PHG4InttDefs::SUPPORTPARAMS);
0075   m_IsSupportActive = par->get_int_param("supportactive");
0076   m_IsEndcapActive = par->get_int_param("endcap_ring_enabled");
0077   std::fill_n(&m_PosZ[0][0], sizeof(m_PosZ) / sizeof(double), std::numeric_limits<double>::quiet_NaN());
0078   std::fill_n(m_SensorRadius, sizeof(m_SensorRadius) / sizeof(double), std::numeric_limits<double>::quiet_NaN());
0079   std::fill_n(m_StripOffsetX, sizeof(m_StripOffsetX) / sizeof(double), std::numeric_limits<double>::quiet_NaN());
0080 }
0081 
0082 //_______________________________________________________________
0083 int PHG4InttDetector::IsInIntt(G4VPhysicalVolume *volume) const
0084 {
0085   // Is this volume one of the sensor strips?
0086   // just checking if the pointer to the logical volume is in the set
0087   // of our active/active absorber ones makes sure we are in an active volume
0088   // name parsing is a bad idea since this is called for all steps
0089   // and we would have to trust that people give different names
0090   // to their volumes
0091   G4LogicalVolume *logvol = volume->GetLogicalVolume();
0092   if (!m_PassiveVolumeTuple.empty() && m_PassiveVolumeTuple.contains(logvol))
0093   {
0094     return -1;
0095   }
0096   if (m_ActiveLogVols.contains(logvol))
0097   {
0098     return 1;
0099   }
0100 
0101   return 0;
0102 }
0103 
0104 void PHG4InttDetector::ConstructMe(G4LogicalVolume *logicWorld)
0105 {
0106   if (Verbosity() > 0)
0107   {
0108     std::cout << "PHG4InttDetector::Construct called for layers " << std::endl;
0109     for (auto layeriter = m_LayerBeginEndIteratorPair.first; layeriter != m_LayerBeginEndIteratorPair.second; ++layeriter)
0110     {
0111       std::cout << "layer " << layeriter->second << std::endl;
0112     }
0113   }
0114   // the tracking layers are placed directly in the world volume, since some layers are (touching) double layers
0115   ConstructIntt(logicWorld);
0116 
0117   // This object provides the strip center locations when given the ladder segment and strip internal cordinates in the sensor
0118   AddGeometryNode();
0119   return;
0120 }
0121 
0122 int PHG4InttDetector::ConstructIntt(G4LogicalVolume *trackerenvelope)
0123 {
0124   recoConsts *rc = recoConsts::instance();  // use for worldmaterial in a few places
0125   // We have an arbitray number of layers (nlayer_) up to 8
0126   // We have 2 types of ladders (vertical strips and horizontal strips)
0127   // We have 2 types of sensors (inner and outer)
0128   std::array<std::array<double, 2>, 8> hdi_z_arr{};
0129   // we loop over layers. All layers have only one laddertype
0130 
0131   // INTT survey map to get the survey information
0132   InttSurveyMap *survey = new InttSurveyMap();
0133   if (useSurvey)
0134   {
0135     std::string url = CDBInterface::instance()->getUrl("InttSurveyMap");
0136     if (!std::filesystem::exists(url))
0137     {
0138       std::cout << PHWHERE << " Could not locate INTT survey geometry " << url << std::endl;
0139       gSystem->Exit(1);
0140     }
0141     std::cout << PHWHERE << "Use the INTT survey geometry. Get the survey map from " << url << std::endl;
0142     if (survey->LoadFromFile(url))
0143     {
0144       std::cout << PHWHERE << "Failed to load INTT survey geometry from the CDB" << std::endl;
0145       gSystem->Exit(1);
0146     }
0147     if (Verbosity() > 0)
0148     {
0149       survey->identify();
0150     }
0151   }
0152   else
0153   {
0154     std::cout << PHWHERE << "Use the default INTT ideal geometry." << std::endl;
0155   }
0156 
0157   for (auto layeriter = m_LayerBeginEndIteratorPair.first; layeriter != m_LayerBeginEndIteratorPair.second; ++layeriter)
0158   {
0159     int inttlayer = layeriter->second;
0160     // get the parameters for this layer
0161     const PHParameters *params1 = m_ParamsContainer->GetParameters(inttlayer);
0162     const int laddertype = params1->get_int_param("laddertype");
0163     const double offsetphi = (params1->get_double_param("offsetphi") * deg) / rad;  // use rad internally
0164     double offsetrot = (params1->get_double_param("offsetrot") * deg) / rad;        // offsetrot is specified in deg, we convert to rad here
0165     m_SensorRadius[inttlayer] = params1->get_double_param("sensor_radius") * cm;
0166     const int nladders_layer = params1->get_int_param("nladder");
0167 
0168     // Look up all remaining parameters by the laddertype for this layer
0169     const PHParameters *params = m_ParamsContainer->GetParameters(laddertype);
0170     const double strip_x = params->get_double_param("strip_x") * cm;
0171     const double strip_y = params->get_double_param("strip_y") * cm;
0172     const int nstrips_phi_sensor = params->get_int_param("nstrips_phi_sensor");
0173     const double sensor_offset_y = params->get_double_param("sensor_offset_y") * cm;
0174     const double hdi_y = params->get_double_param("hdi_y") * cm;
0175     double hdi_kapton_x = params->get_double_param("hdi_kapton_x") * cm;
0176     double hdi_copper_x = params->get_double_param("hdi_copper_x") * cm;
0177     double fphx_x = params->get_double_param("fphx_x") * cm;
0178     double fphx_y = params->get_double_param("fphx_y") * cm;
0179     double fphx_z = params->get_double_param("fphx_z") * cm;
0180     double fphx_offset_z = params->get_double_param("fphx_offset_z") * cm;
0181 
0182     double si_glue_x = params->get_double_param("si_glue_x") * cm;
0183     double fphx_glue_x = params->get_double_param("fphx_glue_x") * cm;
0184     double halfladder_inside_z = params->get_double_param("halfladder_inside_z") * cm;
0185 
0186     if (Verbosity() > 0)
0187     {
0188       std::cout << "Constructing Intt layer: " << std::endl;
0189       std::cout << "  layer " << inttlayer << " laddertype " << laddertype << " nladders_layer " << nladders_layer << " sensor_radius " << m_SensorRadius[inttlayer] << " offsetphi " << offsetphi
0190                 << " rad "
0191                 << " offsetphi " << offsetphi * rad / deg << " deg " << std::endl;
0192     }
0193     // We loop over inner, then outer (wrt the beam-axis), sensors, where  itype specifies the inner or outer sensor
0194     // The rest of this loop will construct and put in place a section of a ladder corresponding to the Z range of this sensor only
0195     for (int itype = 0; itype < 2; ++itype)
0196     {
0197       double strip_z;
0198       int nstrips_z_sensor;
0199       switch (itype)
0200       {
0201       case 0:
0202         strip_z = params->get_double_param("strip_z_0") * cm;
0203         nstrips_z_sensor = params->get_int_param("nstrips_z_sensor_0");
0204         break;
0205       case 1:
0206         strip_z = params->get_double_param("strip_z_1") * cm;
0207         nstrips_z_sensor = params->get_int_param("nstrips_z_sensor_1");
0208         break;
0209       default:
0210         std::cout << "invalid itype " << itype << std::endl;
0211         exit(1);
0212       }
0213 
0214       // ----- Step 1 ------------------------------------------------------
0215       // We make the volumes for Si-sensor, FPHX, HDI, and stave components
0216       // We add them to the ladder later
0217       //====================================================================
0218 
0219       // Create Si-sensor active volume
0220       const double siactive_x = strip_x;
0221       const double siactive_y = strip_y * nstrips_phi_sensor;
0222       const double siactive_z = strip_z * nstrips_z_sensor;
0223       G4VSolid *siactive_box = new G4Box(std::format("siactive_box_{}_{}", inttlayer, itype), siactive_x / 2, siactive_y / 2., siactive_z / 2.);
0224       G4LogicalVolume *siactive_volume = new G4LogicalVolume(siactive_box, GetDetectorMaterial("G4_Si"), std::format("siactive_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0225       if ((m_IsActiveMap.find(inttlayer))->second > 0)
0226       {
0227         m_ActiveLogVols.insert(siactive_volume);
0228       }
0229       m_DisplayAction->AddVolume(siactive_volume, "SiActive");
0230       // We do not subdivide the sensor in G4. We will assign hits to strips in the stepping action, using the geometry object
0231 
0232       // Si-sensor full (active+inactive) area
0233       const double sifull_x = siactive_x;
0234       const double sifull_y = siactive_y + 2.0 * params->get_double_param("sensor_edge_phi") * cm;
0235       const double sifull_z = siactive_z + 2.0 * params->get_double_param("sensor_edge_z") * cm;
0236       G4VSolid *sifull_box = new G4Box(std::format("sifull_box_{}_{}", inttlayer, itype), sifull_x / 2., sifull_y / 2.0, sifull_z / 2.0);
0237 
0238       // Si-sensor inactive area
0239       G4VSolid *siinactive_box = new G4SubtractionSolid(std::format("siinactive_box_{}_{}", inttlayer, itype), sifull_box, siactive_box, nullptr, G4ThreeVector(0, 0, 0));
0240       G4LogicalVolume *siinactive_volume = new G4LogicalVolume(siinactive_box, GetDetectorMaterial("G4_Si"), std::format("siinactive_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0241 
0242       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0243       {
0244         m_PassiveVolumeTuple.insert(std::make_pair(siinactive_volume, std::make_tuple(inttlayer, PHG4InttDefs::SI_INACTIVE)));
0245       }
0246       m_DisplayAction->AddVolume(siinactive_volume, "SiInActive");
0247 
0248       // Glue for Si-sensor full area
0249       G4VSolid *si_glue_box = new G4Box(std::format("si_glue_box_{}_{}", inttlayer, itype), si_glue_x / 2., sifull_y / 2.0, sifull_z / 2.0);
0250 
0251       G4LogicalVolume *si_glue_volume = new G4LogicalVolume(si_glue_box, GetDetectorMaterial("SilverEpoxyGlue_INTT"), std::format("si_glue_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0252 
0253       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0254       {
0255         m_PassiveVolumeTuple.insert(std::make_pair(siinactive_volume, std::make_tuple(inttlayer, PHG4InttDefs::SI_GLUE)));
0256       }
0257       m_DisplayAction->AddVolume(si_glue_volume, "SiGlue");
0258 
0259       // Make the HDI Kapton and copper volumes
0260       // This makes HDI volumes that matche this sensor in Z length
0261       const double hdi_z = sifull_z + params->get_double_param("hdi_edge_z") * cm;
0262       hdi_z_arr[inttlayer][itype] = hdi_z;
0263       G4VSolid *hdi_kapton_box = new G4Box(std::format("hdi_kapton_box_{}_{}", inttlayer, itype), hdi_kapton_x / 2., hdi_y / 2., hdi_z / 2.0);
0264       G4LogicalVolume *hdi_kapton_volume = new G4LogicalVolume(hdi_kapton_box, GetDetectorMaterial("G4_KAPTON"), std::format("hdi_kapton_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0265 
0266       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0267       {
0268         m_PassiveVolumeTuple.insert(std::make_pair(hdi_kapton_volume, std::make_tuple(inttlayer, PHG4InttDefs::HDI_KAPTON)));
0269       }
0270       G4VSolid *hdi_copper_box = new G4Box(std::format("hdi_copper_box_{}_{}", inttlayer, itype), hdi_copper_x / 2., hdi_y / 2., hdi_z / 2.0);
0271       G4LogicalVolume *hdi_copper_volume = new G4LogicalVolume(hdi_copper_box, GetDetectorMaterial("G4_Cu"), std::format("hdi_copper_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0272       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0273       {
0274         m_PassiveVolumeTuple.insert(std::make_pair(hdi_copper_volume, std::make_tuple(inttlayer, PHG4InttDefs::HDI_COPPER)));
0275       }
0276       m_DisplayAction->AddVolume(hdi_kapton_volume, "HdiKapton");
0277       m_DisplayAction->AddVolume(hdi_copper_volume, "HdiCopper");
0278 
0279       // This is the part of the HDI that extends beyond the sensor inside the endcap ring
0280       const double hdiext_z = (itype == 0) ? 0.000001 : halfladder_inside_z - hdi_z_arr[inttlayer][0] - hdi_z;  // need to assign nonzero value for itype=0
0281       // const double hdiext_z = (itype == 0) ? 0.000001 : 0.8;
0282       G4VSolid *hdiext_kapton_box = new G4Box(std::format("hdiext_kapton_box_{}_{}", inttlayer, itype), hdi_kapton_x / 2., hdi_y / 2., hdiext_z / 2.0);
0283       G4LogicalVolume *hdiext_kapton_volume = new G4LogicalVolume(hdiext_kapton_box, GetDetectorMaterial("G4_KAPTON"),  // was "FPC"
0284                                                                   std::format("hdiext_kapton_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0285       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0286       {
0287         m_PassiveVolumeTuple.insert(std::make_pair(hdiext_kapton_volume, std::make_tuple(inttlayer, PHG4InttDefs::HDIEXT_KAPTON)));
0288       }
0289       G4VSolid *hdiext_copper_box = new G4Box(std::format("hdiext_copper_box_{}_{}", inttlayer, itype), hdi_copper_x / 2., hdi_y / 2., hdiext_z / 2.0);
0290       G4LogicalVolume *hdiext_copper_volume = new G4LogicalVolume(hdiext_copper_box, GetDetectorMaterial("G4_Cu"), std::format("hdiext_copper_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0291       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0292       {
0293         m_PassiveVolumeTuple.insert(std::make_pair(hdiext_copper_volume, std::make_tuple(inttlayer, PHG4InttDefs::HDIEXT_COPPER)));
0294       }
0295       m_DisplayAction->AddVolume(hdiext_kapton_volume, "HdiKapton");
0296       m_DisplayAction->AddVolume(hdiext_copper_volume, "HdiCopper");
0297 
0298       // FPHX
0299       G4VSolid *fphx_box = new G4Box(std::format("fphx_box_{}_{}", inttlayer, itype), fphx_x / 2., fphx_y / 2., fphx_z / 2.);
0300       G4LogicalVolume *fphx_volume = new G4LogicalVolume(fphx_box, GetDetectorMaterial("G4_Si"), std::format("fphx_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0301       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0302       {
0303         m_PassiveVolumeTuple.insert(std::make_pair(fphx_volume, std::make_tuple(inttlayer, PHG4InttDefs::FPHX)));
0304       }
0305       m_DisplayAction->AddVolume(fphx_volume, "FPHX");
0306 
0307       const double gap_sensor_fphx = params->get_double_param("gap_sensor_fphx") * cm;
0308 
0309       //  FPHX Container
0310       // make a container for the FPHX chips needed for this sensor, and  then place them in the container
0311       G4VSolid *fphxcontainer_box = new G4Box(std::format("fphxcontainer_box_{}_{}", inttlayer, itype), fphx_x / 2., fphx_y / 2., hdi_z / 2.);
0312       G4LogicalVolume *fphxcontainer_volume =
0313           new G4LogicalVolume(fphxcontainer_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")), std::format("fphxcontainer_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0314       m_DisplayAction->AddVolume(fphxcontainer_volume, "FPHXContainer");
0315 
0316       // Install multiple FPHX volumes in the FPHX container volume
0317       // one FPHX chip per cell - each cell is 128 channels
0318       const double fphx_offsetx = 0.;
0319       const double fphx_offsety = 0.;
0320       int ncopy;
0321       double offsetz;
0322       double cell_length_z;
0323 
0324       if (laddertype == PHG4InttDefs::SEGMENTATION_Z)  // vertical strips
0325       {
0326         // For laddertype 0, we have 5 cells per sensor, but the strips are vertical, so we have to treat it specially
0327         ncopy = nstrips_z_sensor / 128.0;
0328       }
0329       else if (laddertype == PHG4InttDefs::SEGMENTATION_PHI)
0330       {
0331         ncopy = nstrips_z_sensor;
0332       }
0333       else
0334       {
0335         std::cout << PHWHERE << "invalid laddertype " << laddertype << std::endl;
0336         gSystem->Exit(1);
0337         // this is just to make the optimizer happy which otherwise complains about possibly
0338         // uninitialized variables. It doesn't know gSystem->Exit(1) quits,
0339         // this exit here terminates the program for it
0340         exit(1);
0341       }
0342       cell_length_z = strip_z * nstrips_z_sensor / ncopy;
0343       // NOLINTNEXTLINE(bugprone-integer-division)
0344       offsetz = (ncopy % 2 == 0) ? -2. * cell_length_z / 2. * double(ncopy / 2) + cell_length_z / 2. + fphx_offset_z : -2. * cell_length_z / 2. * double(ncopy / 2) + fphx_offset_z;
0345 
0346       G4VPVParameterisation *fphxparam = new PHG4InttFPHXParameterisation(fphx_offsetx, +fphx_offsety, offsetz, 2. * cell_length_z / 2., ncopy);
0347       new G4PVParameterised(std::format("fphxcontainer_{}_{}", inttlayer, itype), fphx_volume, fphxcontainer_volume, kZAxis, ncopy, fphxparam, OverlapCheck());
0348 
0349       // Glue for FPHX, silver powder epoxy, impletemented in the same way as FPHX
0350       G4VSolid *fphx_glue_box = new G4Box(std::format("fphx_glue_box_{}_{}", inttlayer, itype), fphx_glue_x / 2., fphx_y / 2., fphx_z / 2.);
0351 
0352       G4LogicalVolume *fphx_glue_volume =
0353           new G4LogicalVolume(fphx_glue_box, GetDetectorMaterial("SilverEpoxyGlue_INTT"), std::format("fphx_glue_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0354       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0355       {
0356         m_PassiveVolumeTuple.insert(std::make_pair(fphx_glue_volume, std::make_tuple(inttlayer, PHG4InttDefs::FPHX_GLUE)));
0357       }
0358       m_DisplayAction->AddVolume(fphx_glue_volume, "FPHXGlue");
0359 
0360       //  Glue of FPHX Container
0361       // make a container for the glue of FPHX chips, and then place them in the container
0362       G4VSolid *fphx_gluecontainer_box = new G4Box(std::format("fphx_gluecontainer_box_{}_{}", inttlayer, itype), fphx_glue_x / 2., fphx_y / 2., hdi_z / 2.);
0363       G4LogicalVolume *fphx_gluecontainer_volume = new G4LogicalVolume(fphx_gluecontainer_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")),
0364                                                                        std::format("fphx_gluecontainer_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0365 
0366       // Parameters for FPHX glue for G4VPVParameterisation are the same as FPGX's, so reuse them!
0367       G4VPVParameterisation *fphx_glueparam = new PHG4InttFPHXParameterisation(fphx_offsetx, +fphx_offsety, offsetz, 2. * cell_length_z / 2., ncopy);
0368 
0369       new G4PVParameterised(std::format("glue_fphxcontainer_{}_{}", inttlayer, itype), fphx_glue_volume, fphx_gluecontainer_volume, kZAxis, ncopy, fphx_glueparam, OverlapCheck());
0370       m_DisplayAction->AddVolume(fphx_gluecontainer_volume, "FPHXGlueContainer");
0371 
0372       double stave_x = 0.;
0373       double stave_y = 0.;
0374       G4LogicalVolume *stave_volume = nullptr;
0375       G4LogicalVolume *staveext_volume = nullptr;
0376 
0377       // Carbon stave. This consists of the formed sheet, cooling water, the water tube, glue for the tube,
0378       // rohacell foam to fill space around the tube, and the flat CFRP sheet, which completes the outer shell surrounds
0379 
0380       // Rohacel foam and cooling water pipe inside. Formed from straight sections and sections of a tube of
0381       // radius 3.1905 mm. All have wall thickness of 0.1905 mm.
0382       const double stave_thickness = params->get_double_param("stave_straight_cooler_x") * cm;  // stave thickness
0383       const double Rcmin = 0.30 * cm;                                                           // inner radius of curved section, same at both ends
0384       const double Rcmax = Rcmin + stave_thickness;                                             // outer radius of curved section, same at both ends
0385       double Rcavge = (Rcmax + Rcmin) / 2.0;
0386       double dphi_c = 23.19859051 * M_PI / 180.;  // phi of the curved section
0387       const double stave_z = hdi_z;
0388 
0389       // Make CFC structure
0390       //// Make curved sections
0391       const double phic_begin[4] = {M_PI - dphi_c, -dphi_c, 0.0, M_PI};
0392       const double dphic[4] = {dphi_c, dphi_c, dphi_c, dphi_c};
0393 
0394       G4Tubs *stave_curve_cons[4];
0395       G4Tubs *stave_curve_ext_cons[4];
0396       G4LogicalVolume *stave_curve_volume[4];
0397       G4LogicalVolume *stave_curve_ext_volume[4];
0398 
0399       for (int i = 0; i < 4; i++)
0400       {
0401         stave_curve_cons[i] = new G4Tubs(std::format("stave_curve_cons_{}_{}_{}", inttlayer, itype, i), Rcmin, Rcmax, stave_z / 2., phic_begin[i], dphic[i]);
0402         stave_curve_volume[i] =
0403             new G4LogicalVolume(stave_curve_cons[i], GetDetectorMaterial("CFRP_INTT"), std::format("stave_curve_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0404         if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0405         {
0406           m_PassiveVolumeTuple.insert(std::make_pair(stave_curve_volume[i], std::make_tuple(inttlayer, PHG4InttDefs::STAVE_CURVE)));
0407         }
0408         stave_curve_ext_cons[i] = new G4Tubs(std::format("stave_curve_ext_cons_{}_{}_{}", inttlayer, itype, i), Rcmin, Rcmax, hdiext_z / 2., phic_begin[i], dphic[i]);
0409         stave_curve_ext_volume[i] =
0410             new G4LogicalVolume(stave_curve_ext_cons[i], GetDetectorMaterial("CFRP_INTT"), std::format("stave_curve_ext_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0411         if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0412         {
0413           m_PassiveVolumeTuple.insert(std::make_pair(stave_curve_ext_volume[i], std::make_tuple(inttlayer, PHG4InttDefs::STAVEEXT_CURVE)));
0414         }
0415         m_DisplayAction->AddVolume(stave_curve_volume[i], "StaveCurve");
0416         m_DisplayAction->AddVolume(stave_curve_ext_volume[i], "StaveCurve");
0417       }
0418 
0419       // we will need the length in y of the curved section as it is installed in the stave box
0420       double curve_length_y = Rcavge * sin(dphi_c);
0421 
0422       // Make several straight sections for use in making the stave
0423       double stave_straight_outer_y = params->get_double_param("stave_straight_outer_y") * cm;
0424       double stave_straight_cooler_y = params->get_double_param("stave_straight_cooler_y") * cm;
0425       double rohacell_straight_y = params->get_double_param("stave_straight_rohacell_y") * cm;
0426 
0427       // Outer straight sections of stave
0428       G4VSolid *stave_straight_outer_box =
0429           new G4Box(std::format("stave_straight_outer_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_straight_outer_y / 2., stave_z / 2.);
0430       G4LogicalVolume *stave_straight_outer_volume =
0431           new G4LogicalVolume(stave_straight_outer_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_straight_outer_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0432       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0433       {
0434         m_PassiveVolumeTuple.insert(std::make_pair(stave_straight_outer_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVE_STRAIGHT_OUTER)));
0435       }
0436       G4VSolid *stave_straight_outer_ext_box =
0437           new G4Box(std::format("stave_straight_outer_ext_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_straight_outer_y / 2., hdiext_z / 2.);
0438       G4LogicalVolume *stave_straight_outer_ext_volume =
0439           new G4LogicalVolume(stave_straight_outer_ext_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_straight_outer_ext_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0440       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0441       {
0442         m_PassiveVolumeTuple.insert(std::make_pair(stave_straight_outer_ext_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVEEXT_STRAIGHT_OUTER)));
0443       }
0444 
0445       // Top surface of stave
0446       G4VSolid *stave_straight_cooler_box =
0447           new G4Box(std::format("stave_straight_cooler_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_straight_cooler_y / 2., stave_z / 2.);
0448       G4LogicalVolume *stave_straight_cooler_volume =
0449           new G4LogicalVolume(stave_straight_cooler_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_straight_cooler_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0450       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0451       {
0452         m_PassiveVolumeTuple.insert(std::make_pair(stave_straight_cooler_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVE_STRAIGHT_COOLER)));
0453       }
0454       G4VSolid *stave_straight_cooler_ext_box =
0455           new G4Box(std::format("stave_straight_cooler_ext_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_straight_cooler_y / 2., hdiext_z / 2.);
0456       G4LogicalVolume *stave_straight_cooler_ext_volume =
0457           new G4LogicalVolume(stave_straight_cooler_ext_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_straight_cooler_ext_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0458       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0459       {
0460         m_PassiveVolumeTuple.insert(std::make_pair(stave_straight_cooler_ext_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVEEXT_STRAIGHT_COOLER)));
0461       }
0462 
0463       // Slant straight sections of stave
0464       double stave_slant_cooler_y = params->get_double_param("stave_slant_cooler_y") * cm;
0465       G4VSolid *stave_slant_cooler_box = new G4Box(std::format("stave_slant_cooler_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_slant_cooler_y / 2., stave_z / 2.);
0466       G4LogicalVolume *stave_slant_cooler_volume =
0467           new G4LogicalVolume(stave_slant_cooler_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_slant_cooler_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0468       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0469       {
0470         m_PassiveVolumeTuple.insert(std::make_pair(stave_slant_cooler_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVE_STRAIGHT_COOLER)));
0471       }
0472       G4VSolid *stave_slant_cooler_ext_box =
0473           new G4Box(std::format("stave_lant_cooler_ext_box_{}_{}", inttlayer, itype), stave_thickness / 2., stave_slant_cooler_y / 2., hdiext_z / 2.);
0474       G4LogicalVolume *stave_slant_cooler_ext_volume =
0475           new G4LogicalVolume(stave_slant_cooler_ext_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_slant_cooler_ext_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0476       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0477       {
0478         m_PassiveVolumeTuple.insert(std::make_pair(stave_slant_cooler_ext_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVEEXT_STRAIGHT_COOLER)));
0479       }
0480 
0481       // Flat CFRP sheet on the bottom of the stave structure. It was introduced instead of PGS
0482       G4VSolid *stave_bottom_cooler_box = new G4Box(std::format("stave_bottom_cooler_box_{}_{}", inttlayer, itype), stave_thickness / 2., hdi_y / 2., stave_z / 2.);
0483 
0484       G4LogicalVolume *stave_bottom_cooler_volume =
0485           new G4LogicalVolume(stave_bottom_cooler_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_bottom_cooler_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0486       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0487       {
0488         m_PassiveVolumeTuple.insert(std::make_pair(stave_bottom_cooler_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVE_BOTTOM_COOLER)));  // should be changed soon
0489       }
0490 
0491       G4VSolid *stave_bottom_cooler_ext_box = new G4Box(std::format("stave_bottom_cooler_ext_box_{}_{}", inttlayer, itype), stave_thickness / 2., hdi_y / 2., hdiext_z / 2.);
0492       G4LogicalVolume *stave_bottom_cooler_ext_volume =
0493           new G4LogicalVolume(stave_bottom_cooler_ext_box, GetDetectorMaterial("CFRP_INTT"), std::format("stave_bottom_cooler_ext_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0494       if ((m_IsAbsorberActiveMap.find(inttlayer))->second > 0)
0495       {
0496         m_PassiveVolumeTuple.insert(std::make_pair(stave_bottom_cooler_ext_volume, std::make_tuple(inttlayer, PHG4InttDefs::STAVEEXT_BOTTOM_COOLER)));
0497       }
0498 
0499       m_DisplayAction->AddVolume(stave_straight_cooler_volume, "StaveCooler");
0500       m_DisplayAction->AddVolume(stave_straight_cooler_ext_volume, "StaveCooler");
0501       m_DisplayAction->AddVolume(stave_straight_outer_volume, "StaveStraightOuter");
0502       m_DisplayAction->AddVolume(stave_straight_outer_ext_volume, "StaveStraightOuter");
0503       m_DisplayAction->AddVolume(stave_slant_cooler_volume, "StaveCooler");
0504       m_DisplayAction->AddVolume(stave_slant_cooler_ext_volume, "StaveCooler");
0505       m_DisplayAction->AddVolume(stave_bottom_cooler_volume, "StaveCooler");
0506       m_DisplayAction->AddVolume(stave_bottom_cooler_ext_volume, "StaveCooler");
0507 
0508       // cooling pipe + water inside + glue outside
0509       const double Rpmin = 0.10 * cm;  // inner radius of cooling pipe section, same at both ends
0510       const double Rpmax = 0.15 * cm;  // outer radius of cooling pipe section, same at both ends
0511       G4VSolid *stave_glue_box = new G4Box(std::format("stave_glue_box_{}_{}", inttlayer, itype), 3. / 2, 3. / 2., stave_z / 2.);
0512       G4LogicalVolume *stave_glue_volume = new G4LogicalVolume(stave_glue_box, GetDetectorMaterial("Epoxy"), std::format("stave_glue_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0513       G4VSolid *staveext_glue_box = new G4Box(std::format("staveext_glue_box_{}_{}", inttlayer, itype), 3. / 2., 3. / 2., hdiext_z / 2.);
0514       G4LogicalVolume *staveext_glue_volume =
0515           new G4LogicalVolume(staveext_glue_box, GetDetectorMaterial("Epoxy"), std::format("staveext_glue_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0516 
0517       m_DisplayAction->AddVolume(stave_glue_volume, "StaveGlueBox");
0518       m_DisplayAction->AddVolume(staveext_glue_volume, "StaveGlueBox");
0519 
0520       G4VSolid *stave_pipe_cons = new G4Tubs(std::format("stave_pipe_cons_{}_{}", inttlayer, itype), Rpmin, Rpmax, stave_z / 2., 0, 2.0 * M_PI);
0521       G4LogicalVolume *stave_pipe_volume = new G4LogicalVolume(stave_pipe_cons, GetDetectorMaterial("CFRP_INTT"), std::format("stave_pipe_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0522 
0523       G4VSolid *staveext_pipe_cons = new G4Tubs(std::format("staveext_pipe_cons_{}_{}", inttlayer, itype), Rpmin, Rpmax, hdiext_z / 2., 0, 2.0 * M_PI);
0524       G4LogicalVolume *staveext_pipe_volume =
0525           new G4LogicalVolume(staveext_pipe_cons, GetDetectorMaterial("CFRP_INTT"), std::format("staveext_pipe_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0526 
0527       m_DisplayAction->AddVolume(stave_pipe_volume, "StavePipe");
0528       m_DisplayAction->AddVolume(staveext_pipe_volume, "StavePipe");
0529 
0530       G4VSolid *stave_water_cons = new G4Tubs(std::format("stave_water_cons_{}_{}", inttlayer, itype), 0., Rpmin, stave_z / 2., 0, 2.0 * M_PI);
0531       G4LogicalVolume *stave_water_volume =
0532           new G4LogicalVolume(stave_water_cons, GetDetectorMaterial("G4_WATER"), std::format("stave_water_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0533 
0534       G4VSolid *staveext_water_cons = new G4Tubs(std::format("staveext_water_cons_{}_{}", inttlayer, itype), 0., Rpmin, hdiext_z / 2., 0, 2.0 * M_PI);
0535       G4LogicalVolume *staveext_water_volume =
0536           new G4LogicalVolume(staveext_water_cons, GetDetectorMaterial("G4_WATER"), std::format("staveext_water_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0537 
0538       m_DisplayAction->AddVolume(stave_water_volume, "StaveWater");
0539       m_DisplayAction->AddVolume(staveext_water_volume, "StaveWater");
0540 
0541       // rohacell foam
0542       // straight boxes
0543       G4VSolid *rohacell_straight_cons = new G4Box(std::format("rohacell_straight_cons_{}_{}", inttlayer, itype), 3. / 2, rohacell_straight_y / 2., stave_z / 2.);
0544       G4LogicalVolume *rohacell_straight_volume =
0545           new G4LogicalVolume(rohacell_straight_cons, GetDetectorMaterial("ROHACELL_FOAM_51"), std::format("rohacell_straight_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0546 
0547       G4VSolid *rohacellext_straight_cons = new G4Box(std::format("rohacellext_straight_cons_{}_{}", inttlayer, itype), 3. / 2, rohacell_straight_y / 2., hdiext_z / 2.);
0548       G4LogicalVolume *rohacellext_straight_volume =
0549           new G4LogicalVolume(rohacellext_straight_cons, GetDetectorMaterial("ROHACELL_FOAM_51"), std::format("rohacellext_straight_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0550 
0551       // make curved sections for rohacell foam
0552       const double rh_phic_begin[2] = {-dphi_c, 0.0};
0553       const double rh_dphic[2] = {dphi_c, dphi_c};
0554       G4Tubs *rohacell_curve_cons[2];
0555       G4LogicalVolume *rohacell_curve_volume[2];
0556       G4Tubs *rohacellext_curve_cons[2];
0557       G4LogicalVolume *rohacellext_curve_volume[2];
0558       for (int i = 0; i < 2; i++)
0559       {
0560         rohacell_curve_cons[i] = new G4Tubs(std::format("rohacell_curve_cons_{}_{}_{}", inttlayer, itype, i), 0., Rcmin, stave_z / 2., rh_phic_begin[i], rh_dphic[i]);
0561         rohacell_curve_volume[i] =
0562             new G4LogicalVolume(rohacell_curve_cons[i], GetDetectorMaterial("ROHACELL_FOAM_51"), std::format("rohacell_curve_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0563         rohacellext_curve_cons[i] = new G4Tubs(std::format("rohacellext_curve_cons_{}_{}_{}", inttlayer, itype, i), 0., Rcmin, hdiext_z / 2., rh_phic_begin[i], rh_dphic[i]);
0564         rohacellext_curve_volume[i] = new G4LogicalVolume(rohacellext_curve_cons[i], GetDetectorMaterial("ROHACELL_FOAM_51"),
0565                                                           std::format("rohacellext_curve_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0566       }
0567 
0568       // make trapezoidal sections for rohacell foam
0569       G4GenericTrap *rohacell_trap_cons[2];
0570       G4LogicalVolume *rohacell_trap_volume[2];
0571       G4GenericTrap *rohacellext_trap_cons[2];
0572       G4LogicalVolume *rohacellext_trap_volume[2];
0573       for (int i = 0; i < 2; i++)
0574       {
0575         double shift = 1.e-5;  // To mitigate fm order level overlaps reported by GEANT4...
0576         std::vector<G4TwoVector> rohatrap(8);
0577         if (i == 0)
0578         {
0579           rohatrap[0] = G4TwoVector(0. * cm, 0. * cm);
0580           rohatrap[1] = G4TwoVector(Rcmin * cos(dphi_c) - shift, -Rcmin * sin(dphi_c));
0581           rohatrap[2] = G4TwoVector(Rcmin * (1. - cos(dphi_c)) - shift, -stave_slant_cooler_y * cos(dphi_c) - Rcmin * sin(dphi_c));
0582           rohatrap[3] = G4TwoVector(0. * cm, -stave_slant_cooler_y * cos(dphi_c) - Rcmin * sin(dphi_c));
0583         }
0584         else
0585         {
0586           rohatrap[0] = G4TwoVector(0. * cm, +stave_slant_cooler_y * cos(dphi_c) + Rcmin * sin(dphi_c));
0587           rohatrap[1] = G4TwoVector(Rcmax * (1. - cos(dphi_c)) - shift, +stave_slant_cooler_y * cos(dphi_c) + Rcmin * sin(dphi_c));
0588           rohatrap[2] = G4TwoVector(Rcmin * cos(dphi_c) - shift, +Rcmin * sin(dphi_c));
0589           rohatrap[3] = G4TwoVector(0. * cm, 0. * cm);
0590         }
0591         rohatrap[4] = rohatrap[0];
0592         rohatrap[5] = rohatrap[1];
0593         rohatrap[6] = rohatrap[2];
0594         rohatrap[7] = rohatrap[3];
0595 
0596         rohacell_trap_cons[i] = new G4GenericTrap(std::format("rohacell_trap_cons_{}_{}_{}", inttlayer, itype, i), stave_z / 2., rohatrap);
0597         rohacell_trap_volume[i] =
0598             new G4LogicalVolume(rohacell_trap_cons[i], GetDetectorMaterial("ROHACELL_FOAM_51"), std::format("rohacell_trap_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0599 
0600         rohacellext_trap_cons[i] = new G4GenericTrap(std::format("rohacellext_trap_cons_{}_{}_{}", inttlayer, itype, i), hdiext_z / 2., rohatrap);
0601         rohacellext_trap_volume[i] =
0602             new G4LogicalVolume(rohacellext_trap_cons[i], GetDetectorMaterial("ROHACELL_FOAM_51"), std::format("rohacellext_trap_volume_{}_{}_{}", inttlayer, itype, i), nullptr, nullptr, nullptr);
0603       }
0604 
0605       m_DisplayAction->AddVolume(rohacell_straight_volume, "RohaCell");
0606       m_DisplayAction->AddVolume(rohacellext_straight_volume, "RohaCell");
0607       for (int i = 0; i < 2; i++)
0608       {
0609         m_DisplayAction->AddVolume(rohacell_curve_volume[i], "RohaCell");
0610         m_DisplayAction->AddVolume(rohacellext_curve_volume[i], "RohaCell");
0611         m_DisplayAction->AddVolume(rohacell_trap_volume[i], "RohaCell");
0612         m_DisplayAction->AddVolume(rohacellext_trap_volume[i], "RohaCell");
0613       }
0614 
0615       // Now we combine the elements of a stave defined above into a stave
0616       // Create a stave volume to install the stave sections into. The volume has to be big enouigh to contain the cooling tube
0617       double cooler_gap_x = 0.3 * cm;                      // id of cooling tube in cm
0618       double cooler_wall = stave_thickness;                // outer wall thickness of cooling tube
0619       double cooler_x = cooler_gap_x + 2.0 * cooler_wall;  // thickness of the formed sheet, the flat sheet, and the gap b/w the sheets
0620       stave_x = cooler_x;
0621       stave_y = hdi_y;
0622 
0623       // Make stave volume. Drop two corners in positive x to prevent ladder_volume overlapping
0624       // with neighbouring ladders because of small clearance in the latest configuration
0625       G4RotationMatrix *stv_rot_pos = new G4RotationMatrix();
0626       stv_rot_pos->rotateZ(-15. * M_PI / 180.);
0627       G4ThreeVector stvTranspos(stave_x / 2., stave_y / 2., 0.);
0628 
0629       G4RotationMatrix *stv_rot_neg = new G4RotationMatrix();
0630       stv_rot_neg->rotateZ(+15. * M_PI / 180.);
0631       G4ThreeVector stvTransneg(stave_x / 2., -stave_y / 2., 0.);
0632 
0633       G4VSolid *stave_basebox = new G4Box(std::format("stave_basebox_{}_{}", inttlayer, itype), stave_x / 2., stave_y / 2., stave_z / 2.);
0634       G4VSolid *stave_subtbox =
0635           new G4Box(std::format("stave_subtbox_{}_{}", inttlayer, itype), stave_x / 1.5, stave_y / 1.5, stave_z / 1.);  // has to be longer in z to avoid coincident surface
0636 
0637       G4VSolid *stave_box1 = new G4SubtractionSolid(std::format("stave_box1_{}_{}", inttlayer, itype), stave_basebox, stave_subtbox, stv_rot_pos, stvTranspos);
0638 
0639       G4VSolid *stave_box = new G4SubtractionSolid(std::format("stave_box_{}_{}", inttlayer, itype), stave_box1, stave_subtbox, stv_rot_neg, stvTransneg);
0640 
0641       stave_volume = new G4LogicalVolume(stave_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")), std::format("stave_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0642 
0643       G4VSolid *staveext_basebox = new G4Box(std::format("staveext_basebox_{}_{}", inttlayer, itype), stave_x / 2., stave_y / 2., hdiext_z / 2.);
0644       G4VSolid *staveext_subtbox =
0645           new G4Box(std::format("staveext_subtbox_{}_{}", inttlayer, itype), stave_x / 1.5, stave_y / 1.5, hdiext_z / 1.);  // has to be longer in z to avoid coincident surface
0646 
0647       G4VSolid *staveext_box1 = new G4SubtractionSolid(std::format("staveext_box1_{}_{}", inttlayer, itype), staveext_basebox, staveext_subtbox, stv_rot_pos, stvTranspos);
0648 
0649       G4VSolid *staveext_box = new G4SubtractionSolid(std::format("staveext_box_{}_{}", inttlayer, itype), staveext_box1, staveext_subtbox, stv_rot_neg, stvTransneg);
0650 
0651       staveext_volume = new G4LogicalVolume(staveext_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")), std::format("staveext_volume_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0652       // the rotation matrices are just used by G4VSolid, ownership is not taken over
0653       delete stv_rot_pos;
0654       delete stv_rot_neg;
0655 
0656       m_DisplayAction->AddVolume(stave_volume, "StaveBox");
0657       m_DisplayAction->AddVolume(staveext_volume, "StaveBox");
0658 
0659       // Assemble the elements into the stave volume and the stave extension volume
0660       // They are place relative to the center of the stave box. Thus the offset of the center of the segment is relative to the center of the satev box.
0661       // But we want the segment to be located relative to the lowest x limit of the stave box.
0662       if (laddertype == PHG4InttDefs::SEGMENTATION_Z)  // Obsolete!!
0663       {
0664         // only one cooling tube in laddertype 0
0665         // Place the straight sections. We add the middle, then above x axis, then below x axis
0666         double x_off_str[3] = {Rcavge - stave_x / 2., (Rcmax - Rcmin) / 2. - stave_x / 2., (Rcmax - Rcmin) / 2. - stave_x / 2.};
0667         double y_off_str[3] = {0.0, +stave_straight_cooler_y / 2. + 2. * curve_length_y + stave_straight_outer_y / 2.,
0668                                -stave_straight_cooler_y / 2. - 2. * curve_length_y - stave_straight_outer_y / 2.};
0669 
0670         for (int i = 0; i < 3; i++)
0671         {
0672           if (i == 0)
0673           {
0674             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_cooler_volume,
0675                               std::format("stave_straight_cooler_{}_{}_{}", i, inttlayer, itype), stave_volume, false, 0, OverlapCheck());
0676             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_cooler_ext_volume,
0677                               std::format("stave_straight_cooler_ext_{}_{}_{}", i, inttlayer, itype), staveext_volume, false, 0, OverlapCheck());
0678           }
0679           else
0680           {
0681             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_outer_volume,
0682                               std::format("stave_straight_outer_{}_{}_{}", i, inttlayer, itype), stave_volume, false, 0, OverlapCheck());
0683             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_outer_ext_volume,
0684                               std::format("stave_straight_outer_ext_{}_{}_{}", i, inttlayer, itype), staveext_volume, false, 0, OverlapCheck());
0685           }
0686         }
0687         // The cooler curved sections are made using 2 curved sections in a recurve on each side of the cooler straight section
0688         // The tube sections used here have the origin of their volume at their center of rotation. Rcavge
0689         //      Each curve section is moved to the center of the stave volume by a translation of +/- Rcavge
0690         //      Then it is moved to the outside or the inside of the stave volume by a translation of +/-  cooler_gap_x / 2.
0691         // we start at lowest y and work up in y
0692 
0693         double x_off_cooler[4] = {Rcavge - cooler_gap_x / 2., -Rcavge + cooler_gap_x / 2., -Rcavge + cooler_gap_x / 2., Rcavge - cooler_gap_x / 2.};
0694         double y_off_cooler[4] = {-stave_straight_cooler_y / 2. - 2. * curve_length_y, -stave_straight_cooler_y / 2., +stave_straight_cooler_y / 2.,
0695                                   +stave_straight_cooler_y / 2. + 2. * curve_length_y};
0696 
0697         for (int i = 0; i < 4; i++)
0698         {
0699           new G4PVPlacement(nullptr, G4ThreeVector(x_off_cooler[i], y_off_cooler[i], 0.0), stave_curve_volume[i], std::format("stave_curve_{}_{}_{}", inttlayer, itype, i),
0700                             stave_volume, false, 0, OverlapCheck());
0701           new G4PVPlacement(nullptr, G4ThreeVector(x_off_cooler[i], y_off_cooler[i], 0.0), stave_curve_ext_volume[i], std::format("stave_curve_ext_{}_{}_{}", inttlayer, itype, i),
0702                             staveext_volume, false, 0, OverlapCheck());
0703         }
0704       }
0705       else if (laddertype == PHG4InttDefs::SEGMENTATION_PHI)  // The type PHG4InttDefs::SEGMENTATION_PHI ladder
0706       {
0707         // First place the straight sections, do the extension at the same time
0708         // we alternate  positive and negative y values here
0709         double x_off_str[6] = {
0710             (Rcmax + Rcmin) / 2. - stave_x / 2. + stave_thickness,  // inner straight section
0711             (Rcmax + Rcmax) / 4. - stave_x / 2. + stave_thickness,  // slant section
0712             (Rcmax + Rcmax) / 4. - stave_x / 2. + stave_thickness,  // slant section
0713             (Rcmax - Rcmin) / 2. - stave_x / 2. + stave_thickness,  // outer straight section
0714             (Rcmax - Rcmin) / 2. - stave_x / 2. + stave_thickness,  // outer straight section
0715             (Rcmax - Rcmin) / 2. - stave_x / 2.                     // bottom section
0716         };
0717         double y_off_str[6] = {
0718             0.0,                                                                                                                     // inner straight section
0719             -stave_straight_cooler_y / 2. - 1. * curve_length_y - cos(dphi_c) * stave_slant_cooler_y / 2.,                           // slant section
0720             +stave_straight_cooler_y / 2. + 1. * curve_length_y + cos(dphi_c) * stave_slant_cooler_y / 2.,                           // slant section
0721             -stave_straight_cooler_y / 2. - 2. * curve_length_y - cos(dphi_c) * stave_slant_cooler_y - stave_straight_outer_y / 2.,  // outer straight section
0722             +stave_straight_cooler_y / 2. + 2. * curve_length_y + cos(dphi_c) * stave_slant_cooler_y + stave_straight_outer_y / 2.,  // outer straight section
0723             0.0
0724             // bottom straight section
0725         };
0726 
0727         for (int i = 0; i < 6; i++)
0728         {
0729           if (i == 0)  // inner straight section
0730           {
0731             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_cooler_volume,
0732                               std::format("stave_straight_cooler_{}_{}_{}", inttlayer, itype, i), stave_volume, false, 0, OverlapCheck());
0733             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_cooler_ext_volume,
0734                               std::format("stave_straight_cooler_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0735           }
0736           else if (i == 1 || i == 2)  // slant section
0737           {
0738             G4RotationMatrix rotation;
0739             if (i == 1)
0740             {
0741               rotation.rotateZ(-1. * dphi_c);
0742             }
0743             else if (i == 2)
0744             {
0745               rotation.rotateZ(dphi_c);
0746             }
0747             new G4PVPlacement(G4Transform3D(rotation, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0)), stave_slant_cooler_volume,
0748                               std::format("stave_slant_cooler_{}_{}_{}", inttlayer, itype, i), stave_volume, false, 0, OverlapCheck());
0749             new G4PVPlacement(G4Transform3D(rotation, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0)), stave_slant_cooler_ext_volume,
0750                               std::format("stave_slant_cooler_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0751           }
0752           else if (i == 3 || i == 4)  // outer straight section
0753           {
0754             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_outer_volume,
0755                               std::format("stave_straight_outer_{}_{}_{}", inttlayer, itype, i), stave_volume, false, 0, OverlapCheck());
0756             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_straight_outer_ext_volume,
0757                               std::format("stave_straight_outer_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0758           }
0759           else  // bottom
0760           {
0761             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_bottom_cooler_volume, std::format("stave_bottom_cooler_{}_{}_{}", inttlayer, itype, i),
0762                               stave_volume, false, 0, OverlapCheck());
0763             new G4PVPlacement(nullptr, G4ThreeVector(x_off_str[i], y_off_str[i], 0.0), stave_bottom_cooler_ext_volume,
0764                               std::format("stave_bottom_cooler_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0765           }
0766         }
0767 
0768         //// Place the curved sections
0769         //// here we do in order of increasing y
0770 
0771         double x_off_curve[4] = {// increasing in y
0772                                  +Rcavge - cooler_gap_x / 2. + stave_thickness / 2., -Rcavge + cooler_gap_x / 2. + stave_thickness / 2., -Rcavge + cooler_gap_x / 2. + stave_thickness / 2.,
0773                                  +Rcavge - cooler_gap_x / 2. + stave_thickness / 2.};
0774         double y_off_curve[4] = {// increasing in y
0775                                  -stave_straight_cooler_y / 2. - 2. * curve_length_y - cos(dphi_c) * stave_slant_cooler_y, -stave_straight_cooler_y / 2., +stave_straight_cooler_y / 2.,
0776                                  +stave_straight_cooler_y / 2. + 2. * curve_length_y + cos(dphi_c) * stave_slant_cooler_y};
0777 
0778         for (int i = 0; i < 4; i++)
0779         {
0780           new G4PVPlacement(nullptr, G4ThreeVector(x_off_curve[i], y_off_curve[i], 0.0), stave_curve_volume[i], std::format("stave_curve_{}_{}_{}", inttlayer, itype, i), stave_volume,
0781                             false, 0, OverlapCheck());
0782           new G4PVPlacement(nullptr, G4ThreeVector(x_off_curve[i], y_off_curve[i], 0.0), stave_curve_ext_volume[i], std::format("stave_curve_ext_{}_{}_{}", inttlayer, itype, i),
0783                             staveext_volume, false, 0, OverlapCheck());
0784         }
0785 
0786         // Place the rohacell foam
0787         // straight box section
0788         double x_off_roha_str[2] = {// increasing in y
0789                                     -cooler_wall / 2. + stave_thickness / 2., -cooler_wall / 2. + stave_thickness / 2.};
0790         double y_off_roha_str[2] = {// increasing in y
0791                                     -3. / 2. - rohacell_straight_y / 2., +3. / 2. + rohacell_straight_y / 2.};
0792 
0793         for (int i = 0; i < 2; i++)
0794         {
0795           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_str[i], y_off_roha_str[i], 0.0), rohacell_straight_volume,
0796                             std::format("rohacell_straight_{}_{}_{}", inttlayer, itype, i), stave_volume, false, 0, OverlapCheck());
0797           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_str[i], y_off_roha_str[i], 0.0), rohacellext_straight_volume,
0798                             std::format("rohacell_straight_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0799         }
0800 
0801         //// curve section
0802         double x_off_roha_curve[2] = {// increasing in y
0803                                       -Rcavge + cooler_gap_x / 2. + stave_thickness / 2., -Rcavge + cooler_gap_x / 2. + stave_thickness / 2.};
0804         double y_off_roha_curve[2] = {// increasing in y
0805                                       -3. / 2. - rohacell_straight_y, +3. / 2. + rohacell_straight_y};
0806 
0807         for (int i = 0; i < 2; i++)
0808         {
0809           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_curve[i], y_off_roha_curve[i], 0.0), rohacell_curve_volume[i],
0810                             std::format("rohacell_curve_{}_{}_{}", inttlayer, itype, i), stave_volume, false, 0, OverlapCheck());
0811           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_curve[i], y_off_roha_curve[i], 0.0), rohacellext_curve_volume[i],
0812                             std::format("rohacell_curve_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0813         }
0814 
0815         // trapezoidal section
0816         double x_off_roha_trap[2] = {// increasing in y
0817                                      -Rcmin - cooler_wall / 2. + cooler_gap_x / 2. + stave_thickness / 2., -Rcmin - cooler_wall / 2. + cooler_gap_x / 2. + stave_thickness / 2.};
0818         double y_off_roha_trap[2] = {// increasing in y
0819                                      -3. / 2. - rohacell_straight_y, +3. / 2. + rohacell_straight_y};
0820 
0821         for (int i = 0; i < 2; i++)
0822         {
0823           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_trap[i], y_off_roha_trap[i], 0.0), rohacell_trap_volume[i], std::format("rohacell_trap_{}_{}_{}", inttlayer, itype, i),
0824                             stave_volume, false, 0, OverlapCheck());
0825           new G4PVPlacement(nullptr, G4ThreeVector(x_off_roha_trap[i], y_off_roha_trap[i], 0.0), rohacellext_trap_volume[i],
0826                             std::format("rohacell_trap_ext_{}_{}_{}", inttlayer, itype, i), staveext_volume, false, 0, OverlapCheck());
0827         }
0828 
0829         // place glue box, cooling pipe and water inside
0830         new G4PVPlacement(nullptr, G4ThreeVector(0.0, 0.0, 0.0), stave_pipe_volume, std::format("stave_pipe_{}_{}", inttlayer, itype), stave_glue_volume, false, 0, OverlapCheck());
0831         new G4PVPlacement(nullptr, G4ThreeVector(0.0, 0.0, 0.0), staveext_pipe_volume, std::format("stave_pipe_ext_{}_{}", inttlayer, itype), staveext_glue_volume, false, 0,
0832                           OverlapCheck());
0833         new G4PVPlacement(nullptr, G4ThreeVector(0.0, 0.0, 0.0), stave_water_volume, std::format("stave_water_{}_{}", inttlayer, itype), stave_glue_volume, false, 0, OverlapCheck());
0834         new G4PVPlacement(nullptr, G4ThreeVector(0.0, 0.0, 0.0), staveext_water_volume, std::format("stave_water_ext_{}_{}", inttlayer, itype), staveext_glue_volume, false, 0,
0835                           OverlapCheck());
0836 
0837         // place of stave_glue_volume -cooler_wall / 2. + stave_thickness / 2. is actually 0. But I don't put 0 directry to make the origin of the value clear
0838         new G4PVPlacement(nullptr, G4ThreeVector(-cooler_wall / 2. + stave_thickness / 2., 0.0, 0.0), stave_glue_volume, std::format("stave_glue_{}_{}", inttlayer, itype), stave_volume,
0839                           false, 0, OverlapCheck());
0840         new G4PVPlacement(nullptr, G4ThreeVector(-cooler_wall / 2. + stave_thickness / 2., 0.0, 0.0), staveext_glue_volume, std::format("stave_glue_ext_{}_{}", inttlayer, itype),
0841                           staveext_volume, false, 0, OverlapCheck());
0842       }
0843       else
0844       {
0845         std::cout << PHWHERE << "invalid laddertype " << laddertype << std::endl;
0846         gSystem->Exit(1);
0847       }
0848 
0849       // ----- Step 2 ---------------------------------------------------
0850       // We place Si-sensor, FPHX, HDI, and stave in the ladder  volume.
0851       // ================================================================
0852 
0853       // Make the ladder volume first
0854       // We are still in the loop over inner or outer sensors. This is the ladder volume corresponding to this sensor.
0855       // But the thickness of the glue for FPHX is used since it's taller than the sensor in x.
0856 
0857       const double ladder_x = stave_x + hdi_kapton_x + hdi_copper_x + fphx_glue_x + fphx_x;
0858       double ladder_y = hdi_y;
0859 
0860       // Make ladder volume. Drop two corners in positive x as done for stave volume
0861       G4RotationMatrix *lad_box_rotpos = new G4RotationMatrix();
0862       lad_box_rotpos->rotateZ(-15. * M_PI / 180.);
0863       G4ThreeVector ladTranspos(ladder_x / 2., ladder_y / 2., 0.);
0864 
0865       G4RotationMatrix *lad_box_rotneg = new G4RotationMatrix();
0866       lad_box_rotneg->rotateZ(+15. * M_PI / 180.);
0867       G4ThreeVector ladTransneg(ladder_x / 2., -ladder_y / 2., 0.);
0868 
0869       G4VSolid *ladder_basebox = new G4Box(std::format("ladder_basebox_{}_{}", inttlayer, itype), ladder_x / 2., ladder_y / 2., hdi_z / 2.);
0870       G4VSolid *ladder_subtbox =
0871           new G4Box(std::format("ladder_subtbox_{}_{}", inttlayer, itype), stave_x / 1.5, ladder_y / 1.5, hdi_z / 1.);  // has to be longer in z to avoid coincident surface
0872 
0873       G4VSolid *ladder_box1 = new G4SubtractionSolid(std::format("ladder_box1_{}_{}", inttlayer, itype), ladder_basebox, ladder_subtbox, lad_box_rotpos, ladTranspos);
0874 
0875       G4VSolid *ladder_box = new G4SubtractionSolid(std::format("ladder_box_{}_{}", inttlayer, itype), ladder_box1, ladder_subtbox, lad_box_rotneg, ladTransneg);
0876 
0877       G4LogicalVolume *ladder_volume =
0878           new G4LogicalVolume(ladder_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")), std::format("ladder_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0879 
0880       G4VSolid *ladderext_basebox = new G4Box(std::format("ladderext_basebox_{}_{}", inttlayer, itype), ladder_x / 2., ladder_y / 2., hdiext_z / 2.);
0881       G4VSolid *ladderext_subtbox =
0882           new G4Box(std::format("ladderext_subtbox_{}_{}", inttlayer, itype), stave_x / 1.5, ladder_y / 1.5, hdiext_z / 1.);  // has to be longer in z to avoid coincident surface
0883 
0884       G4VSolid *ladderext_box1 = new G4SubtractionSolid(std::format("ladderext_box1_{}_{}", inttlayer, itype), ladderext_basebox, ladderext_subtbox, lad_box_rotpos, ladTranspos);
0885 
0886       G4VSolid *ladderext_box = new G4SubtractionSolid(std::format("ladderext_box_{}_{}", inttlayer, itype), ladderext_box1, ladderext_subtbox, lad_box_rotneg, ladTransneg);
0887 
0888       G4LogicalVolume *ladderext_volume =
0889           new G4LogicalVolume(ladderext_box, GetDetectorMaterial(rc->get_StringFlag("WorldMaterial")), std::format("ladderext_{}_{}", inttlayer, itype), nullptr, nullptr, nullptr);
0890       // the rotation matrices are just used by G4VSolid, ownership is not taken over
0891       delete lad_box_rotpos;
0892       delete lad_box_rotneg;
0893       m_DisplayAction->AddVolume(ladder_volume, "Ladder");
0894       m_DisplayAction->AddVolume(ladderext_volume, "Ladder");
0895 
0896       // Now add the components of the ladder to the ladder volume
0897       // The sensor is closest to the beam pipe, the stave cooler is furthest away
0898       // Note that the cooler has been assembled in the stave volume with the top at larger x, so the sensor will be at smaller x
0899       // That will be the configuration when the ladder is at phi = 0 degrees, the positive x direction
0900 
0901       // We start at the most positive x value and add the stave first
0902 
0903       // Carbon stave
0904       double TVstave_y = 0.0;
0905       const double TVstave_x = ladder_x / 2. - stave_x / 2.;
0906 
0907       new G4PVPlacement(nullptr, G4ThreeVector(TVstave_x, TVstave_y, 0.0), stave_volume, std::format("stave_{}_{}", inttlayer, itype), ladder_volume, false, 0, OverlapCheck());
0908       new G4PVPlacement(nullptr, G4ThreeVector(TVstave_x, TVstave_y, 0.0), staveext_volume, std::format("staveext_{}_{}", inttlayer, itype), ladderext_volume, false, 0, OverlapCheck());
0909 
0910       // HDI Kapton
0911       const double TVhdi_kapton_x = TVstave_x - stave_x / 2. - hdi_kapton_x / 2.;
0912       new G4PVPlacement(nullptr, G4ThreeVector(TVhdi_kapton_x, TVstave_y, 0.0), hdi_kapton_volume, std::format("hdikapton_{}_{}", inttlayer, itype), ladder_volume, false, 0,
0913                         OverlapCheck());
0914       new G4PVPlacement(nullptr, G4ThreeVector(TVhdi_kapton_x, TVstave_y, 0.0), hdiext_kapton_volume, std::format("hdiextkapton_{}_{}", inttlayer, itype), ladderext_volume, false, 0,
0915                         OverlapCheck());
0916 
0917       // HDI copper
0918       const double TVhdi_copper_x = TVhdi_kapton_x - hdi_kapton_x / 2. - hdi_copper_x / 2.;
0919       new G4PVPlacement(nullptr, G4ThreeVector(TVhdi_copper_x, TVstave_y, 0.0), hdi_copper_volume, std::format("hdicopper_{}_{}", inttlayer, itype), ladder_volume, false, 0,
0920                         OverlapCheck());
0921       new G4PVPlacement(nullptr, G4ThreeVector(TVhdi_copper_x, TVstave_y, 0.0), hdiext_copper_volume, std::format("hdiextcopper_{}_{}", inttlayer, itype), ladderext_volume, false, 0,
0922                         OverlapCheck());
0923 
0924       // Glue for Si-sensor
0925       const double TVsi_glue_x = TVhdi_copper_x - hdi_copper_x / 2. - si_glue_x / 2.;
0926       new G4PVPlacement(nullptr, G4ThreeVector(TVsi_glue_x, TVstave_y, 0.0), si_glue_volume, std::format("si_glue_{}_{}", inttlayer, itype), ladder_volume, false, 0, OverlapCheck());
0927 
0928       // Si-sensor
0929       double TVSi_y = 0.0;
0930       // sensor is not centered in y in the ladder volume for the Z sensitive ladders
0931       if (laddertype == PHG4InttDefs::SEGMENTATION_Z)
0932       {
0933         TVSi_y = +sensor_offset_y;
0934       }
0935       // const double TVSi_x = TVhdi_copper_x - hdi_copper_x / 2. - siactive_x / 2.;
0936       const double TVSi_x = TVsi_glue_x - si_glue_x / 2. - siactive_x / 2.;
0937       new G4PVPlacement(nullptr, G4ThreeVector(TVSi_x, TVSi_y, 0.0), siinactive_volume, std::format("siinactive_{}_{}", inttlayer, itype), ladder_volume, false, 0, OverlapCheck());
0938       new G4PVPlacement(nullptr, G4ThreeVector(TVSi_x, TVSi_y, 0.0), siactive_volume, std::format("siactive_{}_{}", inttlayer, itype), ladder_volume, false, 0, OverlapCheck());
0939 
0940       // FPHX glue
0941       const double TVfphx_glue_x = TVhdi_copper_x - hdi_copper_x / 2. - fphx_glue_x / 2.;
0942       double TVfphx_glue_y = sifull_y / 2. + gap_sensor_fphx + fphx_y / 2.;
0943       if (laddertype == PHG4InttDefs::SEGMENTATION_Z)
0944       {
0945         TVfphx_glue_y -= sensor_offset_y;
0946       }
0947 
0948       // laddertype PHG4InttDefs::SEGMENTATION_Z has only one FPHX, laddertype PHG4InttDefs::SEGMENTATION_PHI has two
0949       if (laddertype == PHG4InttDefs::SEGMENTATION_PHI)
0950       {
0951         new G4PVPlacement(nullptr, G4ThreeVector(TVfphx_glue_x, +TVfphx_glue_y, 0.0), fphx_gluecontainer_volume, std::format("fphx_gluecontainerp_{}_{}", inttlayer, itype),
0952                           ladder_volume, false, 0, OverlapCheck());
0953       }
0954       new G4PVPlacement(nullptr, G4ThreeVector(TVfphx_glue_x, -TVfphx_glue_y, 0.0), fphx_gluecontainer_volume, std::format("fphx_gluecontainerm_{}_{}", inttlayer, itype), ladder_volume,
0955                         false, 0, OverlapCheck());
0956 
0957       // FPHX container
0958       const double TVfphx_x = TVfphx_glue_x - fphx_glue_x / 2. - fphx_x / 2.;
0959       double TVfphx_y = sifull_y / 2. + gap_sensor_fphx + fphx_y / 2.;
0960       if (laddertype == PHG4InttDefs::SEGMENTATION_Z)
0961       {
0962         TVfphx_y -= sensor_offset_y;
0963       }
0964 
0965       // laddertype PHG4InttDefs::SEGMENTATION_Z has only one FPHX, laddertype PHG4InttDefs::SEGMENTATION_PHI has two
0966       if (laddertype == PHG4InttDefs::SEGMENTATION_PHI)
0967       {
0968         // new G4PVPlacement(0, G4ThreeVector(TVfphx_x, +TVfphx_y+100., 100.0), fphxcontainer_volume, std::format("fphxcontainerp_{}_{}", inttlayer, itype), ladder_volume, false,
0969         // 0, OverlapCheck());
0970         new G4PVPlacement(nullptr, G4ThreeVector(TVfphx_x, +TVfphx_y, 0.0), fphxcontainer_volume, std::format("fphxcontainerp_{}_{}", inttlayer, itype), ladder_volume, false, 0,
0971                           OverlapCheck());
0972       }
0973       new G4PVPlacement(nullptr, G4ThreeVector(TVfphx_x, -TVfphx_y, 0.0), fphxcontainer_volume, std::format("fphxcontainerm_{}_{}", inttlayer, itype), ladder_volume, false, 0,
0974                         OverlapCheck());
0975 
0976       // ----- Step 3 -----------------------------------------------------------------------------------------------
0977       // We install the section of ladder for this sensor at all requested phi values and at positive and negative Z
0978       //=============================================================================================================
0979 
0980       // Distribute Ladders in phi
0981       // We are still in the loops over layer and sensor type, we will place copies of the ladder section for this sensor
0982       //  at all ladder phi values, and at both positive and negative Z values.
0983 
0984       // given radius values are for the center of the sensor, we need the x offset from center of ladder to center of sensor so we can place the ladder
0985       double sensor_offset_x_ladder = 0.0 - TVSi_x;  // ladder center is at x = 0.0 by construction. Sensor is at lower x, so TVSi_x is negative
0986 
0987       const double dphi = 2 * M_PI / nladders_layer;
0988 
0989       m_PosZ[inttlayer][itype] = (itype == 0) ? hdi_z / 2. : hdi_z_arr[inttlayer][0] + hdi_z / 2.;  // location of center of ladder in Z
0990       m_StripOffsetX[inttlayer] = sensor_offset_x_ladder;
0991 
0992       // The sensors have no tilt in the new design
0993       //    The type 1 ladders have the sensor at the center of the ladder in phi, so that is easy
0994       //    The type 0 ladders are more complicated because the sensor center is perpendicular to the radial vector and the sensor is not at the ladder center
0995       //         We made the stave box symmetric in y around the sensor center to simplify things
0996 
0997       for (int icopy = 0; icopy < nladders_layer; icopy++)
0998       {
0999         // sensor center
1000         const double phi = offsetphi + dphi * icopy;  // if offsetphi is zero we start at zero
1001 
1002         double radius;
1003         // Make each layer at a single radius - i.e. what was formerly a sub-layer is now considered a layer
1004         radius = m_SensorRadius[inttlayer];
1005         radius += sensor_offset_x_ladder;
1006 
1007         InttNameSpace::Offline_s k_ladder;
1008         Eigen::Affine3d const *abs_transform_ladder;
1009         Eigen::Matrix3d rot_ladder;
1010         Eigen::Matrix3d sphnx2geant;
1011         sphnx2geant << 0, 1, 0,  //
1012             -1, 0, 0,            //
1013             0, 0, 1;             //
1014         Eigen::Matrix3d rot_ladder_geant;
1015         double rotation_z_survey = 0.0;
1016         Eigen::Vector3d translation_ladder;
1017         double posx_ladder_survey = -1E10;
1018         double posy_ladder_survey = -1E10;
1019         double posz_ladder_survey = -1E10;
1020 
1021         if (useSurvey)
1022         {
1023           // set up the key for the survey map
1024           k_ladder.layer = inttlayer + 3;
1025           k_ladder.ladder_phi = icopy;
1026           k_ladder.ladder_z = InttSurveyMap::Wildcard;  // Note: positive ladder_z = itype + 2, negative ladder_z = itype
1027           k_ladder.strip_x = InttSurveyMap::Wildcard;
1028           k_ladder.strip_y = InttSurveyMap::Wildcard;
1029 
1030           abs_transform_ladder = survey->GetAbsoluteTransform(k_ladder);
1031           if (!abs_transform_ladder)
1032           {
1033             std::cout << "Eigen::Affine3d absolute transform for (layer,ladder_phi,ladder_z,strip_phi,strip_z) = (" << k_ladder.layer << "," << k_ladder.ladder_phi << ","
1034                       << k_ladder.ladder_z << "," << k_ladder.strip_x << "," << k_ladder.strip_y << ") not valid" << std::endl;
1035             continue;
1036           }
1037 
1038           // Rotation
1039           rot_ladder = abs_transform_ladder->rotation();
1040           rot_ladder_geant = sphnx2geant * rot_ladder;
1041 
1042           // Rotation
1043           Eigen::Vector3d angles_ladder = rot_ladder_geant.eulerAngles(2, 0, 1);  // Tait-Bryan angles (z -> x' -> y'')
1044           // std::cout << "angles_ladder = " << angles_ladder << std::endl;
1045           if (inttlayer < 2)  // layer = 3 and 4
1046           {
1047             if (inttlayer == 0 && icopy == 0)
1048             {
1049               rotation_z_survey = angles_ladder[0] - M_PI;
1050             }
1051             else if ((inttlayer == 0 && icopy >= 7) || (inttlayer == 1 && icopy >= 6))
1052             {
1053               rotation_z_survey = angles_ladder[0] + M_PI;
1054             }
1055             else
1056             {
1057               rotation_z_survey = angles_ladder[0];
1058             }
1059           }
1060           else  // layer 5 and 6
1061           {
1062             if (inttlayer == 2 && icopy == 0)
1063             {
1064               rotation_z_survey = angles_ladder[0] - M_PI;
1065             }
1066             else if (icopy >= 9)
1067             {
1068               rotation_z_survey = angles_ladder[0] + M_PI;
1069             }
1070             else
1071             {
1072               rotation_z_survey = angles_ladder[0];
1073             }
1074           }
1075 
1076           // Translation
1077           translation_ladder = abs_transform_ladder->translation();
1078           posx_ladder_survey = translation_ladder.x();
1079           posy_ladder_survey = translation_ladder.y();
1080           posz_ladder_survey = translation_ladder.z();  //! Q: is this the whole ladder shift in the z direction?
1081         }
1082 
1083         double p = 0.0;
1084         if (laddertype == PHG4InttDefs::SEGMENTATION_Z)
1085         {
1086           // The Z sensitive ladders have the sensors offset in y relative to the ladder center
1087           // We have to slightly rotate the ladder in its own frame to make the radial vector to the sensor center normal to the sensor face
1088           p = atan(sensor_offset_y / radius);
1089           // then we adjust the distance to the center of the ladder to put the sensor at the requested distance from the center of the barrel
1090           radius /= cos(p);
1091         }
1092 
1093         // these describe the center of the ladder volume, placing it so that the center of the sensor is at phi = dphi * icopy, and at the correct radius
1094         const double posx_ladder_ideal = radius * cos(phi - p);
1095         const double posy_ladder_ideal = radius * sin(phi - p);
1096         // rotate in its own frame to make sensor perp to radial vector (p), then additionally rotate to account for ladder phi
1097         const double rotation_z_ideal = p + (phi - p) + offsetrot;
1098 
1099         const double posx = (useSurvey) ? posx_ladder_survey : posx_ladder_ideal;
1100         const double posy = (useSurvey) ? posy_ladder_survey : posy_ladder_ideal;
1101         const double laddercenter_shiftz = (useSurvey) ? params->get_double_param("ladder_center_avgshift_z") * cm : 0.0;
1102         const double fRotate = (useSurvey) ? rotation_z_survey : rotation_z_ideal;
1103         G4RotationMatrix ladderrotation;
1104         ladderrotation.rotateZ(fRotate);
1105 
1106         if (Verbosity() > 100 && useSurvey)
1107         {
1108           std::cout << "(dphi,layer,icopy,itype,p,phi,offsetrot,rotation_z_ideal,rotation_z_survey)=(" << dphi << "," << inttlayer << "," << icopy << "," << itype << "," << p << "," << phi
1109                     << "," << offsetrot << "," << rotation_z_ideal << "," << rotation_z_survey << ")" << std::endl;
1110         }
1111 
1112         // this placement version rotates the ladder in its own frame by fRotate, then translates the center to posx, posy, +/- m_PosZ
1113         auto *pointer_negz = new G4PVPlacement(G4Transform3D(ladderrotation, G4ThreeVector(posx, posy, -m_PosZ[inttlayer][itype] + laddercenter_shiftz)), ladder_volume,
1114                            std::format("ladder_{}_{}_{}_negz", inttlayer, itype, icopy), trackerenvelope, false, 0, OverlapCheck());
1115         auto *pointer_posz = new G4PVPlacement(G4Transform3D(ladderrotation, G4ThreeVector(posx, posy, +m_PosZ[inttlayer][itype] + laddercenter_shiftz)), ladder_volume,
1116                            std::format("ladder_{}_{}_{}_posz", inttlayer, itype, icopy), trackerenvelope, false, 0, OverlapCheck());
1117 
1118         std::string ladder_id = std::format("(inttlayer,itype,icopy)=({},{},{})", inttlayer, itype, icopy);
1119         if (Verbosity() > 100 && useSurvey && itype == 0)
1120         {
1121           std::cout << " G4PVPlacement for " << ladder_id << " ladder volume" << std::endl
1122                     << ladder_id << ", Rotation about the Z-axis (survey): " << rotation_z_survey << std::endl
1123                     << ladder_id << ", Rotation about the Z-axis (ideal): " << rotation_z_ideal << std::endl
1124                     << ladder_id << ", Translation (survey): " << G4ThreeVector(posx_ladder_survey, posy_ladder_survey, posz_ladder_survey) << std::endl
1125                     << ladder_id << ", Translation (ideal): " << G4ThreeVector(posx_ladder_ideal, posy_ladder_ideal, 0.0) << std::endl;
1126         }
1127 
1128         if (m_IsActiveMap.contains(inttlayer))
1129         {
1130           m_ActiveVolumeTuple.insert(std::make_pair(pointer_negz, std::make_tuple(inttlayer, itype, icopy, -1)));
1131           m_ActiveVolumeTuple.insert(std::make_pair(pointer_posz, std::make_tuple(inttlayer, itype, icopy, 1)));
1132         }
1133 
1134         // The net effect of the above manipulations for the Z sensitive ladders is that the center of the sensor is at dphi * icopy and at the requested radius
1135         // That us all that the geometry object needs to know, so no changes to that are necessary
1136 
1137         if (itype != 0)
1138         {
1139           // We have added the outer sensor above, now we add the HDI extension tab to the end of the outer sensor HDI
1140           const double posz_ext = (hdi_z_arr[inttlayer][0] + hdi_z) + hdiext_z / 2.;
1141 
1142           new G4PVPlacement(G4Transform3D(ladderrotation, G4ThreeVector(posx, posy, -posz_ext + laddercenter_shiftz)), ladderext_volume,
1143                             std::format("ladderext_{}_{}_{}_negz", inttlayer, itype, icopy), trackerenvelope, false, 0, OverlapCheck());
1144           new G4PVPlacement(G4Transform3D(ladderrotation, G4ThreeVector(posx, posy, +posz_ext + laddercenter_shiftz)), ladderext_volume,
1145                             std::format("ladderext_{}_{}_{}_posz", inttlayer, itype, icopy), trackerenvelope, false, 0, OverlapCheck());
1146         }
1147 
1148         // if (Verbosity() > 100)
1149         // {
1150         //     std::cout << "   Ladder copy " << icopy << " radius " << radius << " phi " << phi << " itype " << itype << " posz " << m_PosZ[inttlayer][itype] << " fRotate " << fRotate
1151         //               << " posx " << posx << " posy " << posy << std::endl
1152         //               << "Transformation for the positive ladder z (layer,ladder_phi,ladder_z,strip_phi,strip_z)=(" << k_ladder.layer << "," << k_ladder.ladder_phi << "," << k_ladder.ladder_z
1153         //               << "," << k_ladder.strip_phi << "," << k_ladder.strip_z << ")" << std::endl
1154         //               << "Ideal case, fRotate (radian) around the Z axis = " << fRotate << std::endl
1155         //               << "G4RotationMatrix ladderrotation (default,ideal)" << std::endl
1156         //               << ladderrotation << std::endl
1157         //               << "Eigen::Affine3D transform matrix (from InttSurveyMap)" << std::endl
1158         //               << abs_transform_ladder->matrix() << std::endl
1159         //               << "Translation vector (absolute)" << std::endl
1160         //               << translation_ladder << std::endl
1161         //               << "(inttlayer,icopy,itype)=(" << k_ladder.layer << "," << icopy << "," << itype << "), difference between default and survey (dX, dY)=(" << (posx_ladder - posx) << ","
1162         //               << (posy_ladder - posy) << ")" << std::endl
1163         //               << "Rotation matrix (absolute, before the temporary transform)" << std::endl
1164         //               << abs_transform_ladder->rotation()
1165         //               << std::endl
1166         //               << "Rotation matrix (absolute, after the temporary transform)" << std::endl
1167         //               << ladderrotation_survey << std::endl
1168         //               << "------------------------------------------------------------------------------------------" << std::endl;
1169         // }
1170       }  // end loop over ladder copy placement in phi and positive and negative Z
1171     }    // end loop over inner or outer sensor
1172   }      // end loop over layers
1173 
1174   // Finally, we add some support material for the silicon detectors
1175 
1176   //
1177   /*
1178     4 rails, which are 12mm OD and 9mm ID tubes at a radius of 168.5 mm.  They are spaced equidistantly in phi.
1179     The rails run along the entire length of the TPC and even stick out of the TPC, but I think for the moment you don't have to put the parts that stick out in the simulation.
1180     An inner skin with a ID at 62.416 mm and a thickness of 0.250 mm.
1181     An outer skin with a ID at 120.444 mm and a sandwich of 0.25 mm cfc, 1.5 mm foam and 0.25 mm cfc.
1182 
1183     All of the above are carbon fiber.
1184   */
1185   const PHParameters *supportparams = m_ParamsContainer->GetParameters(PHG4InttDefs::SUPPORTPARAMS);
1186 
1187   // rails
1188   G4Tubs *rail_tube = new G4Tubs(std::format("si_support_rail"), supportparams->get_double_param("rail_inner_radius") * cm, supportparams->get_double_param("rail_outer_radius") * cm,
1189                                  supportparams->get_double_param("rail_length") * cm / 2.0, 0, 2.0 * M_PI);
1190   G4LogicalVolume *rail_volume = new G4LogicalVolume(rail_tube, GetDetectorMaterial("CFRP_INTT"), "rail_volume", nullptr, nullptr, nullptr);
1191   if (m_IsSupportActive > 0)
1192   {
1193     m_PassiveVolumeTuple.insert(std::make_pair(rail_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::SUPPORT_RAIL)));
1194   }
1195   m_DisplayAction->AddVolume(rail_volume, "Rail");
1196 
1197   double rail_dphi = supportparams->get_double_param("rail_dphi") * deg / rad;
1198   double rail_phi_start = supportparams->get_double_param("rail_phi_start") * deg / rad;
1199   double rail_radius = supportparams->get_double_param("rail_radius") * cm;
1200   for (int i = 0; i < 4; i++)
1201   {
1202     double phi = rail_phi_start + i * rail_dphi;
1203 
1204     // place a copy at each rail phi value
1205     const double posx = rail_radius * cos(phi);
1206     const double posy = rail_radius * sin(phi);
1207 
1208     new G4PVPlacement(nullptr, G4ThreeVector(posx, posy, 0.0), rail_volume, std::format("si_support_rail_{}", i), trackerenvelope, false, 0, OverlapCheck());
1209   }
1210 
1211   // Outer skin
1212   // G4Tubs *outer_skin_cfcin_tube = new G4Tubs("si_outer_skin_cfcin",
1213   //                                            supportparams->get_double_param("outer_skin_cfcin_inner_radius") * cm,
1214   //                                            supportparams->get_double_param("outer_skin_cfcin_outer_radius") * cm,
1215   //                                            supportparams->get_double_param("outer_skin_cfcin_length") * cm / 2.,
1216   //                                            0, 2.0 * M_PI);
1217   // G4LogicalVolume *outer_skin_cfcin_volume = new G4LogicalVolume(outer_skin_cfcin_tube, GetDetectorMaterial("CFRP_INTT"),
1218   //                                                                "outer_skin_cfcin_volume", 0, 0, 0);
1219 
1220   // G4Tubs *outer_skin_foam_tube = new G4Tubs("si_outer_skin_foam",
1221   //                                           supportparams->get_double_param("outer_skin_foam_inner_radius") * cm,
1222   //                                           supportparams->get_double_param("outer_skin_foam_outer_radius") * cm,
1223   //                                           supportparams->get_double_param("outer_skin_foam_length") * cm / 2.,
1224   //                                           0, 2.0 * M_PI);
1225   // G4LogicalVolume *outer_skin_foam_volume = new G4LogicalVolume(outer_skin_foam_tube, GetDetectorMaterial("ROHACELL_FOAM_110"),
1226   //                                                               "outer_skin_foam_volume", 0, 0, 0);
1227 
1228   // G4Tubs *outer_skin_cfstd::cout_tube = new G4Tubs("si_outer_skin_cfstd::cout",
1229   //                                             supportparams->get_double_param("outer_skin_cfstd::cout_inner_radius") * cm,
1230   //                                             supportparams->get_double_param("outer_skin_cfstd::cout_outer_radius") * cm,
1231   //                                             supportparams->get_double_param("outer_skin_cfstd::cout_length") * cm / 2.,
1232   //                                             0, 2.0 * M_PI);
1233   // G4LogicalVolume *outer_skin_cfstd::cout_volume = new G4LogicalVolume(outer_skin_cfstd::cout_tube, GetDetectorMaterial("CFRP_INTT"),
1234   //                                                                 "outer_skin_cfstd::cout_volume", 0, 0, 0);
1235 
1236   G4Tubs *outer_skin_tube = new G4Tubs("si_outer_skin", supportparams->get_double_param("outer_skin_inner_radius") * cm, supportparams->get_double_param("outer_skin_outer_radius") * cm,
1237                                        supportparams->get_double_param("outer_skin_length") * cm / 2., 0, 2.0 * M_PI);
1238   G4LogicalVolume *outer_skin_volume = new G4LogicalVolume(outer_skin_tube, GetDetectorMaterial("CFRP_INTT"), "outer_skin_volume", nullptr, nullptr, nullptr);
1239 
1240   if (m_IsSupportActive > 0)
1241   {
1242     // m_PassiveVolumeTuple.insert(std::make_pair(outer_skin_cfcin_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::INTT_OUTER_SKIN)));
1243     // m_PassiveVolumeTuple.insert(std::make_pair(outer_skin_foam_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::INTT_OUTER_SKIN)));
1244     // m_PassiveVolumeTuple.insert(std::make_pair(outer_skin_cfstd::cout_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::INTT_OUTER_SKIN)));
1245     m_PassiveVolumeTuple.insert(std::make_pair(outer_skin_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::INTT_OUTER_SKIN)));
1246   }
1247   // m_DisplayAction->AddVolume(outer_skin_cfcin_volume, "Skin");
1248   // m_DisplayAction->AddVolume(outer_skin_foam_volume, "Skin");
1249   // m_DisplayAction->AddVolume(outer_skin_cfstd::cout_volume, "Skin");
1250   // new G4PVPlacement(0, G4ThreeVector(0, 0.0, 0), outer_skin_cfcin_volume,
1251   //                   "si_support_outer_skin_cfcin", trackerenvelope, false, 0, OverlapCheck());
1252   // new G4PVPlacement(0, G4ThreeVector(0, 0.0), outer_skin_foam_volume,
1253   //                   "si_support_outer_skin_foam", trackerenvelope, false, 0, OverlapCheck());
1254   // new G4PVPlacement(0, G4ThreeVector(0, 0.0), outer_skin_cfstd::cout_volume,
1255   //                   "si_support_outer_skin_cfstd::cout", trackerenvelope, false, 0, OverlapCheck());
1256   m_DisplayAction->AddVolume(outer_skin_volume, "Skin");
1257 
1258   double si_support_skin_center_x = (useSurvey) ? supportparams->get_double_param("endcap_ring_x") * cm : 0.0;
1259   double si_support_skin_center_y = (useSurvey) ? supportparams->get_double_param("endcap_ring_y") * cm : 0.0;
1260   new G4PVPlacement(nullptr, G4ThreeVector(si_support_skin_center_x, si_support_skin_center_y, 0), outer_skin_volume,
1261                     "si_support_outer_skin_cfcin", trackerenvelope, false, 0, OverlapCheck());
1262 
1263   // Inner skin
1264   G4Tubs *inner_skin_tube = new G4Tubs("si_inner_skin", supportparams->get_double_param("inner_skin_inner_radius") * cm, supportparams->get_double_param("inner_skin_outer_radius") * cm,
1265                                        supportparams->get_double_param("inner_skin_length") * cm / 2., 0, 2.0 * M_PI);
1266   G4LogicalVolume *inner_skin_volume = new G4LogicalVolume(inner_skin_tube, GetDetectorMaterial("CFRP_INTT"), "inner_skin_volume", nullptr, nullptr, nullptr);
1267   if (m_IsSupportActive > 0)
1268   {
1269     m_PassiveVolumeTuple.insert(std::make_pair(inner_skin_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::INTT_INNER_SKIN)));
1270   }
1271   m_DisplayAction->AddVolume(inner_skin_volume, "Skin");
1272 
1273   new G4PVPlacement(nullptr, G4ThreeVector(si_support_skin_center_x, si_support_skin_center_y), inner_skin_volume, "si_support_inner_skin",
1274                     trackerenvelope, false, 0, OverlapCheck());
1275 
1276   // Service barrel outer  ////////////////////////////////////////////////////////////////////////////////////
1277   G4Tubs *service_barrel_outer_tube =
1278       new G4Tubs("si_service_barrel_outer", supportparams->get_double_param("service_barrel_outer_inner_radius") * cm, supportparams->get_double_param("service_barrel_outer_outer_radius") * cm,
1279                  supportparams->get_double_param("service_barrel_outer_length") * cm / 2., 0, 2.0 * M_PI);
1280   G4LogicalVolume *service_barrel_outer_volume = new G4LogicalVolume(service_barrel_outer_tube, GetDetectorMaterial("CFRP_INTT"), "service_barrel_outer_volume", nullptr, nullptr, nullptr);
1281   if (m_IsSupportActive > 0)
1282   {
1283     m_PassiveVolumeTuple.insert(std::make_pair(service_barrel_outer_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::SERVICE_BARREL_OUTER)));
1284   }
1285   m_DisplayAction->AddVolume(service_barrel_outer_volume, "Skin");
1286 
1287   new G4PVPlacement(nullptr, G4ThreeVector(0, 0.0), service_barrel_outer_volume, "si_support_service_barrel_outer", trackerenvelope, false, 0, OverlapCheck());
1288 
1289   // Support Tube  ////////////////////////////////////////////////////////////////////////////////////
1290   G4Tubs *support_tube_tube = new G4Tubs("si_support_tube", supportparams->get_double_param("support_tube_inner_radius") * cm, supportparams->get_double_param("support_tube_outer_radius") * cm,
1291                                          supportparams->get_double_param("support_tube_length") * cm / 2., 0, 2.0 * M_PI);
1292   G4LogicalVolume *support_tube_volume = new G4LogicalVolume(support_tube_tube, GetDetectorMaterial("CFRP_INTT"), "support_tube_volume", nullptr, nullptr, nullptr);
1293   if (m_IsSupportActive > 0)
1294   {
1295     m_PassiveVolumeTuple.insert(std::make_pair(support_tube_volume, std::make_tuple(PHG4InttDefs::SUPPORT_DETID, PHG4InttDefs::SUPPORT_TUBE)));
1296   }
1297   m_DisplayAction->AddVolume(support_tube_volume, "Skin");
1298 
1299   new G4PVPlacement(nullptr, G4ThreeVector(0, 0.0), support_tube_volume, "si_support_support_tube", trackerenvelope, false, 0, OverlapCheck());
1300 
1301   // Endcap ring in simulations = Endcap rings + endcap staves
1302   int inttlayer = (m_LayerBeginEndIteratorPair.first)->second;
1303   const PHParameters *params1 = m_ParamsContainer->GetParameters(inttlayer);
1304   const int laddertype = params1->get_int_param("laddertype");
1305   const PHParameters *params = m_ParamsContainer->GetParameters(laddertype);
1306 
1307   // Carbon PEEK ring ////////////////////////////////////////////////////////////////////////////////////
1308   G4Tubs *endcap_CP_ring = new G4Tubs("endcap_CP_ring", supportparams->get_double_param("endcap_CPring_inner_radius") * cm, supportparams->get_double_param("endcap_CPring_outer_radius") * cm,
1309                                       supportparams->get_double_param("endcap_CPring_length") * cm / 2., 0, 2.0 * M_PI);
1310 
1311   G4LogicalVolume *endcap_CP_ring_volume = new G4LogicalVolume(endcap_CP_ring, GetDetectorMaterial("CF30_PEEK70"), "endcap_CP_ring_volume", nullptr, nullptr, nullptr);
1312   m_DisplayAction->AddVolume(endcap_CP_ring_volume, "EndcapCPRing");
1313 
1314   // new Al-PEEK ring from Jan/2021    //////////////////////////////////////////////////////////////////////////
1315   // outermost part (Al)
1316   G4Tubs *endcap_AlPEEK_Alring_1 =
1317       new G4Tubs("endcap_AlPEEK_Alring_1", supportparams->get_double_param("endcap_AlPEEK_Alring_1_inner_radius") * cm, supportparams->get_double_param("endcap_AlPEEK_Alring_1_outer_radius") * cm,
1318                  supportparams->get_double_param("endcap_AlPEEK_Alring_length") * cm / 2., 0, 2.0 * M_PI);
1319 
1320   G4LogicalVolume *endcap_AlPEEK_Alring_1_volume = new G4LogicalVolume(endcap_AlPEEK_Alring_1, GetDetectorMaterial("G4_Al"), "endcap_AlPEEK_Alring_1_volume", nullptr, nullptr, nullptr);
1321   m_DisplayAction->AddVolume(endcap_AlPEEK_Alring_1_volume, "EndcapAlPEEK_Al1");
1322 
1323   // 2nd outermost part (Carbon PEEK)
1324   G4Tubs *endcap_AlPEEK_Cring_1 =
1325       new G4Tubs("endcap_AlPEEK_Cring_1", supportparams->get_double_param("endcap_AlPEEK_Cring_1_inner_radius") * cm, supportparams->get_double_param("endcap_AlPEEK_Cring_1_outer_radius") * cm,
1326                  supportparams->get_double_param("endcap_AlPEEK_Cring_length") * cm / 2., 0, 2.0 * M_PI);
1327 
1328   G4LogicalVolume *endcap_AlPEEK_Cring_1_volume = new G4LogicalVolume(endcap_AlPEEK_Cring_1, GetDetectorMaterial("CF30_PEEK70"), "endcap_AlPEEK_Cring_1_volume", nullptr, nullptr, nullptr);
1329   m_DisplayAction->AddVolume(endcap_AlPEEK_Cring_1_volume, "EndcapAlPEEK_C1");
1330 
1331   // 3rd outermost part (Al)
1332   G4Tubs *endcap_AlPEEK_Alring_2 =
1333       new G4Tubs("endcap_AlPEEK_Alring_2", supportparams->get_double_param("endcap_AlPEEK_Alring_2_inner_radius") * cm, supportparams->get_double_param("endcap_AlPEEK_Alring_2_outer_radius") * cm,
1334                  supportparams->get_double_param("endcap_AlPEEK_Alring_length") * cm / 2., 0, 2.0 * M_PI);
1335 
1336   G4LogicalVolume *endcap_AlPEEK_Alring_2_volume = new G4LogicalVolume(endcap_AlPEEK_Alring_2, GetDetectorMaterial("G4_Al"), "endcap_AlPEEK_Alring_2_volume", nullptr, nullptr, nullptr);
1337   m_DisplayAction->AddVolume(endcap_AlPEEK_Alring_2_volume, "EndcapAlPEEK_Al2");
1338 
1339   // 4th outermost part (Carbon PEEK)
1340   G4Tubs *endcap_AlPEEK_Cring_2 =
1341       new G4Tubs("endcap_AlPEEK_Cring_2", supportparams->get_double_param("endcap_AlPEEK_Cring_2_inner_radius") * cm, supportparams->get_double_param("endcap_AlPEEK_Cring_2_outer_radius") * cm,
1342                  supportparams->get_double_param("endcap_AlPEEK_Cring_length") * cm / 2., 0, 2.0 * M_PI);
1343 
1344   G4LogicalVolume *endcap_AlPEEK_Cring_2_volume = new G4LogicalVolume(endcap_AlPEEK_Cring_2, GetDetectorMaterial("CF30_PEEK70"), "endcap_AlPEEK_Cring_2_volume", nullptr, nullptr, nullptr);
1345   m_DisplayAction->AddVolume(endcap_AlPEEK_Cring_2_volume, "EndcapAlPEEK_C2");
1346 
1347   // 5th outermost part (innermost) (Al)
1348   G4Tubs *endcap_AlPEEK_Alring_3 =
1349       new G4Tubs("endcap_AlPEEK_Alring_3", supportparams->get_double_param("endcap_AlPEEK_Alring_3_inner_radius") * cm, supportparams->get_double_param("endcap_AlPEEK_Alring_3_outer_radius") * cm,
1350                  supportparams->get_double_param("endcap_AlPEEK_Alring_length") * cm / 2., 0, 2.0 * M_PI);
1351 
1352   G4LogicalVolume *endcap_AlPEEK_Alring_3_volume = new G4LogicalVolume(endcap_AlPEEK_Alring_3, GetDetectorMaterial("G4_Al"), "endcap_AlPEEK_Alring_3_volume", nullptr, nullptr, nullptr);
1353   m_DisplayAction->AddVolume(endcap_AlPEEK_Alring_3_volume, "EndcapAlPEEK_Al3");
1354 
1355   if (m_IsEndcapActive)
1356   {
1357     double endcap_outer_edge_z = 0.0;                           // absolute z-coordinate of outer edge (furthest side from the origin) of the endcap, used for bus extender
1358     if (supportparams->get_int_param("endcap_ring_type") == 0)  // Place Al endcap rings
1359     {
1360       // Aluminum ring
1361       G4Tubs *endcap_Al_ring = new G4Tubs("endcap_Al_ring", supportparams->get_double_param("endcap_Alring_inner_radius") * cm,
1362                                           supportparams->get_double_param("endcap_Alring_outer_radius") * cm, supportparams->get_double_param("endcap_Alring_length") * cm / 2., 0, 2.0 * M_PI);
1363 
1364       G4LogicalVolume *endcap_Al_ring_volume = new G4LogicalVolume(endcap_Al_ring, GetDetectorMaterial("Al6061T6"), "endcap_Al_ring_volume", nullptr, nullptr, nullptr);
1365 
1366       // Stainlees steal ring
1367       G4Tubs *endcap_SS_ring = new G4Tubs("endcap_SS_ring", supportparams->get_double_param("endcap_SSring_inner_radius") * cm,
1368                                           supportparams->get_double_param("endcap_SSring_outer_radius") * cm, supportparams->get_double_param("endcap_SSring_length") * cm / 2., 0, 2.0 * M_PI);
1369 
1370       G4LogicalVolume *endcap_SS_ring_volume = new G4LogicalVolume(endcap_SS_ring, GetDetectorMaterial("SS316"), "endcap_SS_ring_volume", nullptr, nullptr, nullptr);
1371 
1372       // Water Glycol ring
1373       G4Tubs *endcap_WG_ring = new G4Tubs("endcap_WG_ring", supportparams->get_double_param("endcap_WGring_inner_radius") * cm,
1374                                           supportparams->get_double_param("endcap_WGring_outer_radius") * cm, supportparams->get_double_param("endcap_WGring_length") * cm / 2., 0, 2.0 * M_PI);
1375 
1376       G4LogicalVolume *endcap_WG_ring_volume = new G4LogicalVolume(endcap_WG_ring, GetDetectorMaterial("WaterGlycol_INTT"), "endcap_WG_ring_volume", nullptr, nullptr, nullptr);
1377 
1378       double endcap_ring_z = supportparams->get_double_param("endcap_ring_z") * cm;
1379       for (int i = 0; i < 2; i++)  // i=0 : positive z, i=1 negative z
1380       {
1381         endcap_ring_z = (i == 0) ? endcap_ring_z : -1.0 * endcap_ring_z;
1382 
1383         double width_WGring_z = supportparams->get_double_param("endcap_WGring_length") * cm;
1384         double width_SSring_z = supportparams->get_double_param("endcap_SSring_length") * cm;
1385         double width_Alring_z = supportparams->get_double_param("endcap_Alring_length") * cm;
1386 
1387         for (int j = 0; j < 2; j++)  // j=0 : positive side z, j=1 negative side z
1388         {
1389           width_WGring_z = (j == 0) ? width_WGring_z : -1.0 * width_WGring_z;
1390           width_SSring_z = (j == 0) ? width_SSring_z : -1.0 * width_SSring_z;
1391           width_Alring_z = (j == 0) ? width_Alring_z : -1.0 * width_Alring_z;
1392 
1393           double cent_WGring_z = endcap_ring_z + width_WGring_z / 2.;
1394           double cent_SSring_z = endcap_ring_z + width_WGring_z + width_SSring_z / 2.;
1395           double cent_Alring_z = endcap_ring_z + width_WGring_z + width_SSring_z + width_Alring_z / 2.;
1396           endcap_outer_edge_z = fabs(endcap_ring_z) + fabs(width_WGring_z + width_SSring_z + width_Alring_z / 2.);  // absolute distance from origin
1397 
1398           new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_WGring_z), endcap_WG_ring_volume, std::format("endcap_WG_ring_pv_{}_{}", i, j), trackerenvelope, false, 0,
1399                             OverlapCheck());
1400 
1401           new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_SSring_z), endcap_SS_ring_volume, std::format("endcap_SS_ring_pv_{}_{}", i, j), trackerenvelope, false, 0,
1402                             OverlapCheck());
1403 
1404           new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_Alring_z), endcap_Al_ring_volume, std::format("endcap_Al_ring_pv_{}_{}", i, j), trackerenvelope, false, 0,
1405                             OverlapCheck());
1406         }                                                            // end of the loop over positive/negative sides with j
1407       }                                                              // end of the loop over positive/negative sides with i
1408     }                                                                // end of endcap_ring_type == 0
1409     else if (supportparams->get_int_param("endcap_ring_type") == 1)  // Place CP endcap rings
1410     {
1411       double endcap_ring_z = supportparams->get_double_param("endcap_CPring_z") * cm;
1412 
1413       for (int i = 0; i < 2; i++)  // i=0 : positive z, i=1 negative z
1414       {
1415         endcap_ring_z = (i == 0) ? endcap_ring_z : -1.0 * endcap_ring_z;
1416         endcap_outer_edge_z = fabs(endcap_ring_z);
1417 
1418         new G4PVPlacement(nullptr, G4ThreeVector(0, 0, endcap_ring_z), endcap_CP_ring_volume, std::format("endcap_CP_ring_pv_{}",  i), trackerenvelope, false, 0, OverlapCheck());
1419 
1420       }                                                              // end of the loop over positive/negative sides
1421     }                                                                // end of endcap_ring_type == 1
1422     else if (supportparams->get_int_param("endcap_ring_type") == 2)  // the new endcap introduced in Jan/2021
1423     {
1424       double si_0_width =
1425           params->get_double_param("strip_z_0") * params->get_int_param("nstrips_z_sensor_0") * cm + 2 * params->get_double_param("sensor_edge_z") * cm;  // length of the smaller cells
1426       double si_1_width =
1427           params->get_double_param("strip_z_1") * params->get_int_param("nstrips_z_sensor_1") * cm + 2 * params->get_double_param("sensor_edge_z") * cm;  // length of the larger cells
1428       double sifull_width = si_0_width + si_1_width;                                                                                                      // length of the Si module
1429       double hdi_width = sifull_width + params->get_double_param("hdi_edge_z") * cm;
1430       double hdiext_width = params->get_double_param("halfladder_inside_z") * cm - sifull_width;
1431       double hdifull_width = hdi_width + hdiext_width;  // length of a half lader
1432       double endcap_ring_z = hdifull_width + supportparams->get_double_param("endcap_AlPEEK_Cring_length") / 2.0 * cm;
1433 
1434       for (int i = 0; i < 2; i++)  // i=0 : positive z, i=1 negative z
1435       {
1436         endcap_ring_z = (i == 0) ? endcap_ring_z : -1.0 * endcap_ring_z;
1437         endcap_outer_edge_z = fabs(endcap_ring_z) + supportparams->get_double_param("endcap_AlPEEK_Alring_length") * cm / 2.0;
1438         double endcap_ring_x = (useSurvey) ? supportparams->get_double_param("endcap_ring_x") * cm : 0.0;
1439         double endcap_ring_y = (useSurvey) ? supportparams->get_double_param("endcap_ring_y") * cm : 0.0;
1440         double endcap_ring_z_shift = (useSurvey) ? params->get_double_param("ladder_center_avgshift_z") * cm : 0.0;
1441 
1442         new G4PVPlacement(nullptr, G4ThreeVector(endcap_ring_x, endcap_ring_y, endcap_ring_z + endcap_ring_z_shift), endcap_AlPEEK_Alring_1_volume, std::format("endcap_AlPEEK_Alring_1_pv_{}",  i), trackerenvelope, false, 0,
1443                           OverlapCheck());
1444 
1445         new G4PVPlacement(nullptr, G4ThreeVector(endcap_ring_x, endcap_ring_y, endcap_ring_z + endcap_ring_z_shift), endcap_AlPEEK_Cring_1_volume, std::format("endcap_AlPEEK_Cring_1_pv_{}",  i), trackerenvelope, false, 0,
1446                           OverlapCheck());
1447 
1448         new G4PVPlacement(nullptr, G4ThreeVector(endcap_ring_x, endcap_ring_y, endcap_ring_z + endcap_ring_z_shift), endcap_AlPEEK_Alring_2_volume, std::format("endcap_AlPEEK_Alring_2_pv_{}",  i), trackerenvelope, false, 0,
1449                           OverlapCheck());
1450 
1451         new G4PVPlacement(nullptr, G4ThreeVector(endcap_ring_x, endcap_ring_y, endcap_ring_z + endcap_ring_z_shift), endcap_AlPEEK_Cring_2_volume, std::format("endcap_AlPEEK_Cring_2_pv_{}",  i), trackerenvelope, false, 0,
1452                           OverlapCheck());
1453 
1454         new G4PVPlacement(nullptr, G4ThreeVector(endcap_ring_x, endcap_ring_y, endcap_ring_z + endcap_ring_z_shift), endcap_AlPEEK_Alring_3_volume, std::format("endcap_AlPEEK_Alring_3_pv_{}",  i), trackerenvelope, false, 0,
1455                           OverlapCheck());
1456       }
1457     }
1458 
1459     /////////////////////////////////////////////////////////////////////
1460     // Mimic cylinders for the bus extender (very simplified for the moment)
1461     double bus_extender_radius_innermost = supportparams->get_double_param("bus_extender_radius") * cm;
1462     double bus_extender_length = supportparams->get_double_param("bus_extender_length") * cm;
1463     double bus_extender_copper_x = supportparams->get_double_param("bus_extender_copper_x") * cm;
1464     double bus_extender_kapton_x = supportparams->get_double_param("bus_extender_kapton_x") * cm;
1465 
1466     // copper layer of the bus extender for the inner barrel
1467     double inner_radius = bus_extender_radius_innermost;
1468     G4Tubs *bus_extender_copper_inner = new G4Tubs("bus_extender_coppe_innerr", inner_radius, inner_radius + bus_extender_copper_x, bus_extender_length / 2.0, 0, 2.0 * M_PI);
1469 
1470     G4LogicalVolume *bus_extender_copper_inner_volume = new G4LogicalVolume(bus_extender_copper_inner, GetDetectorMaterial("G4_Cu"), "bus_extender_copper_inner_volume", nullptr, nullptr, nullptr);
1471     m_DisplayAction->AddVolume(bus_extender_copper_inner_volume, "BusExtenderCopperInner");
1472 
1473     // kapton layer of the bus extender for the inner barrel
1474     inner_radius += bus_extender_copper_x;
1475     G4Tubs *bus_extender_kapton_inner = new G4Tubs("bus_extender_kapton_inner", inner_radius, inner_radius + bus_extender_kapton_x, bus_extender_length / 2.0, 0, 2.0 * M_PI);
1476 
1477     G4LogicalVolume *bus_extender_kapton_inner_volume = new G4LogicalVolume(bus_extender_kapton_inner, GetDetectorMaterial("G4_KAPTON"), "bus_extender_kapton_inner_volume", nullptr, nullptr, nullptr);
1478     m_DisplayAction->AddVolume(bus_extender_kapton_inner_volume, "BusExtenderKaptonInner");
1479 
1480     // copper layer of the bus extender for the outer outerbarrel
1481     inner_radius += bus_extender_kapton_x;
1482     G4Tubs *bus_extender_copper_outer = new G4Tubs("bus_extender_copper_outer", inner_radius, inner_radius + bus_extender_copper_x, bus_extender_length / 2.0, 0, 2.0 * M_PI);
1483 
1484     G4LogicalVolume *bus_extender_copper_outer_volume = new G4LogicalVolume(bus_extender_copper_outer, GetDetectorMaterial("G4_Cu"), "bus_extender_copper_outer_volume", nullptr, nullptr, nullptr);
1485     m_DisplayAction->AddVolume(bus_extender_copper_outer_volume, "BusExtenderCopperOuter");
1486 
1487     // kapton layer of the bus extender for the outer barrel
1488     inner_radius += bus_extender_copper_x;
1489     G4Tubs *bus_extender_kapton_outer = new G4Tubs("bus_extender_kapton_outer", inner_radius, inner_radius + bus_extender_kapton_x, bus_extender_length / 2.0, 0, 2.0 * M_PI);
1490 
1491     G4LogicalVolume *bus_extender_kapton_outer_volume = new G4LogicalVolume(bus_extender_kapton_outer, GetDetectorMaterial("G4_KAPTON"), "bus_extender_kapton_outer_volume", nullptr, nullptr, nullptr);
1492     m_DisplayAction->AddVolume(bus_extender_kapton_outer_volume, "BusExtenderKaptonOuter");
1493 
1494     double bus_extender_z = endcap_outer_edge_z;
1495     for (int i = 0; i < 2; i++)  // loop over both positive and negative sides to put the cylinders,  i=0 : positive z, i=1 negative z
1496     {
1497       // copper layer of bus extender for the inner barrel
1498       double cent_bus_extender_z = bus_extender_z + bus_extender_length / 2.0;
1499       cent_bus_extender_z *= (i == 0) ? 1.0 : -1.0;
1500       new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_bus_extender_z), bus_extender_copper_inner_volume, std::format("bus_extender_copper_inner_layer_pv_{}",  i), trackerenvelope, false,
1501                         0, OverlapCheck());
1502 
1503       // kapton layer of bus extender for the inner barrel
1504       new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_bus_extender_z), bus_extender_kapton_inner_volume, std::format("bus_extender_kapton_inner_layer_pv_{}",  i), trackerenvelope, false,
1505                         0, OverlapCheck());
1506 
1507       // copper layer of bus extender for the outer barrel
1508       new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_bus_extender_z), bus_extender_copper_outer_volume, std::format("bus_extender_copper_outer_layer_pv_{}",  i), trackerenvelope, false,
1509                         0, OverlapCheck());
1510 
1511       // kapton layer of bus extender for the outer barrel
1512       new G4PVPlacement(nullptr, G4ThreeVector(0, 0, cent_bus_extender_z), bus_extender_kapton_outer_volume, std::format("bus_extender_kapton_outer_layer_pv_{}",  i), trackerenvelope, false,
1513                         0, OverlapCheck());
1514     }
1515   }
1516 
1517   delete survey;
1518 
1519   return 0;
1520 }
1521 
1522 void PHG4InttDetector::AddGeometryNode()
1523 {
1524   int active = 0;
1525   //  std::map<int, int>::const_iterator iter;
1526   for (auto &iter : m_IsActiveMap)
1527   {
1528     if (iter.second > 0)
1529     {
1530       active = 1;
1531       break;
1532     }
1533   }
1534   if (active)
1535   {
1536     std::string geonode = (m_SuperDetector != "NONE") ? std::format("CYLINDERGEOM_{}", m_SuperDetector) : std::format("CYLINDERGEOM_{}", m_DetectorType);
1537 
1538     PHG4CylinderGeomContainer *geo = findNode::getClass<PHG4CylinderGeomContainer>(topNode(), geonode);
1539     if (!geo)
1540     {
1541       geo = new PHG4CylinderGeomContainer();
1542       PHNodeIterator iter(topNode());
1543       PHCompositeNode *runNode = dynamic_cast<PHCompositeNode *>(iter.findFirst("PHCompositeNode", "RUN"));
1544       PHNodeIterator runiter(runNode);
1545       PHCompositeNode *geomNode = dynamic_cast<PHCompositeNode *>(runiter.findFirst("PHCompositeNode", "RECO_TRACKING_GEOMETRY"));
1546       if(!geomNode)
1547       {
1548         geomNode = new PHCompositeNode("RECO_TRACKING_GEOMETRY");
1549         runNode->addNode(geomNode);
1550       }
1551       PHIODataNode<PHObject> *newNode = new PHIODataNode<PHObject>(geo, geonode, "PHObject");
1552       geomNode->addNode(newNode);
1553     }
1554 
1555     for (auto layeriter = m_LayerBeginEndIteratorPair.first; layeriter != m_LayerBeginEndIteratorPair.second; ++layeriter)
1556     {
1557       const int sphxlayer = layeriter->first;
1558       const int inttlayer = layeriter->second;
1559       int ilayer = inttlayer;
1560       const PHParameters *params_layer = m_ParamsContainer->GetParameters(inttlayer);
1561       const int laddertype = params_layer->get_int_param("laddertype");
1562       // parameters are stored in cm per our convention
1563       const PHParameters *params = m_ParamsContainer->GetParameters(laddertype);
1564       CylinderGeomIntt *mygeom = new CylinderGeomIntt(sphxlayer, params->get_double_param("strip_x"), params->get_double_param("strip_y"), params->get_double_param("strip_z_0"),
1565                                                       params->get_double_param("strip_z_1"), params->get_int_param("nstrips_z_sensor_0"), params->get_int_param("nstrips_z_sensor_1"),
1566                                                       params->get_int_param("nstrips_phi_sensor"), params_layer->get_int_param("nladder"),
1567                                                       m_PosZ[ilayer][0] / cm,  // m_PosZ uses G4 internal units, needs to be converted to cm
1568                                                       m_PosZ[ilayer][1] / cm, m_SensorRadius[ilayer] / cm, 0.0,
1569                                                       params_layer->get_double_param("offsetphi") * deg / rad,  // expects radians
1570                                                       params_layer->get_double_param("offsetrot") * deg / rad);
1571       geo->AddLayerGeom(sphxlayer, mygeom);
1572       if (Verbosity() > 0)
1573       {
1574         geo->identify();
1575       }
1576     }
1577   }
1578 }
1579 
1580 std::map<G4VPhysicalVolume *, std::tuple<int, int, int, int>>::const_iterator PHG4InttDetector::get_ActiveVolumeTuple(G4VPhysicalVolume *physvol) const
1581 {
1582   auto iter = m_ActiveVolumeTuple.find(physvol);
1583   if (iter == m_ActiveVolumeTuple.end())
1584   {
1585     std::cout << PHWHERE << " Volume " << physvol->GetName() << " not in active volume tuple" << std::endl;
1586     gSystem->Exit(1);
1587   }
1588   return iter;
1589 }
1590 
1591 std::map<G4LogicalVolume *, std::tuple<int, int>>::const_iterator PHG4InttDetector::get_PassiveVolumeTuple(G4LogicalVolume *logvol) const
1592 {
1593   auto iter = m_PassiveVolumeTuple.find(logvol);
1594   if (iter == m_PassiveVolumeTuple.end())
1595   {
1596     std::cout << PHWHERE << " Volume " << logvol->GetName() << " not in passive volume tuple" << std::endl;
1597     gSystem->Exit(1);
1598   }
1599   return iter;
1600 }