EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
PHG4EPDDetector.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file PHG4EPDDetector.cc
1 /* vim: set sw=2 ft=cpp: */
2 
3 #include "PHG4EPDDetector.h"
4 
5 #include "PHG4EPDDisplayAction.h"
6 
7 #include <g4main/PHG4Detector.h>
8 #include <g4main/PHG4DisplayAction.h> // for PHG4DisplayAction
9 #include <g4main/PHG4Subsystem.h>
10 
11 #include <phparameter/PHParameters.h>
12 
13 #include <Geant4/G4ExtrudedSolid.hh>
14 #include <Geant4/G4LogicalVolume.hh>
15 #include <Geant4/G4Material.hh>
16 #include <Geant4/G4PVPlacement.hh>
17 #include <Geant4/G4RotationMatrix.hh>
18 #include <Geant4/G4SystemOfUnits.hh>
19 #include <Geant4/G4ThreeVector.hh>
20 #include <Geant4/G4TwoVector.hh> // for G4TwoVector
21 #include <Geant4/G4VPhysicalVolume.hh> // for G4VPhysicalVolume
22 
23 #include <algorithm> // for max
24 #include <cmath>
25 #include <vector> // for vector
26 
28  PHCompositeNode* node,
30  std::string const& name)
31  : PHG4Detector(subsys, node, name)
32  , m_DisplayAction(dynamic_cast<PHG4EPDDisplayAction*>(subsys->GetDisplayAction()))
33  , m_Params(parameters)
34  , m_ActiveFlag(m_Params->get_int_param("active"))
35  , m_SupportActiveFlag(m_Params->get_int_param("supportactive"))
36 {
37 }
38 
39 void PHG4EPDDetector::ConstructMe(G4LogicalVolume* world)
40 {
41  G4Material* material = GetDetectorMaterial("G4_PLASTIC_SC_VINYLTOLUENE");
42 
43  G4ThreeVector positive(0., 0., m_Params->get_double_param("place_z") * cm);
44  G4ThreeVector negative(0., 0., -m_Params->get_double_param("place_z") * cm);
45 
46  constexpr int32_t ntiles = 31;
47  constexpr int32_t nslices = 12;
48 
49  for (int32_t i = 0; i < ntiles; ++i)
50  {
51  std::string label = "EPD_tile_" + std::to_string(i);
52 
53  G4ExtrudedSolid* block = construct_block(i);
54  G4LogicalVolume* volume = new G4LogicalVolume(block, material, label, 0, 0, 0);
55 
56  GetDisplayAction()->AddVolume(volume, volume->GetName());
57  m_ActiveLogVolSet.insert(volume);
58 
59  for (int32_t k = 0; k < nslices; ++k)
60  {
61  G4RotationMatrix* rotate = new G4RotationMatrix();
62 
63  rotate->rotateZ(k * 2 * M_PI / nslices);
64 
65  m_volumes.emplace(
66  new G4PVPlacement(
67  rotate, positive, volume, label, world,
68  false, 2 * k + 0, OverlapCheck()),
69  module_id_for(i, k, 0));
70 
71  m_volumes.emplace(
72  new G4PVPlacement(
73  rotate, negative, volume, label, world,
74  false, 2 * k + 1, OverlapCheck()),
75  module_id_for(i, k, 1));
76  }
77  }
78 }
79 
80 int PHG4EPDDetector::IsInDetector(G4VPhysicalVolume* volume) const
81 {
82  G4LogicalVolume* mylogvol = volume->GetLogicalVolume();
83  if (m_ActiveFlag)
84  {
85  if (m_ActiveLogVolSet.find(mylogvol) != m_ActiveLogVolSet.end())
86  {
87  return 1;
88  }
89  }
91  {
92  if (m_SupportLogVolSet.find(mylogvol) != m_SupportLogVolSet.end())
93  {
94  return -2;
95  }
96  }
97  return 0;
98 }
99 
100 uint32_t PHG4EPDDetector::module_id_for(int32_t index, int32_t slice,
101  int32_t side)
102 {
103  return (side << 9 | slice << 5 | index) & 0x3FF;
104 }
105 
107 {
108  return m_volumes[volume];
109 }
110 
111 static constexpr double dz = 6.;
112 
113 static constexpr double coordinates[31][6][2] = {
114  {{-20.5, 82.4}, {0.0, 85.1}, {20.5, 82.4}, {9.7, 41.9}, {-9.7, 41.9}, {0.0, 0.0}},
115  {{31.9, 124.9}, {1.9, 128.8}, {0.8, 129.0}, {0.8, 86.6}, {0.0, 0.0}, {20.9, 84.0}},
116  {{-31.9, 124.9}, {-1.9, 128.8}, {-0.8, 129.0}, {-0.8, 86.6}, {0.0, 0.0}, {-20.9, 84.0}},
117  {{43.3, 167.4}, {3.0, 172.7}, {0.8, 173.0}, {0.8, 130.6}, {1.9, 130.5}, {32.3, 126.5}},
118  {{-43.3, 167.4}, {-3.0, 172.7}, {-0.8, 173.0}, {-0.8, 130.6}, {-1.9, 130.5}, {-32.3, 126.5}},
119  {{57.6, 220.8}, {4.4, 227.8}, {0.8, 228.3}, {0.8, 174.6}, {3.1, 174.3}, {43.7, 169.0}},
120  {{-57.6, 220.8}, {-4.4, 227.8}, {-0.8, 228.3}, {-0.8, 174.6}, {-3.1, 174.3}, {-43.7, 169.0}},
121  {{71.9, 274.2}, {5.9, 282.9}, {0.8, 283.6}, {0.8, 229.9}, {4.5, 229.4}, {58.0, 222.4}},
122  {{-71.9, 274.2}, {-5.9, 282.9}, {-0.8, 283.6}, {-0.8, 229.9}, {-4.5, 229.4}, {-58.0, 222.4}},
123  {{86.2, 327.6}, {7.3, 338.0}, {0.8, 338.9}, {0.8, 285.2}, {5.9, 284.6}, {72.3, 275.8}},
124  {{-86.2, 327.6}, {-7.3, 338.0}, {-0.8, 338.9}, {-0.8, 285.2}, {-5.9, 284.6}, {-72.3, 275.8}},
125  {{100.5, 381.1}, {8.3, 394.2}, {0.8, 394.2}, {0.8, 340.5}, {7.3, 339.7}, {86.7, 329.2}},
126  {{-100.5, 381.1}, {-8.3, 394.2}, {-0.8, 394.2}, {-0.8, 340.5}, {-7.3, 339.7}, {-86.7, 329.2}},
127  {{114.9, 434.5}, {8.3, 448.5}, {0.8, 449.5}, {0.8, 395.9}, {8.3, 394.9}, {101.0, 382.7}},
128  {{-114.9, 434.5}, {-8.3, 448.5}, {-0.8, 449.5}, {-0.8, 395.9}, {-8.3, 394.9}, {-101.0, 382.7}},
129  {{129.2, 487.9}, {8.3, 503.8}, {0.8, 504.8}, {0.8, 451.2}, {8.3, 450.2}, {115.3, 436.1}},
130  {{-129.2, 487.9}, {-8.3, 503.8}, {-0.8, 504.8}, {-0.8, 451.2}, {-8.3, 450.2}, {-115.3, 436.1}},
131  {{143.5, 541.3}, {8.3, 559.1}, {0.8, 560.1}, {0.8, 506.5}, {8.3, 505.5}, {129.6, 489.5}},
132  {{-143.5, 541.3}, {-8.3, 559.1}, {-0.8, 560.1}, {-0.8, 506.5}, {-8.3, 505.5}, {-129.6, 489.5}},
133  {{157.8, 594.8}, {8.3, 614.4}, {0.8, 615.4}, {0.8, 561.8}, {8.3, 560.8}, {143.9, 542.9}},
134  {{-157.8, 594.8}, {-8.3, 614.4}, {-0.8, 615.4}, {-0.8, 561.8}, {-8.3, 560.8}, {-143.9, 542.9}},
135  {{172.1, 648.2}, {8.3, 669.7}, {0.8, 670.7}, {0.8, 617.1}, {8.3, 616.1}, {158.2, 596.4}},
136  {{-172.1, 648.2}, {-8.3, 669.7}, {-0.8, 670.7}, {-0.8, 617.1}, {-8.3, 616.1}, {-158.2, 596.4}},
137  {{186.4, 701.6}, {8.3, 725.0}, {0.8, 726.0}, {0.8, 672.4}, {8.3, 671.4}, {172.5, 649.8}},
138  {{-186.4, 701.6}, {-8.3, 725.0}, {-0.8, 726.0}, {-0.8, 672.4}, {-8.3, 671.4}, {-172.5, 649.8}},
139  {{200.7, 755.0}, {8.3, 780.4}, {0.8, 781.3}, {0.8, 727.7}, {8.3, 726.7}, {186.8, 703.2}},
140  {{-200.7, 755.0}, {-8.3, 780.4}, {-0.8, 781.3}, {-0.8, 727.7}, {-8.3, 726.7}, {-186.8, 703.2}},
141  {{215.0, 808.4}, {8.3, 835.7}, {0.8, 836.6}, {0.8, 783.0}, {8.3, 782.0}, {201.2, 756.6}},
142  {{-215.0, 808.4}, {-8.3, 835.7}, {-0.8, 836.6}, {-0.8, 783.0}, {-8.3, 782.0}, {-201.2, 756.6}},
143  {{229.4, 861.8}, {17.5, 890.6}, {0.8, 892.7}, {0.8, 838.3}, {8.3, 837.3}, {215.5, 810.0}},
144  {{-229.4, 861.8}, {-17.5, 890.6}, {-0.8, 892.7}, {-0.8, 838.3}, {-8.3, 837.3}, {-215.5, 810.0}},
145 };
146 
147 G4ExtrudedSolid* PHG4EPDDetector::construct_block(int32_t index)
148 {
149  std::string label("tile_" + std::to_string(index));
150 
151  const double(*coords)[6][2] = &coordinates[index];
152 
153  std::vector<G4TwoVector> vertices;
154 
155  for (int32_t i = 0; i < 6; ++i)
156  {
157  double x = (*coords)[i][0];
158  double y = (*coords)[i][1];
159 
160  if (x == 0. && y == 0.)
161  {
162  continue;
163  }
164  vertices.emplace_back(x * mm, y * mm);
165  }
166 
167  std::vector<G4ExtrudedSolid::ZSection> zsections = {
168  G4ExtrudedSolid::ZSection(-dz, G4TwoVector(), 1.),
169  G4ExtrudedSolid::ZSection(dz, G4TwoVector(), 1.),
170  };
171 
172  return new G4ExtrudedSolid(label, vertices, zsections);
173 }