Class Reference for E1039 Core & Analysis Software
PHG4EMCalDetector.cc
Go to the documentation of this file.
1 #include "PHG4EMCalDetector.h"
2 
3 #include <phparameter/PHParameters.h>
4 
7 
8 #include <g4main/PHG4Detector.h> // for PHG4Detector
9 #include <g4main/PHG4Subsystem.h>
10 
11 #include <Geant4/G4Box.hh>
12 #include <Geant4/G4LogicalVolume.hh>
13 #include <Geant4/G4Material.hh>
14 #include <Geant4/G4PVPlacement.hh>
15 #include <Geant4/G4PhysicalConstants.hh>
16 #include <Geant4/G4RotationMatrix.hh> // for G4RotationMatrix
17 #include <Geant4/G4String.hh> // for G4String
18 #include <Geant4/G4SystemOfUnits.hh>
19 #include <Geant4/G4ThreeVector.hh> // for G4ThreeVector
20 #include <Geant4/G4Transform3D.hh> // for G4Transform3D
21 #include <Geant4/G4Tubs.hh>
22 #include <Geant4/G4Types.hh> // for G4double, G4int
23 #include <Geant4/G4VPhysicalVolume.hh> // for G4VPhysicalVolume
24 
25 #include <TSystem.h>
26 
27 #include <cassert>
28 #include <cmath>
29 #include <fstream>
30 #include <iostream>
31 #include <sstream>
32 #include <utility> // for pair, make_pair
33 
34 class G4VSolid;
35 class PHCompositeNode;
36 
37 using namespace std;
38 
39 //_______________________________________________________________________
40 PHG4EMCalDetector::PHG4EMCalDetector(PHCompositeNode* node, PHParameters* parameters, const std::string& dnam, const int lyr)
41  : PHG4Detector(node, dnam)
42  , m_Params(parameters)
43  , m_scintLogical(nullptr)
44  , m_absorberLogical(nullptr)
45  , m_tower_size_x(m_Params->get_double_param("tower_x")*cm)
46  , m_tower_size_y(m_Params->get_double_param("tower_y")*cm)
47  , m_tower_size_z(m_Params->get_double_param("tower_z")*cm)
48  , m_ntowers_x(m_Params->get_int_param("n_towers_x"))
49  , m_ntowers_y(m_Params->get_int_param("n_towers_y"))
50  , m_npbsc_layers(m_Params->get_int_param("n_layers"))
51  , m_ActiveFlag(m_Params->get_int_param("active"))
52  , m_AbsorberActiveFlag(m_Params->get_int_param("absorberactive"))
53  , m_SuperDetector("NONE")
54  , m_Layer(lyr)
55 {}
56 
57 //_______________________________________________________________________
58 int PHG4EMCalDetector::IsInForwardEcal(G4VPhysicalVolume* volume) const
59 {
60  G4LogicalVolume* mylogvol = volume->GetLogicalVolume();
61  if(m_ActiveFlag > 0)
62  {
63  if(m_scintLogical == mylogvol) return 1;
64  }
65  if(m_AbsorberActiveFlag > 0)
66  {
67  if(m_absorberLogical == mylogvol) return -1;
68  }
69  return 0;
70 }
71 
72 //_______________________________________________________________________
73 void PHG4EMCalDetector::Construct(G4LogicalVolume* logicWorld)
74 {
75  if(Verbosity() > 0)
76  {
77  std::cout << "PHG4EMCalDetector: Begin Construction for " << name << std::endl;
78  }
79 
80  // Overal envelop volume for the entire detector
81  double xLength_enve = m_tower_size_x*m_ntowers_x*cm;
82  double yLength_enve = m_tower_size_y*m_ntowers_y*cm;
83  double zLength_enve = m_tower_size_z*cm;
84  G4VSolid* emcal_enve_solid = new G4Box(name+"_enve_solid", xLength_enve/2., yLength_enve/2., zLength_enve/2.);
85 
86  G4Material* mat_Air = G4Material::GetMaterial("G4_AIR");
87  G4LogicalVolume* emcal_enve_logic = new G4LogicalVolume(emcal_enve_solid, mat_Air, name+"_enve_logic");
88 
89  G4RotationMatrix rot_enve;
90  rot_enve.rotateX(m_Params->get_double_param("rot_x")*rad);
91  rot_enve.rotateY(m_Params->get_double_param("rot_y")*rad);
92  rot_enve.rotateZ(m_Params->get_double_param("rot_z")*rad);
93 
94  G4ThreeVector loc_enve;
95  loc_enve.setX(m_Params->get_double_param("place_x")*cm);
96  loc_enve.setY(m_Params->get_double_param("place_y")*cm);
97  loc_enve.setZ(m_Params->get_double_param("place_z")*cm);
98 
99  new G4PVPlacement(G4Transform3D(rot_enve, loc_enve), emcal_enve_logic, name+"_enve_phys", logicWorld, false, 0, overlapcheck);
100 
101  // Construct a single tower
102  G4LogicalVolume* singleTower_logic = ConstructSingleTower();
103 
104  // Place calorimeter towers within envelope
105  PlaceTower(emcal_enve_logic, singleTower_logic);
106 }
107 
108 G4LogicalVolume* PHG4EMCalDetector::ConstructSingleTower()
109 {
110  if(Verbosity() > 0)
111  {
112  std::cout << "PHG4EMCalDetector: Build logical volume for a single tower ..." << std::endl;
113  }
114 
115  // create logical volume for single tower
116  G4Material* mat_Air = G4Material::GetMaterial("G4_AIR");
117  G4VSolid* singleTower_solid = new G4Box(name+"_singleTower_solid", m_tower_size_x/2., m_tower_size_y/2., m_tower_size_z/2.);
118  G4LogicalVolume* singleTower_logic = new G4LogicalVolume(singleTower_solid, mat_Air, name+"_singleTower_logic");
119 
120  /* create geometry volumes for scintillator and absorber plates to place inside single_tower */
121  // PHENIX EMCal JGL 3/27/2016
122  double thickness_layer = m_tower_size_z/(float)m_npbsc_layers;
123  double thickness_absorber = 1.5*mm;
124  double thickness_scint = 4.0*mm;
125  double thickness_airgap = thickness_layer - thickness_absorber - thickness_scint;
126 
127  G4VSolid* absorber_solid = new G4Box(name+"_absorberplate_solid", m_tower_size_x/2., m_tower_size_y/2., thickness_absorber/2.);
128  G4VSolid* scint_solid = new G4Box(name+"_scintplate_solid", m_tower_size_x/2., m_tower_size_y/2., thickness_scint/2.);
129 
130  G4Material* mat_absorber = G4Material::GetMaterial("G4_Pb");
131  G4Material* mat_scint = G4Material::GetMaterial("G4_POLYSTYRENE");
132  m_absorberLogical = new G4LogicalVolume(absorber_solid, mat_absorber, name+"_absorberplate_logic");
133  m_scintLogical = new G4LogicalVolume(scint_solid, mat_scint, name+"_scintplace_logic");
134 
135  /* place physical volumes for absorber and scintillator plates */
136  double zpos_i = -m_tower_size_z/2. + thickness_absorber/2.;
137  for(int i = 0; i < m_npbsc_layers; ++i)
138  {
139  //Absorber plate first
140  new G4PVPlacement(0, G4ThreeVector(0., 0., zpos_i), m_absorberLogical, TString::Format("%s_absplate_%d", name.c_str(), i).Data(), singleTower_logic, false, i, overlapcheck);
141  //std::cout << " -- creating absorber plate at z = " << zpos_i/cm << std::endl;
142 
143  //Scintilator plate second
144  zpos_i += (thickness_absorber/2. + thickness_airgap/2. + thickness_scint/2.);
145  new G4PVPlacement(0, G4ThreeVector(0., 0., zpos_i), m_scintLogical, TString::Format("%s_sciplate_%d", name.c_str(), i).Data(), singleTower_logic, false, i, overlapcheck);
146  //std::cout << " -- creating scintillator plate at z = " << zpos_i/cm << std::endl;
147 
148  //move to the next layer
149  zpos_i += (thickness_scint/2. + thickness_airgap/2. + thickness_absorber/2.);
150  }
151 
152  return singleTower_logic;
153 }
154 
155 void PHG4EMCalDetector::PlaceTower(G4LogicalVolume* envelope_Logic, G4LogicalVolume* singleTower_logic)
156 {
157  int copyNo = 0;
158  for(int i = 0; i < m_ntowers_x; ++i)
159  {
160  double x_pos = -0.5*m_ntowers_x*m_tower_size_x + (i + 0.5)*m_tower_size_x;
161  for(int j = 0; j < m_ntowers_y; ++j)
162  {
163  double y_pos = -0.5*m_ntowers_y*m_tower_size_y + (j + 0.5)*m_tower_size_y;
165  {
166  std::cout << " -- Placing tower (" << i << ", " << j << ") at x = " << x_pos/cm << " cm, y = " << y_pos/cm << " cm" << std::endl;
167  }
168 
169  new G4PVPlacement(0, G4ThreeVector(x_pos, y_pos, 0.), singleTower_logic, TString::Format("%s_tower_%d_%d", name.c_str(), i, j).Data(), envelope_Logic, false, copyNo++, overlapcheck);
170  }
171  }
172 }
@ VERBOSITY_SOME
Output some useful messages during manual command line running.
Definition: Fun4AllBase.h:39
base class for phenix detector creation
Definition: PHG4Detector.h:14
virtual int Verbosity() const
Definition: PHG4Detector.h:47
bool overlapcheck
Definition: PHG4Detector.h:62
std::string name
Definition: PHG4Detector.h:61
int IsInForwardEcal(G4VPhysicalVolume *) const
virtual void Construct(G4LogicalVolume *world)
construct
PHG4EMCalDetector(PHCompositeNode *Node, PHParameters *parameters, const std::string &dnam="EMCal", const int lyr=0)
constructor
double get_double_param(const std::string &name) const