EIC Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4JLeicVTXDetector.cc
Go to the documentation of this file. Or view the newest version in sPHENIX GitHub for file G4JLeicVTXDetector.cc
1 #include "G4JLeicVTXDetector.h"
2 
3 #include <phparameter/PHParameters.h>
4 #include <phparameter/PHParametersContainer.h>
5 
6 #include <g4main/PHG4Detector.h> // for PHG4Detector
7 
8 #include <Geant4/G4Box.hh>
9 #include <Geant4/G4Color.hh>
10 #include <Geant4/G4LogicalVolume.hh>
11 #include <Geant4/G4Material.hh>
12 #include <Geant4/G4PVPlacement.hh>
13 #include <Geant4/G4RotationMatrix.hh> // for G4RotationMatrix
14 #include <Geant4/G4String.hh> // for G4String
15 #include <Geant4/G4SystemOfUnits.hh>
16 #include <Geant4/G4ThreeVector.hh> // for G4ThreeVector
17 #include <Geant4/G4Transform3D.hh> // for G4Transform3D
18 #include <Geant4/G4VisAttributes.hh>
19 
20 #include <cmath>
21 #include <iostream> // for operator<<, endl, bas...
22 #include <utility> // for pair
23 
24 class G4VSolid;
25 class PHCompositeNode;
26 
27 using namespace std;
28 
30  : PHG4Detector(subsys, Node, dnam)
31  , m_ParamsContainer(params)
32 {
34  m_IsActiveFlag = par->get_int_param("active");
35  m_IsAbsorberActiveFlag = par->get_int_param("absorberactive");
36  m_Layers = par->get_int_param("layers");
37 }
38 
39 //_______________________________________________________________
40 //_______________________________________________________________
41 int G4JLeicVTXDetector::IsInVTX(G4VPhysicalVolume *volume) const
42 {
43  map<G4VPhysicalVolume *, int>::const_iterator iter = m_PhysicalVolumesMap.find(volume);
44 
45  if (iter != m_PhysicalVolumesMap.end())
46  {
47  return iter->second;
48  }
49 
50  return 0;
51 }
52 
53 void G4JLeicVTXDetector::ConstructMe(G4LogicalVolume *logicWorld)
54 {
55  for (int ilayer = 0; ilayer < m_Layers; ilayer++)
56  {
57  // get parameters for this layer
58  const PHParameters *par = m_ParamsContainer->GetParameters(ilayer);
59  double cb_VTX_ladder_DZ = par->get_double_param("Dz") * cm;
60  double cb_VTX_ladder_DY = par->get_double_param("Dy") * cm;
61  double cb_VTX_ladder_Thickness = par->get_double_param("Dx") * cm;
62  double dR = par->get_double_param("Rin") * cm;
63  double myL = 2 * M_PI * dR;
64  int NUM = myL / cb_VTX_ladder_DY;
65 
66  for (int i = 0; i < 2; i++)
67  {
68  double LN = cb_VTX_ladder_DY * NUM;
69  double LN1 = cb_VTX_ladder_DY * (NUM + 1 + i);
70  if (LN / LN1 > 0.8) NUM = NUM + 1;
71  }
72 
73  double cb_VTX_ladder_deltaphi = 2 * M_PI / NUM;
74  string solidname = "cb_VTX_ladder_Solid_" + to_string(ilayer);
75  G4VSolid *solid = new G4Box(solidname, cb_VTX_ladder_Thickness / 2., cb_VTX_ladder_DY / 2., cb_VTX_ladder_DZ / 2.);
76  string logical_name = "cb_VTX_ladder_Logic_" + to_string(ilayer);
77  G4LogicalVolume *logical = new G4LogicalVolume(solid, GetDetectorMaterial("G4_Si"), logical_name);
78  G4VisAttributes *attr_cb_VTX_ladder = nullptr;
79  switch (ilayer)
80  {
81  case 0:
82  case 1:
83  attr_cb_VTX_ladder = new G4VisAttributes(G4Color(0.0, 0.2, 0.8, 2.0));
84  break;
85  case 2:
86  attr_cb_VTX_ladder = new G4VisAttributes(G4Color(0.0, 0.2, 0.8, 0.7));
87  break;
88  default:
89  attr_cb_VTX_ladder = new G4VisAttributes(G4Color(0.0 + 0.1 * double(ilayer - 3), 1., 1. - 0.1 * double(ilayer - 3), 1.0));
90  break;
91  }
92  attr_cb_VTX_ladder->SetForceSolid(true);
93  logical->SetVisAttributes(attr_cb_VTX_ladder);
94  for (int ia = 0; ia < NUM; ia++)
95  {
96  double phi = (ia * (cb_VTX_ladder_deltaphi));
97  double x = -dR * cos(phi);
98  double y = -dR * sin(phi);
99  G4RotationMatrix rot;
100  rot.rotateZ(cb_VTX_ladder_deltaphi * ia);
101  rot.rotateZ(-7. * deg);
102  ostringstream physname;
103  physname << "cb_VTX_ladder_Phys_" << ilayer << "_" << ia;
104  G4VPhysicalVolume *phy = new G4PVPlacement(G4Transform3D(rot, G4ThreeVector(x, y, 0)),
105  logical, physname.str(),
106  logicWorld, 0, false, OverlapCheck());
107  // layer starts at zero but needs to be positive to mark active volume
108  m_PhysicalVolumesMap.insert(make_pair(phy, ilayer + 1));
109  }
110  }
111 
112  return;
113 }
114 
115 void G4JLeicVTXDetector::Print(const std::string &what) const
116 {
117  cout << "JLeic VTX Detector:" << endl;
118  if (what == "ALL" || what == "VOLUME")
119  {
120  cout << "Version 0.1" << endl;
121  }
122  return;
123 }