EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4MvtxDetector.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4MvtxDetector.cc
1 #include "PHG4MvtxDetector.h"
2 
3 #include "PHG4MvtxDefs.h"
5 
7 
9 
10 #include <phparameter/PHParameters.h>
11 #include <phparameter/PHParametersContainer.h>
12 
13 #include <g4main/PHG4Detector.h> // for PHG4Detector
14 #include <g4main/PHG4DisplayAction.h> // for PHG4DisplayAction
15 #include <g4main/PHG4Subsystem.h> // for PHG4Subsystem
16 
17 #include <phool/PHCompositeNode.h>
18 #include <phool/PHIODataNode.h>
19 #include <phool/PHNode.h> // for PHNode
20 #include <phool/PHNodeIterator.h> // for PHNodeIterator
21 #include <phool/PHObject.h> // for PHObject
22 #include <phool/getClass.h>
23 
24 #include <Geant4/G4AssemblyVolume.hh>
25 #include <Geant4/G4GDMLParser.hh>
26 #include <Geant4/G4GDMLReadStructure.hh> // for G4GDMLReadStructure
27 #include <Geant4/G4LogicalVolume.hh>
28 #include <Geant4/G4Material.hh>
29 #include <Geant4/G4RotationMatrix.hh> // for G4RotationMatrix
30 #include <Geant4/G4String.hh> // for G4String
31 #include <Geant4/G4SystemOfUnits.hh>
32 #include <Geant4/G4ThreeVector.hh> // for G4ThreeVector
33 #include <Geant4/G4Transform3D.hh> // for G4Transform3D
34 #include <Geant4/G4Types.hh> // for G4double
35 #include <Geant4/G4VPhysicalVolume.hh> // for G4VPhysicalVolume
36 #include <Geant4/G4PVPlacement.hh>
37 #include <Geant4/G4Tubs.hh>
38 
39 #include <cmath>
40 #include <cstdio> // for sprintf
41 #include <iostream> // for operator<<, basic...
42 #include <memory>
43 #include <sstream>
44 #include <utility> // for pair, make_pair
45 #include <vector> // for vector, vector<>:...
46 
47 using namespace std;
48 
49 namespace mvtxGeomDef
50 {
51  double mvtx_shell_inner_radius = 4.8 * cm;
52  double skin_thickness = 0.01 * cm;
53  double foam_core_thickness = 0.18 * cm;
54  double mvtx_shell_length = 50. * cm;
56 
57  double wrap_rmin = 2.1 * cm;
60 };
61 
62 PHG4MvtxDetector::PHG4MvtxDetector(PHG4Subsystem* subsys, PHCompositeNode* Node, const PHParametersContainer* _paramsContainer, const std::string& dnam)
63  : PHG4Detector(subsys, Node, dnam)
64  , m_DisplayAction(dynamic_cast<PHG4MvtxDisplayAction*>(subsys->GetDisplayAction()))
65  , m_ParamsContainer(_paramsContainer)
66  , m_StaveGeometryFile(_paramsContainer->GetParameters(PHG4MvtxDefs::GLOBAL)->get_string_param("stave_geometry_file"))
67  , m_EndWheelsSideS(_paramsContainer->GetParameters(PHG4MvtxDefs::GLOBAL)->get_string_param("end_wheels_sideS"))
68  , m_EndWheelsSideN(_paramsContainer->GetParameters(PHG4MvtxDefs::GLOBAL)->get_string_param("end_wheels_sideN"))
69 {
70  // Verbosity(2);
71 
72  if (Verbosity() > 0)
73  cout << "PHG4MvtxDetector constructor called" << endl;
74 
75  if (Verbosity() > 10)
76  cout << " cm " << cm << " mm " << mm << endl;
77  for (int ilayer = 0; ilayer < n_Layers; ++ilayer)
78  {
79  const PHParameters* params = m_ParamsContainer->GetParameters(ilayer);
80  m_IsLayerActive[ilayer] = params->get_int_param("active");
81  m_IsLayerAbsorberActive[ilayer] = params->get_int_param("absorberactive");
82  m_IsBlackHole[ilayer] = params->get_int_param("blackhole");
83  m_N_staves[ilayer] = params->get_int_param("N_staves");
84  m_nominal_radius[ilayer] = params->get_double_param("layer_nominal_radius");
85  m_nominal_phitilt[ilayer] = params->get_double_param("phitilt");
86  m_nominal_phi0[ilayer] = params->get_double_param("phi0");
87  }
88  /*
89  const PHParameters* alpide_params = m_ParamsContainer->GetParameters(PHG4MvtxDefs::ALPIDE_SEGMENTATION);
90  m_PixelX = alpide_params->get_double_param("pixel_x");
91  m_PixelZ = alpide_params->get_double_param("pixel_z");
92  m_PixelThickness = alpide_params->get_double_param("pixel_thickness");
93  */
94  if (Verbosity() > 0)
95  {
96  cout << "PHG4MvtxDetector constructor: making Mvtx detector. " << endl;
97  }
98 }
99 
100 //_______________________________________________________________
101 //_______________________________________________________________
102 int PHG4MvtxDetector::IsSensor(G4VPhysicalVolume* volume) const
103 {
104  // Is this volume one of the sensors?
105  // Checks if pointer matches one of our stored sensors for this layer
106  if (m_SensorPV.find(volume) != m_SensorPV.end())
107  {
108  if (Verbosity() > 0)
109  {
110  cout << " -- PHG4MvtxTDetector::IsSensor --" << endl;
111  cout << " volume Name : " << volume->GetName() << endl;
112  cout << " -----------------------------------------" << endl;
113  }
114  return 1;
115  }
116 
117  return 0;
118 }
119 
120 int PHG4MvtxDetector::IsInMvtx(G4VPhysicalVolume* volume, int& layer, int& stave) const
121 {
122  // Does this stave belong to this layer?
123  // Since the Assembly volume read from GDML does not give unique pointers
124  // to sensors, we need to check the stave, which is unique
125  auto iter = m_StavePV.find(volume);
126  if (iter != m_StavePV.end())
127  {
128  tie(layer, stave) = iter->second;
129  if (Verbosity() > 0)
130  {
131  cout << " -- PHG4MvtxDetector::IsInMvtx --" << endl;
132  cout << " layer: " << layer << endl;
133  cout << " stave: " << stave << endl;
134  cout << " volume Name : " << volume->GetName() << endl;
135  cout << " stave Name : " << iter->first->GetName() << endl;
136  cout << " -----------------------------------------" << endl;
137  }
138  return 1;
139  }
140 
141  return 0;
142 }
143 
144 int PHG4MvtxDetector::get_layer(int index) const
145 {
146  // Get Mvtx layer from stave index in the Mvtx
147  // Mvtx stave index start from 0, i.e index = 0 for stave 0 in layer 0
148  int lay = 0;
149  while (!(index < m_N_staves[lay]))
150  {
151  index -= m_N_staves[lay];
152  lay++;
153  }
154  return lay;
155 }
156 
157 int PHG4MvtxDetector::get_stave(int index) const
158 {
159  // Get stave index in the layer from stave index in the Mvtx
160  int lay = 0;
161  while (!(index < m_N_staves[lay]))
162  {
163  index -= m_N_staves[lay];
164  lay++;
165  }
166  return index;
167 }
168 
169 void PHG4MvtxDetector::ConstructMe(G4LogicalVolume* logicWorld)
170 {
171  // This is called from PHG4PhenixDetector::Construct()
172 
173  if (Verbosity() > 0)
174  {
175  cout << endl
176  << "PHG4MvtxDetector::Construct called for Mvtx " << endl;
177  }
178 
179  //Create a wrapper volume
180  auto tube = new G4Tubs("sol_MVTX_Wrapper", mvtxGeomDef::wrap_rmin, mvtxGeomDef::wrap_rmax,
181  mvtxGeomDef::wrap_zlen / 2.0, -M_PI, 2.0 * M_PI);
182  auto world_mat = logicWorld->GetMaterial();
183  auto logicMVTX = new G4LogicalVolume(tube, world_mat, "log_MVTX_Wrapper");
184  new G4PVPlacement(new G4RotationMatrix(), G4ThreeVector(), logicMVTX, "MVTX_Wrapper", logicWorld, false, 0, false);
185 
186  // the tracking layers are placed directly in the world volume, since some layers are (touching) double layers
187  // this reads in the ITS stave geometry from a file and constructs the layer from it
188  ConstructMvtx(logicMVTX);
189  ConstructMvtxPassiveVol(logicMVTX);
190 
191  AddGeometryNode();
192  return;
193 }
194 
195 int PHG4MvtxDetector::ConstructMvtx(G4LogicalVolume* trackerenvelope)
196 {
197  if (Verbosity() > 0)
198  {
199  cout << " PHG4MvtxDetector::ConstructMvtx:" << endl;
200  cout << endl;
201  }
202  //===================================
203  // Import the stave physical volume here
204  //===================================
205 
206  // import the staves from the gemetry file
207  std::unique_ptr<G4GDMLReadStructure> reader(new G4GDMLReadStructure());
208  G4GDMLParser gdmlParser(reader.get());
209  gdmlParser.Read(m_StaveGeometryFile, false);
210 
211  // figure out which assembly we want
212  char assemblyname[500];
213  sprintf(assemblyname, "MVTXStave");
214 
215  if (Verbosity() > 0)
216  {
217  cout << "Geting the stave assembly named " << assemblyname << endl;
218  }
219  G4AssemblyVolume* av_ITSUStave = reader->GetAssembly(assemblyname);
220 
221  for (unsigned short ilayer = 0; ilayer < n_Layers; ++ilayer)
222  {
223  if (m_IsLayerActive[ilayer])
224  {
225  if (Verbosity() > 0)
226  {
227  cout << endl;
228  cout << " Constructing Layer " << ilayer << endl;
229  }
230  ConstructMvtx_Layer(ilayer, av_ITSUStave, trackerenvelope);
231  }
232  }
233  FillPVArray(av_ITSUStave);
234  SetDisplayProperty(av_ITSUStave);
235 
236  return 0;
237 }
238 
239 int PHG4MvtxDetector::ConstructMvtx_Layer(int layer, G4AssemblyVolume* av_ITSUStave, G4LogicalVolume*& trackerenvelope)
240 {
241  //=========================================
242  // Now we populate the whole layer with the staves
243  //==========================================
244 
245  int N_staves = m_N_staves[layer];
246  G4double layer_nominal_radius = m_nominal_radius[layer];
247  G4double phitilt = m_nominal_phitilt[layer];
248  G4double phi0 = m_nominal_phi0[layer]; //YCM: azimuthal offset for the first stave
249 
250  if (N_staves < 0)
251  {
252  // The user did not specify how many staves to use for this layer, so we calculate the best value
253  // We choose a phistep that fits N_staves into this radius, but with an arclength separation of AT LEAST arcstep
254  // ideally, the radius would be chosen so that numstaves = N_staves exactly, to get the closest spacing of staves possible without overlaps
255  double arcstep = 12.25;
256  double numstaves = 2.0 * M_PI * layer_nominal_radius / arcstep; // this is just to print out
257  N_staves = int(2.0 * M_PI * layer_nominal_radius / arcstep); // this is the number of staves used
258 
259  // Also suggest the ideal radius for this layer
260  if (Verbosity() > 0)
261  {
262  cout << " Calculated N_staves for layer " /*<< layer*/
263  << " layer_nominal_radius " << layer_nominal_radius
264  << " ITS arcstep " << arcstep
265  << " circumference divided by arcstep " << numstaves
266  << " N_staves " << N_staves
267  << endl;
268  cout << "A radius for this layer of " << (double) N_staves * arcstep / (2.0 * M_PI) + 0.01 << " or "
269  << (double) (N_staves + 1) * arcstep / (2.0 * M_PI) + 0.01 << " would produce perfect stave spacing" << endl;
270  }
271  }
272 
273  G4double phistep = get_phistep(layer); // this produces even stave spacing
274  double z_location = 0.0;
275 
276  if (Verbosity() > 0)
277  {
278  cout << " layer " /*<< layer*/
279  << " layer_nominal_radius " << layer_nominal_radius
280  << " N_staves " << N_staves
281  << " phistep " << phistep
282  << " phitilt " << phitilt
283  << " phi0 " << phi0
284  << endl;
285  }
286 
287  // The stave starts out at (0,0,0) and oriented so that the sensors face upward in y
288  // So we need to rotate the sensor 90 degrees before placing it using phi_offset
289  // note that the gdml file uses a negative phi_offset - different coord system, apparently - the following works
290  double phi_offset = M_PI / 2.0;
291 
292  for (int iphi = 0; iphi < N_staves; iphi++)
293  {
294  // Place the ladder segment envelopes at the correct z and phi
295  // This is the azimuthal angle at which we place the stave
296  // phi0 is the azimuthal offset for the first stave
297  G4double phi_rotation = phi0 + (double) iphi * phistep;
298 
299  G4RotationMatrix Ra;
300  G4ThreeVector Ta;
301 
302  if (Verbosity() > 0)
303  {
304  cout << "phi_offset = " << phi_offset << " iphi " << iphi << " phi_rotation = " << phi_rotation << " phitilt " << phitilt << endl;
305  }
306  // It is first rotated in phi by the azimuthal angle phi_rotation, plus the 90 degrees needed to point the face of the sensor at the origin, plus the tilt (if a tilt is appropriate)
307 
308  // note - if this is layer 0-2, phitilt is the additional tilt for clearance. Otherwise it is zero
309  Ra.rotateZ(phi_rotation + phi_offset + phitilt);
310  // Then translated as follows
311 
312  Ta.setX(layer_nominal_radius * cos(phi_rotation));
313  Ta.setY(layer_nominal_radius * sin(phi_rotation));
314  Ta.setZ(z_location);
315 
316  if (Verbosity() > 0)
317  {
318  cout << " iphi " << iphi << " phi_rotation " << phi_rotation
319  << " x " << layer_nominal_radius * cos(phi_rotation)
320  << " y " << layer_nominal_radius * sin(phi_rotation)
321  << " z " << z_location
322  << endl;
323  }
324  G4Transform3D Tr(Ra, Ta);
325 
326  av_ITSUStave->MakeImprint(trackerenvelope, Tr, 0, OverlapCheck());
327  }
328 
329  if (Verbosity() > 0)
330  {
331  cout << "This layer has a total of " << N_staves << " staves" << endl;
332  }
333  return 0;
334 }
335 
337 {
338  if (Verbosity() > 0)
339  {
340  cout << " PHG4MvtxDetector::ConstructMvtxServices:" << endl;
341  cout << endl;
342  }
343 
344  //=======================================================
345  // Add an outer shell for the MVTX - moved it from INTT PHG4InttDetector.cc
346  //=======================================================
347  G4LogicalVolume *mvtx_shell_outer_skin_volume = GetMvtxOuterShell(lv);
348  new G4PVPlacement(0, G4ThreeVector(0, 0.0), mvtx_shell_outer_skin_volume,
349  "mvtx_shell_outer_skin_volume", lv, false, 0, OverlapCheck());
350 
351  //===================================
352  // Construct Services geometry
353  //===================================
354  if ( (! m_EndWheelsSideN.empty()) && (!m_EndWheelsSideS.empty()) )
355  {
356  // import the end_wheels from the geometry file
357  std::unique_ptr<G4GDMLReadStructure> reader(new G4GDMLReadStructure());
358  G4GDMLParser gdmlParser(reader.get());
359 
360  G4RotationMatrix Ra;
361  Ra.rotateY(M_PI);
362  G4ThreeVector Ta;
363 
364  Ta.setZ(-30);
365  G4Transform3D TrS(Ra, Ta);
366 
367  gdmlParser.Read(m_EndWheelsSideS, false);
368  G4AssemblyVolume* av_EW_S = reader->GetAssembly("EndWheelsSideA");
369  av_EW_S->MakeImprint(lv, TrS, 0, OverlapCheck());
370 
371  Ta.setZ(20);
372  G4Transform3D TrN(Ra, Ta);
373 
374  gdmlParser.Read(m_EndWheelsSideN, false);
375  G4AssemblyVolume* av_EW_N = reader->GetAssembly("EndWheelsSideC");
376  av_EW_N->MakeImprint(lv, TrN, 0, OverlapCheck());
377  }
378  return 0;
379 }
380 
381 G4LogicalVolume* PHG4MvtxDetector::GetMvtxOuterShell(G4LogicalVolume*& trackerenvelope)
382 {
383  // A Rohacell foam sandwich made of 0.1 mm thick CFRP skin and 1.8 mm Rohacell 110 foam core, it has a density of 110 kg/m**3.
384  //mvtx_outer_shell
385  G4Tubs* mvtx_outer_shell_tube = new G4Tubs("mvtx_outer_shell",
388  mvtxGeomDef::mvtx_shell_length / 2.0, -M_PI, 2.0 * M_PI);
389 
390  G4LogicalVolume *mvtx_outer_shell_volume = new G4LogicalVolume(mvtx_outer_shell_tube,
391  trackerenvelope->GetMaterial(),
392  "mvtx_outer_shell_volume", 0, 0, 0);
393 
394 
395  double mvtx_shell_inner_skin_inner_radius = mvtxGeomDef::mvtx_shell_inner_radius;
396  double mvtx_shell_foam_core_inner_radius = mvtx_shell_inner_skin_inner_radius + mvtxGeomDef::skin_thickness;
397  double mvtx_shell_outer_skin_inner_radius = mvtx_shell_foam_core_inner_radius + mvtxGeomDef::foam_core_thickness;
398 
399  G4Tubs *mvtx_shell_inner_skin_tube = new G4Tubs("mvtx_shell_inner_skin",
400  mvtx_shell_inner_skin_inner_radius,
401  mvtx_shell_inner_skin_inner_radius + mvtxGeomDef::skin_thickness,
402  mvtxGeomDef::mvtx_shell_length / 2.0, -M_PI, 2.0 * M_PI);
403 
404  G4LogicalVolume *mvtx_shell_inner_skin_volume = new G4LogicalVolume(mvtx_shell_inner_skin_tube,
405  GetDetectorMaterial("CFRP_INTT"),
406  "mvtx_shell_inner_skin_volume", 0, 0, 0);
407 
408  new G4PVPlacement(0, G4ThreeVector(0, 0.0), mvtx_shell_inner_skin_volume,
409  "mvtx_shell_inner_skin", mvtx_outer_shell_volume, false, 0, OverlapCheck());
410  //m_DisplayAction->AddVolume(mvtx_shell_inner_skin_volume, "Rail");
411 
412  G4Tubs *mvtx_shell_foam_core_tube = new G4Tubs("mvtx_shell_foam_core",
413  mvtx_shell_foam_core_inner_radius,
414  mvtx_shell_foam_core_inner_radius + mvtxGeomDef::foam_core_thickness,
415  mvtxGeomDef::mvtx_shell_length / 2.0, -M_PI, 2.0 * M_PI);
416 
417  G4LogicalVolume *mvtx_shell_foam_core_volume = new G4LogicalVolume(mvtx_shell_foam_core_tube,
418  GetDetectorMaterial("ROHACELL_FOAM_110"),
419  "mvtx_shell_foam_core_volume", 0, 0, 0);
420 
421  new G4PVPlacement(0, G4ThreeVector(0, 0.0), mvtx_shell_foam_core_volume,
422  "mvtx_shell_foam_core", mvtx_outer_shell_volume, false, 0, OverlapCheck());
423  //m_DisplayAction->AddVolume(mvtx_shell_foam_core_volume, "Rail");
424 
425  G4Tubs *mvtx_shell_outer_skin_tube = new G4Tubs("mvtx_shell_outer_skin",
426  mvtx_shell_outer_skin_inner_radius,
427  mvtx_shell_outer_skin_inner_radius + mvtxGeomDef::skin_thickness,
428  mvtxGeomDef::mvtx_shell_length / 2.0, -M_PI, 2.0 * M_PI);
429 
430  G4LogicalVolume *mvtx_shell_outer_skin_volume = new G4LogicalVolume(mvtx_shell_outer_skin_tube,
431  GetDetectorMaterial("CFRP_INTT"),
432  "mvtx_shell_outer_skin_volume", 0, 0, 0);
433 
434  new G4PVPlacement(0, G4ThreeVector(0, 0.0), mvtx_shell_outer_skin_volume,
435  "mvtx_shell_outer_skin", mvtx_outer_shell_volume, false, 0, OverlapCheck());
436  //m_DisplayAction->AddVolume(mvtx_shell_outer_skin_volume, "");
437 
438  return mvtx_outer_shell_volume;
439 }
440 
441 
442 void PHG4MvtxDetector::SetDisplayProperty(G4AssemblyVolume* av)
443 {
444  // cout <<"SetDisplayProperty - G4AssemblyVolume w/ TotalImprintedVolumes "<<av->TotalImprintedVolumes()
445  // <<"/"<<av->GetImprintsCount()<<endl;
446 
447  std::vector<G4VPhysicalVolume*>::iterator it = av->GetVolumesIterator();
448 
449  int nDaughters = av->TotalImprintedVolumes();
450  for (int i = 0; i < nDaughters; ++i, ++it)
451  {
452  // cout <<"SetDisplayProperty - AV["<<i<<"] = "<<(*it)->GetName()<<endl;
453  G4VPhysicalVolume* pv = (*it);
454 
455  G4LogicalVolume* worldLogical = pv->GetLogicalVolume();
456  SetDisplayProperty(worldLogical);
457  }
458 }
459 
460 void PHG4MvtxDetector::SetDisplayProperty(G4LogicalVolume* lv)
461 {
462  string material_name(lv->GetMaterial()->GetName());
463 
464  if (Verbosity() >= 50)
465  {
466  cout << "SetDisplayProperty - LV " << lv->GetName() << " built with "
467  << material_name << endl;
468  }
469  vector<string> matname = {"SI", "KAPTON", "ALUMINUM", "Carbon", "M60J3K", "WATER"};
470  bool found = false;
471  for (string nam : matname)
472  {
473  if (material_name.find(nam) != std::string::npos)
474  {
475  m_DisplayAction->AddVolume(lv, nam);
476  if (Verbosity() >= 50)
477  {
478  cout << "SetDisplayProperty - LV " << lv->GetName() << " display with " << nam << endl;
479  }
480  found = true;
481  break;
482  }
483  }
484  if (!found)
485  {
486  m_DisplayAction->AddVolume(lv, "ANYTHING_ELSE");
487  }
488  int nDaughters = lv->GetNoDaughters();
489  for (int i = 0; i < nDaughters; ++i)
490  {
491  G4VPhysicalVolume* pv = lv->GetDaughter(i);
492 
493  // cout <<"SetDisplayProperty - PV["<<i<<"] = "<<pv->GetName()<<endl;
494 
495  G4LogicalVolume* worldLogical = pv->GetLogicalVolume();
496  SetDisplayProperty(worldLogical);
497  }
498 }
499 
501 {
502  int active = 0;
503  for (auto& isAct : m_IsLayerActive)
504  {
505  active |= isAct;
506  }
507  if (active) // At least one layer is active
508  {
509  ostringstream geonode;
510  geonode << "CYLINDERGEOM_" << ((m_SuperDetector != "NONE") ? m_SuperDetector : m_Detector);
511  PHG4CylinderGeomContainer* geo = findNode::getClass<PHG4CylinderGeomContainer>(topNode(), geonode.str().c_str());
512  if (!geo)
513  {
514  geo = new PHG4CylinderGeomContainer();
515  PHNodeIterator iter(topNode());
516  PHCompositeNode* runNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "RUN"));
517  PHIODataNode<PHObject>* newNode = new PHIODataNode<PHObject>(geo, geonode.str().c_str(), "PHObject");
518  runNode->addNode(newNode);
519  }
520  // here in the detector class we have internal units(mm), convert to cm
521  // before putting into the geom object
522  for (unsigned short ilayer = 0; ilayer < n_Layers; ++ilayer)
523  {
524  CylinderGeom_Mvtx* mygeom = new CylinderGeom_Mvtx(ilayer,
525  m_N_staves[ilayer],
526  m_nominal_radius[ilayer] / cm,
527  get_phistep(ilayer) / rad,
528  m_nominal_phitilt[ilayer] / rad,
529  m_nominal_phi0[ilayer] / rad
530  );
531  geo->AddLayerGeom(ilayer, mygeom);
532  } // loop per layers
533  if (Verbosity())
534  {
535  geo->identify();
536  }
537  } //is active
538 } // AddGeometryNode
539 
540 void PHG4MvtxDetector::FillPVArray(G4AssemblyVolume* av)
541 {
542  if (Verbosity() > 0)
543  {
544  cout << "-- FillPVArray --" << endl;
545  }
546  std::vector<G4VPhysicalVolume*>::iterator it = av->GetVolumesIterator();
547 
548  int nDaughters = av->TotalImprintedVolumes();
549  for (int i = 0; i < nDaughters; ++i, ++it)
550  {
551  G4VPhysicalVolume* pv = (*it);
552 
553  G4LogicalVolume* worldLogical = pv->GetLogicalVolume();
554  // we only care about the staves, which contain the sensors, not the structures
555  if (pv->GetName().find("MVTXHalfStave_pv") != string::npos)
556  {
557  int layer = get_layer(m_StavePV.size());
558  int stave = get_stave(m_StavePV.size());
559 
560  m_StavePV.insert(make_pair(pv, make_tuple(layer, stave)));
561 
562  if (Verbosity() > 0)
563  {
564  cout << "Mvtx layer id " << layer << endl;
565  cout << "Stave in layer id " << stave << endl;
566  cout << "Mvtx stave count " << m_StavePV.size() << endl;
567  cout << "FillPVArray - AV[" << i << "] = " << (*it)->GetName() << endl;
568  cout << " LV[" << i << "] = " << worldLogical->GetName() << endl;
569  }
570 
571  FindSensor(worldLogical);
572  }
573  else // in case of stave structure
574  {
575  if (Verbosity() > 0)
576  {
577  cout << "FillPVArray - AV[" << i << "] = " << (*it)->GetName() << endl;
578  cout << " LV[" << i << "] = " << worldLogical->GetName() << endl;
579  }
580  }
581  }
582 }
583 
584 void PHG4MvtxDetector::FindSensor(G4LogicalVolume* lv)
585 {
586  int nDaughters = lv->GetNoDaughters();
587  for (int i = 0; i < nDaughters; ++i)
588  {
589  G4VPhysicalVolume* pv = lv->GetDaughter(i);
590  if (Verbosity() > 0)
591  {
592  cout << " PV[" << i << "]: " << pv->GetName() << endl;
593  }
594  if (pv->GetName().find("MVTXSensor_") != string::npos)
595  {
596  m_SensorPV.insert(pv);
597 
598  if (Verbosity() > 0)
599  {
600  cout << " Adding Sensor Vol <" << pv->GetName() << " (" << m_SensorPV.size() << ")>" << endl;
601  }
602  }
603 
604  G4LogicalVolume* worldLogical = pv->GetLogicalVolume();
605 
606  if (Verbosity() > 0)
607  {
608  cout << " LV[" << i << "]: " << worldLogical->GetName() << endl;
609  }
610  FindSensor(worldLogical);
611  }
612 }