Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-06 08:22:07

0001 #include "PHG4Prototype3InnerHcalDetector.h"
0002 
0003 #include <phparameter/PHParameters.h>
0004 
0005 #include <g4main/PHG4Detector.h>  // for PHG4Detector
0006 #include <g4main/PHG4Units.h>
0007 
0008 #include <TSystem.h>
0009 
0010 #include <Geant4/G4AssemblyVolume.hh>
0011 #include <Geant4/G4Box.hh>
0012 #include <Geant4/G4Colour.hh>
0013 #include <Geant4/G4ExtrudedSolid.hh>
0014 #include <Geant4/G4LogicalVolume.hh>
0015 #include <Geant4/G4Material.hh>
0016 #include <Geant4/G4PVPlacement.hh>
0017 #include <Geant4/G4RotationMatrix.hh>  // for G4RotationMatrix
0018 #include <Geant4/G4String.hh>          // for G4String
0019 #include <Geant4/G4SystemOfUnits.hh>   // for mm, deg, cm, cm3, rad
0020 #include <Geant4/G4ThreeVector.hh>     // for G4ThreeVector
0021 #include <Geant4/G4TwoVector.hh>
0022 #include <Geant4/G4VPhysicalVolume.hh>  // for G4VPhysicalVolume
0023 #include <Geant4/G4VSolid.hh>           // for G4VSolid
0024 #include <Geant4/G4VisAttributes.hh>
0025 
0026 #include <boost/format.hpp>
0027 
0028 #include <cmath>
0029 #include <iostream>  // for operator<<, endl
0030 #include <sstream>
0031 #include <utility>  // for pair, make_pair
0032 #include <vector>   // for vector, vector<>::...
0033 
0034 class PHCompositeNode;
0035 
0036 using namespace std;
0037 
0038 static const string scintimothername = "InnerHcalScintiMother";
0039 static const string steelplatename = "InnerHcalSteelPlate";
0040 
0041 PHG4Prototype3InnerHcalDetector::PHG4Prototype3InnerHcalDetector(PHG4Subsystem *subsys, PHCompositeNode *Node, PHParameters *parameters, const std::string &dnam)
0042   : PHG4Detector(subsys, Node, dnam)
0043   , m_params(parameters)
0044   , m_InnerHcalSteelPlate(nullptr)
0045   , m_InnerHcalAssembly(nullptr)
0046   , m_scintibox(nullptr)
0047   , m_SteelPlateCornerUpperLeft(1154.49 * mm, -189.06 * mm)
0048   , m_SteelPlateCornerUpperRight(1297.94 * mm, -349.22 * mm)
0049   , m_SteelPlateCornerLowerRight(1288.18 * mm, -357.8 * mm)
0050   , m_SteelPlateCornerLowerLeft(1157.3 * mm, -205.56 * mm)
0051   //
0052   , m_ScintiTile9DistanceToCorner(26.44 * mm)
0053   , m_ScintiTile9FrontSize(140.3 * mm)
0054   , m_ScintiTile9CornerUpperLeft(0 * mm, 0 * mm)
0055   , m_ScintiTile9CornerUpperRight(191. * mm, -134.4 * 191. / 198.1 * mm)
0056   , m_ScintiTile9CornerLowerRight(191. * mm, -191. * mm / tan(52.02 / 180. * M_PI) - m_ScintiTile9FrontSize)
0057   , m_ScintiTile9CornerLowerLeft(0 * mm, -m_ScintiTile9FrontSize)
0058   //
0059   , m_ScintiTile10FrontSize(149.2 * mm)
0060   , m_ScintiTile10CornerUpperLeft(0 * mm, 0 * mm)
0061   , m_ScintiTile10CornerUpperRight(191. * mm, -154.6 * 191. / 198.1 * mm)
0062   , m_ScintiTile10CornerLowerRight(191. * mm, -191. * mm / tan(48.34 / 180. * M_PI) - m_ScintiTile10FrontSize)
0063   , m_ScintiTile10CornerLowerLeft(0 * mm, -m_ScintiTile10FrontSize)
0064   //
0065   , m_ScintiTile11FrontSize(144.3 * mm)
0066   , m_ScintiTile11CornerUpperLeft(0 * mm, 0 * mm)
0067   , m_ScintiTile11CornerUpperRight(191. * mm, -176.2 * 191. / 198.1 * mm)
0068   , m_ScintiTile11CornerLowerRight(191. * mm, -191. * mm / tan(45.14 / 180. * M_PI) - m_ScintiTile11FrontSize)
0069   , m_ScintiTile11CornerLowerLeft(0 * mm, -m_ScintiTile11FrontSize)
0070   //
0071   , m_ScintiTile12FrontSize(186.6 * mm)
0072   , m_ScintiTile12CornerUpperLeft(0 * mm, 0 * mm)
0073   , m_ScintiTile12CornerUpperRight(191. * mm, -197.11 * 191. / 198.1 * mm)
0074   , m_ScintiTile12CornerLowerRight(191. * mm, -191. * mm / tan(41.47 / 180. * M_PI) - m_ScintiTile12FrontSize)
0075   , m_ScintiTile12CornerLowerLeft(0 * mm, -m_ScintiTile12FrontSize)
0076   //
0077   , m_ScintiX(198.1 * mm)
0078   , m_SteelZ(901.7 * mm)
0079   , m_ScintiTileZ(m_SteelZ)
0080   , m_ScintiTileThickness(7 * mm)
0081   , m_GapBetweenTiles(1 * mm)
0082   , m_ScintiGap(10 * mm)
0083   , m_DeltaPhi(2 * M_PI / 256.)
0084   , m_VolumeSteel(NAN)
0085   , m_VolumeScintillator(0)
0086   , m_ScintiCornerLowerLeft(45.573 * inch, -7.383 * inch)
0087   , m_ScintiCornerLowerRight(50.604 * inch, -12.972 * inch)
0088   // leave this in in case we ever need those coordinates
0089   // m_ScintiCornerUpperRight(50.809*inch,-12.787*inch),
0090   // m_ScintiCornerUpperLeft(45.778*inch,-7.198*inch),
0091   , m_NumScintiPlates(16)
0092   , m_NumSteelPlates(m_NumScintiPlates + 1)
0093   , m_Active(m_params->get_int_param("active"))
0094   , m_AbsorberActive(m_params->get_int_param("absorberactive"))
0095   , m_Layer(0)
0096 {
0097   // patch to get the upper steel plate location right within less than a mm
0098   m_DeltaPhi += 0.00125 * m_DeltaPhi;
0099 }
0100 
0101 PHG4Prototype3InnerHcalDetector::~PHG4Prototype3InnerHcalDetector()
0102 {
0103   delete m_InnerHcalAssembly;
0104 }
0105 
0106 //_______________________________________________________________
0107 //_______________________________________________________________
0108 int PHG4Prototype3InnerHcalDetector::IsInPrototype3InnerHcal(G4VPhysicalVolume *volume) const
0109 {
0110   G4LogicalVolume *logvol = volume->GetLogicalVolume();
0111   if (m_AbsorberActive && logvol == m_InnerHcalSteelPlate)
0112   {
0113     return -1;
0114   }
0115   if (m_Active && m_ActiveVolumeSet.find(logvol) != m_ActiveVolumeSet.end())
0116   {
0117     return 1;
0118   }
0119   return 0;
0120 }
0121 
0122 G4LogicalVolume *
0123 PHG4Prototype3InnerHcalDetector::ConstructSteelPlate(G4LogicalVolume */*hcalenvelope*/)
0124 {
0125   if (!m_InnerHcalSteelPlate)
0126   {
0127     G4VSolid *steel_plate;
0128 
0129     std::vector<G4TwoVector> vertexes;
0130     vertexes.push_back(m_SteelPlateCornerUpperLeft);
0131     vertexes.push_back(m_SteelPlateCornerUpperRight);
0132     vertexes.push_back(m_SteelPlateCornerLowerRight);
0133     vertexes.push_back(m_SteelPlateCornerLowerLeft);
0134     G4TwoVector zero(0, 0);
0135     steel_plate = new G4ExtrudedSolid("InnerHcalSteelPlateSolid",
0136                                       vertexes,
0137                                       m_SteelZ / 2.0,
0138                                       zero, 1.0,
0139                                       zero, 1.0);
0140 
0141     m_VolumeSteel = steel_plate->GetCubicVolume() * m_NumSteelPlates;
0142     m_InnerHcalSteelPlate = new G4LogicalVolume(steel_plate, G4Material::GetMaterial(m_params->get_string_param("material")), steelplatename, 0, 0, 0);
0143     G4VisAttributes visattchk;
0144     visattchk.SetVisibility(true);
0145     visattchk.SetForceSolid(false);
0146     visattchk.SetColour(G4Colour::Blue());
0147     m_InnerHcalSteelPlate->SetVisAttributes(visattchk);
0148   }
0149   return m_InnerHcalSteelPlate;
0150 }
0151 
0152 G4LogicalVolume *
0153 PHG4Prototype3InnerHcalDetector::ConstructScintillatorBoxHiEta(G4LogicalVolume *hcalenvelope)
0154 {
0155   int copynum = 0;
0156   G4VSolid *scintiboxsolid = new G4Box(scintimothername, m_ScintiX / 2., (m_ScintiGap) / 2., m_ScintiTileZ / 2.);
0157   //DisplayVolume(scintiboxsolid,hcalenvelope);
0158 
0159   G4LogicalVolume *scintiboxlogical = new G4LogicalVolume(scintiboxsolid, G4Material::GetMaterial("G4_AIR"), G4String(scintimothername), 0, 0, 0);
0160   G4VisAttributes hcalVisAtt;
0161   hcalVisAtt.SetVisibility(true);
0162   hcalVisAtt.SetForceSolid(false);
0163   hcalVisAtt.SetColour(G4Colour::Magenta());
0164   G4LogicalVolume *scintit9_logic = ConstructScintiTile9(hcalenvelope);
0165   scintit9_logic->SetVisAttributes(hcalVisAtt);
0166   double distance_to_corner = -m_SteelZ / 2. + m_ScintiTile9DistanceToCorner;
0167   G4RotationMatrix *Rot;
0168   Rot = new G4RotationMatrix();
0169   Rot->rotateX(90 * deg);
0170   new G4PVPlacement(Rot, G4ThreeVector(-m_ScintiX / 2., 0, distance_to_corner), scintit9_logic, (boost::format("InnerScinti_%d") % copynum).str(), scintiboxlogical, false, copynum, OverlapCheck());
0171 
0172   hcalVisAtt.SetVisibility(true);
0173   hcalVisAtt.SetForceSolid(false);
0174   hcalVisAtt.SetColour(G4Colour::Blue());
0175   G4LogicalVolume *scintit10_logic = ConstructScintiTile10(hcalenvelope);
0176   scintit10_logic->SetVisAttributes(hcalVisAtt);
0177 
0178   distance_to_corner += m_ScintiTile9FrontSize + m_GapBetweenTiles;
0179   Rot = new G4RotationMatrix();
0180   Rot->rotateX(90 * deg);
0181   copynum++;
0182   new G4PVPlacement(Rot, G4ThreeVector(-m_ScintiX / 2., 0, distance_to_corner), scintit10_logic, (boost::format("InnerScinti_%d") % copynum).str(), scintiboxlogical, false, copynum, OverlapCheck());
0183 
0184   hcalVisAtt.SetVisibility(true);
0185   hcalVisAtt.SetForceSolid(false);
0186   hcalVisAtt.SetColour(G4Colour::Yellow());
0187   G4LogicalVolume *scintit11_logic = ConstructScintiTile11(hcalenvelope);
0188   scintit11_logic->SetVisAttributes(hcalVisAtt);
0189 
0190   distance_to_corner += m_ScintiTile10FrontSize + m_GapBetweenTiles;
0191   Rot = new G4RotationMatrix();
0192   Rot->rotateX(90 * deg);
0193   copynum++;
0194   new G4PVPlacement(Rot, G4ThreeVector(-m_ScintiX / 2., 0, distance_to_corner), scintit11_logic, (boost::format("InnerScinti_%d") % copynum).str(), scintiboxlogical, false, copynum, OverlapCheck());
0195 
0196   hcalVisAtt.SetVisibility(true);
0197   hcalVisAtt.SetForceSolid(false);
0198   hcalVisAtt.SetColour(G4Colour::Cyan());
0199   G4LogicalVolume *scintit12_logic = ConstructScintiTile12(hcalenvelope);
0200   scintit12_logic->SetVisAttributes(hcalVisAtt);
0201 
0202   distance_to_corner += m_ScintiTile11FrontSize + m_GapBetweenTiles;
0203   Rot = new G4RotationMatrix();
0204   Rot->rotateX(90 * deg);
0205   copynum++;
0206   new G4PVPlacement(Rot, G4ThreeVector(-m_ScintiX / 2., 0, distance_to_corner), scintit12_logic, (boost::format("InnerScinti_%d") % copynum).str(), scintiboxlogical, false, copynum, OverlapCheck());
0207   //    DisplayVolume(scintiboxlogical,hcalenvelope);
0208   return scintiboxlogical;
0209 }
0210 
0211 G4LogicalVolume *
0212 PHG4Prototype3InnerHcalDetector::ConstructScintiTile9(G4LogicalVolume */*hcalenvelope*/)
0213 {
0214   std::vector<G4TwoVector> vertexes;
0215   vertexes.push_back(m_ScintiTile9CornerUpperLeft);
0216   vertexes.push_back(m_ScintiTile9CornerUpperRight);
0217   vertexes.push_back(m_ScintiTile9CornerLowerRight);
0218   vertexes.push_back(m_ScintiTile9CornerLowerLeft);
0219   G4TwoVector zero(0, 0);
0220   G4VSolid *scintit9 = new G4ExtrudedSolid("InnerHcalScintiT9",
0221                                            vertexes,
0222                                            m_ScintiTileThickness / 2.0,
0223                                            zero, 1.0,
0224                                            zero, 1.0);
0225 
0226   m_VolumeScintillator += scintit9->GetCubicVolume() * m_NumScintiPlates;
0227   G4LogicalVolume *scintit9_logic = new G4LogicalVolume(scintit9, G4Material::GetMaterial("G4_POLYSTYRENE"), "InnerHcalScintiT9", nullptr, nullptr, nullptr);
0228   //     DisplayVolume(scintit9,hcalenvelope);
0229   m_ActiveVolumeSet.insert(scintit9_logic);
0230   return scintit9_logic;
0231 }
0232 
0233 G4LogicalVolume *
0234 PHG4Prototype3InnerHcalDetector::ConstructScintiTile10(G4LogicalVolume */*hcalenvelope*/)
0235 {
0236   std::vector<G4TwoVector> vertexes;
0237   vertexes.push_back(m_ScintiTile10CornerUpperLeft);
0238   vertexes.push_back(m_ScintiTile10CornerUpperRight);
0239   vertexes.push_back(m_ScintiTile10CornerLowerRight);
0240   vertexes.push_back(m_ScintiTile10CornerLowerLeft);
0241   G4TwoVector zero(0, 0);
0242   G4VSolid *scintit10 = new G4ExtrudedSolid("InnerHcalScintiT10",
0243                                             vertexes,
0244                                             m_ScintiTileThickness / 2.0,
0245                                             zero, 1.0,
0246                                             zero, 1.0);
0247 
0248   m_VolumeScintillator += scintit10->GetCubicVolume() * m_NumScintiPlates;
0249   G4LogicalVolume *scintit10_logic = new G4LogicalVolume(scintit10, G4Material::GetMaterial("G4_POLYSTYRENE"), "InnerHcalScintiT10", nullptr, nullptr, nullptr);
0250   //     DisplayVolume(scintit10,hcalenvelope);
0251   m_ActiveVolumeSet.insert(scintit10_logic);
0252   return scintit10_logic;
0253 }
0254 
0255 G4LogicalVolume *
0256 PHG4Prototype3InnerHcalDetector::ConstructScintiTile11(G4LogicalVolume */*hcalenvelope*/)
0257 {
0258   std::vector<G4TwoVector> vertexes;
0259   vertexes.push_back(m_ScintiTile11CornerUpperLeft);
0260   vertexes.push_back(m_ScintiTile11CornerUpperRight);
0261   vertexes.push_back(m_ScintiTile11CornerLowerRight);
0262   vertexes.push_back(m_ScintiTile11CornerLowerLeft);
0263   G4TwoVector zero(0, 0);
0264   G4VSolid *scintit11 = new G4ExtrudedSolid("InnerHcalScintiT11",
0265                                             vertexes,
0266                                             m_ScintiTileThickness / 2.0,
0267                                             zero, 1.0,
0268                                             zero, 1.0);
0269 
0270   m_VolumeScintillator += scintit11->GetCubicVolume() * m_NumScintiPlates;
0271   G4LogicalVolume *scintit11_logic = new G4LogicalVolume(scintit11, G4Material::GetMaterial("G4_POLYSTYRENE"), "InnerHcalScintiT11", nullptr, nullptr, nullptr);
0272   //     DisplayVolume(scintit11,hcalenvelope);
0273   m_ActiveVolumeSet.insert(scintit11_logic);
0274   return scintit11_logic;
0275 }
0276 
0277 G4LogicalVolume *
0278 PHG4Prototype3InnerHcalDetector::ConstructScintiTile12(G4LogicalVolume */*hcalenvelope*/)
0279 {
0280   std::vector<G4TwoVector> vertexes;
0281   vertexes.push_back(m_ScintiTile12CornerUpperLeft);
0282   vertexes.push_back(m_ScintiTile12CornerUpperRight);
0283   vertexes.push_back(m_ScintiTile12CornerLowerRight);
0284   vertexes.push_back(m_ScintiTile12CornerLowerLeft);
0285   G4TwoVector zero(0, 0);
0286   G4VSolid *scintit12 = new G4ExtrudedSolid("InnerHcalScintiT12",
0287                                             vertexes,
0288                                             m_ScintiTileThickness / 2.0,
0289                                             zero, 1.0,
0290                                             zero, 1.0);
0291 
0292   m_VolumeScintillator += scintit12->GetCubicVolume() * m_NumScintiPlates;
0293   G4LogicalVolume *scintit12_logic = new G4LogicalVolume(scintit12, G4Material::GetMaterial("G4_POLYSTYRENE"), "InnerHcalScintiT12", nullptr, nullptr, nullptr);
0294   //       DisplayVolume(scintit12,hcalenvelope);
0295   m_ActiveVolumeSet.insert(scintit12_logic);
0296   return scintit12_logic;
0297 }
0298 
0299 // Construct the envelope and the call the
0300 // actual inner hcal construction
0301 void PHG4Prototype3InnerHcalDetector::ConstructMe(G4LogicalVolume *logicWorld)
0302 {
0303   G4ThreeVector g4vec(m_params->get_double_param("place_x") * cm,
0304                       m_params->get_double_param("place_y") * cm,
0305                       m_params->get_double_param("place_z") * cm);
0306   G4RotationMatrix Rot;
0307   Rot.rotateX(m_params->get_double_param("rot_x") * deg);
0308   Rot.rotateY(m_params->get_double_param("rot_y") * deg);
0309   Rot.rotateZ(m_params->get_double_param("rot_z") * deg);
0310   m_InnerHcalAssembly = new G4AssemblyVolume();
0311   ConstructInnerHcal(logicWorld);
0312   m_InnerHcalAssembly->MakeImprint(logicWorld, g4vec, &Rot, 0, OverlapCheck());
0313   // this is rather pathetic - there is no way to extract the name when a volume is added
0314   // to the assembly. The only thing we can do is get an iterator over the placed volumes
0315   // in the order in which they were placed. Since this code does not install the scintillators
0316   // for the Al version, parsing the volume names to get the id does not work since it changes
0317   // So now we loop over all volumes and store them in a map for fast lookup of the row
0318   int isteel = 0;
0319   int iscinti = 0;
0320   vector<G4VPhysicalVolume *>::iterator it = m_InnerHcalAssembly->GetVolumesIterator();
0321   for (unsigned int i = 0; i < m_InnerHcalAssembly->TotalImprintedVolumes(); i++)
0322   {
0323     string volname = (*it)->GetName();
0324     if (volname.find(steelplatename) != string::npos)
0325     {
0326       m_SteelPlateIdMap.insert(make_pair(volname, isteel));
0327       ++isteel;
0328     }
0329     else if (volname.find(scintimothername) != string::npos)
0330     {
0331       m_ScintillatorIdMap.insert(make_pair(volname, iscinti));
0332       ++iscinti;
0333     }
0334     ++it;
0335   }
0336   // print out volume names and their assigned id
0337   // map<string,int>::const_iterator iter;
0338   // for (iter = m_SteelPlateIdMap.begin(); iter != m_SteelPlateIdMap.end(); ++iter)
0339   // {
0340   //   cout << iter->first << ", " << iter->second << endl;
0341   // }
0342   // for (iter = m_ScintillatorIdMap.begin(); iter != m_ScintillatorIdMap.end(); ++iter)
0343   // {
0344   //   cout << iter->first << ", " << iter->second << endl;
0345   // }
0346   return;
0347 }
0348 
0349 int PHG4Prototype3InnerHcalDetector::ConstructInnerHcal(G4LogicalVolume *hcalenvelope)
0350 {
0351   G4LogicalVolume *steel_plate = ConstructSteelPlate(hcalenvelope);  // bottom steel plate
0352   if (m_params->get_int_param("scintillators"))
0353   {
0354     if (m_params->get_int_param("hi_eta"))
0355     {
0356       m_scintibox = ConstructScintillatorBoxHiEta(hcalenvelope);
0357     }
0358     else
0359     {
0360       cout << "midrapidity scintillator not implemented" << endl;
0361       gSystem->Exit(1);
0362     }
0363   }
0364   double phi = 0.;
0365   double phislat = 0.;
0366   // the coordinate of the center of the bottom of the bottom steel plate
0367   // to get the radius of the circle which is the center of the scintillator box
0368 
0369   double bottom_xmiddle_steel_tile = (m_SteelPlateCornerLowerRight.x() + m_SteelPlateCornerLowerLeft.x()) / 2.;
0370   double bottom_ymiddle_steel_tile = (m_SteelPlateCornerLowerLeft.y() + m_SteelPlateCornerLowerRight.y()) / 2.;
0371   // the math is not exact, need to move the middle radius by 14mm to
0372   // get the upper steel plate right
0373   double middlerad = sqrt(bottom_xmiddle_steel_tile * bottom_xmiddle_steel_tile + bottom_ymiddle_steel_tile * bottom_ymiddle_steel_tile) - 0.14 * cm;
0374   // another fudge factor to get the upper steel location right
0375   double philow = atan((bottom_ymiddle_steel_tile - (m_ScintiGap * 25. / 32.)) / bottom_xmiddle_steel_tile);
0376 
0377   double scintiangle = GetScintiAngle();
0378 
0379   // need to shift in x to get the upper steel plate close to right
0380   double xstart = 0;
0381   double xoff = 0.015 * cm;
0382   for (int i = 0; i < m_NumSteelPlates; i++)
0383   {
0384     G4RotationMatrix Rot;
0385     Rot.rotateZ(phi * rad);
0386     G4ThreeVector g4vec(xstart, 0, 0);
0387     m_InnerHcalAssembly->AddPlacedVolume(steel_plate, g4vec, &Rot);
0388     if (m_scintibox && i > 0)
0389     {
0390       double ypos = sin(phi + philow) * middlerad;
0391       double xpos = cos(phi + philow) * middlerad;
0392       G4RotationMatrix Rot1;
0393       Rot1.rotateZ(scintiangle + phislat);
0394       G4ThreeVector g4vecsc(xpos + xstart, ypos, 0);
0395       m_InnerHcalAssembly->AddPlacedVolume(m_scintibox, g4vecsc, &Rot1);
0396       phislat += m_DeltaPhi;
0397     }
0398     phi += m_DeltaPhi;
0399     xstart += xoff;
0400   }
0401   return 0;
0402 }
0403 
0404 // calculate the angle of the bottom scintillator. It is the angle of the top edge
0405 // of the steel plate
0406 double
0407 PHG4Prototype3InnerHcalDetector::GetScintiAngle()
0408 {
0409   double xlen = m_ScintiCornerLowerRight.x() - m_ScintiCornerLowerLeft.x();
0410   double ylen = m_ScintiCornerLowerRight.y() - m_ScintiCornerLowerLeft.y();
0411   double angle = atan(ylen / xlen);
0412   return angle;
0413 }
0414 
0415 void PHG4Prototype3InnerHcalDetector::Print(const string &what) const
0416 {
0417   cout << "Inner Hcal Detector:" << endl;
0418   if (what == "ALL" || what == "VOLUME")
0419   {
0420     cout << "Volume Steel: " << m_VolumeSteel / cm3 << " cm^3" << endl;
0421     cout << "Volume Scintillator: " << m_VolumeScintillator / cm3 << " cm^3" << endl;
0422   }
0423   return;
0424 }
0425 
0426 int PHG4Prototype3InnerHcalDetector::get_scinti_row_id(const string &volname)
0427 {
0428   int id = -9999;
0429   auto it = m_ScintillatorIdMap.find(volname);
0430   if (it != m_ScintillatorIdMap.end())
0431   {
0432     id = it->second;
0433   }
0434   else
0435   {
0436     cout << "unknown scintillator volume name: " << volname << endl;
0437   }
0438 
0439   return id;
0440 }
0441 
0442 int PHG4Prototype3InnerHcalDetector::get_steel_plate_id(const string &volname)
0443 {
0444   int id = -9999;
0445   auto it = m_SteelPlateIdMap.find(volname);
0446   if (it != m_SteelPlateIdMap.end())
0447   {
0448     id = it->second;
0449   }
0450   else
0451   {
0452     cout << "unknown steel volume name: " << volname << endl;
0453   }
0454   return id;
0455 }