Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:18:05

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