Belle II Software  release-05-01-25
CaveCreator.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter, Igal Jaegle *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <beast/cave/geometry/CaveCreator.h>
12 #include <beast/cave/simulation/SensitiveDetector.h>
13 
14 #include <geometry/CreatorFactory.h>
15 #include <framework/gearbox/GearDir.h>
16 #include <framework/logging/Logger.h>
17 
18 #include <boost/format.hpp>
19 #include <boost/foreach.hpp>
20 #include <boost/algorithm/string.hpp>
21 
22 #include <G4LogicalVolume.hh>
23 #include <G4PVPlacement.hh>
24 
25 //Shapes
26 #include <G4Box.hh>
27 #include <G4UserLimits.hh>
28 
29 using namespace std;
30 using namespace boost;
31 
32 namespace Belle2 {
39  namespace cave {
40 
41  // Register the creator
43  geometry::CreatorFactory<CaveCreator> CaveFactory("CAVECreator");
44 
45  CaveCreator::CaveCreator(): m_sensitive(0)
46  {
48  }
49 
50  CaveCreator::~CaveCreator()
51  {
52  if (m_sensitive) delete m_sensitive;
53  }
54 
55  void CaveCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
56  {
57  G4double density;
58  G4double A;
59  G4int Z;
60 
61  G4String name, symbol;
62  G4double fractionmass;
63 
64  A = 1.01 * CLHEP::g / CLHEP::mole;
65  G4Element* elH = new G4Element(name = "Hydrogen", symbol = "H" , Z = 1, A);
66 
67  A = 12.01 * CLHEP::g / CLHEP::mole;
68  G4Element* elC = new G4Element(name = "Carbon" , symbol = "C" , Z = 6, A);
69 
70  A = 16.00 * CLHEP::g / CLHEP::mole;
71  G4Element* elO = new G4Element(name = "Oxygen" , symbol = "O" , Z = 8, A);
72 
73  A = 22.99 * CLHEP::g / CLHEP::mole;
74  G4Element* elNa = new G4Element(name = "Natrium" , symbol = "Na" , Z = 11 , A);
75 
76  A = 200.59 * CLHEP::g / CLHEP::mole;
77  G4Element* elHg = new G4Element(name = "Hg" , symbol = "Hg" , Z = 80, A);
78 
79  A = 26.98 * CLHEP::g / CLHEP::mole;
80  G4Element* elAl = new G4Element(name = "Aluminium" , symbol = "Al" , Z = 13, A);
81 
82  A = 28.09 * CLHEP::g / CLHEP::mole;
83  G4Element* elSi = new G4Element(name = "Silicon", symbol = "Si", Z = 14, A);
84 
85  A = 39.1 * CLHEP::g / CLHEP::mole;
86  G4Element* elK = new G4Element(name = "K" , symbol = "K" , Z = 19 , A);
87 
88  A = 69.72 * CLHEP::g / CLHEP::mole;
89  G4Element* elCa = new G4Element(name = "Calzium" , symbol = "Ca" , Z = 31 , A);
90 
91  A = 55.85 * CLHEP::g / CLHEP::mole;
92  G4Element* elFe = new G4Element(name = "Iron" , symbol = "Fe", Z = 26, A);
93 
94  density = 2.03 * CLHEP::g / CLHEP::cm3;
95  G4Material* Concrete = new G4Material("Concrete", density, 10);
96  Concrete->AddElement(elH , fractionmass = 0.01);
97  Concrete->AddElement(elO , fractionmass = 0.529);
98  Concrete->AddElement(elNa , fractionmass = 0.016);
99  Concrete->AddElement(elHg , fractionmass = 0.002);
100  Concrete->AddElement(elAl , fractionmass = 0.034);
101  Concrete->AddElement(elSi , fractionmass = 0.337);
102  Concrete->AddElement(elK , fractionmass = 0.013);
103  Concrete->AddElement(elCa , fractionmass = 0.044);
104  Concrete->AddElement(elFe , fractionmass = 0.014);
105  Concrete->AddElement(elC , fractionmass = 0.001);
106 
107 
108  //lets get the stepsize parameter with a default value of 5 µm
109  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
110 
111  //no get the array. Notice that the default framework unit is cm, so the
112  //values will be automatically converted
113  vector<double> bar = content.getArray("bar");
114  B2INFO("Contents of bar: ");
115  BOOST_FOREACH(double value, bar) {
116  B2INFO("value: " << value);
117  }
118  int detID = 0;
119  //Lets loop over all the Active nodes
120  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
121 
122  //create cave volume
123  G4Box* s_CAVE = new G4Box("s_CAVE",
124  activeParams.getLength("px")*CLHEP::cm ,
125  activeParams.getLength("py")*CLHEP::cm ,
126  activeParams.getLength("pz")*CLHEP::cm);
127 
128  //G4LogicalVolume* l_CAVE = new G4LogicalVolume(s_CAVE, geometry::Materials::get("CAVE"), "l_CAVE", 0, m_sensitive);
129  G4LogicalVolume* l_CAVE = new G4LogicalVolume(s_CAVE, Concrete, "l_CAVE", 0, 0);
130 
131  //Lets limit the Geant4 stepsize inside the volume
132  l_CAVE->SetUserLimits(new G4UserLimits(stepSize));
133 
134  //position cave volume
135  G4ThreeVector CAVEpos = G4ThreeVector(
136  activeParams.getLength("x_cave") * CLHEP::cm,
137  activeParams.getLength("y_cave") * CLHEP::cm,
138  activeParams.getLength("z_cave") * CLHEP::cm
139  );
140 
141  G4RotationMatrix* rot_cave = new G4RotationMatrix();
142  rot_cave->rotateX(activeParams.getLength("AngleX"));
143  rot_cave->rotateY(activeParams.getLength("AngleY"));
144  rot_cave->rotateZ(activeParams.getLength("AngleZ"));
145  //geometry::setColor(*l_CAVE, "#006699");
146  //double angle = activeParams.getDouble("angle");
147  //double rx = activeParams.getDouble("rx");
148  //double ry = activeParams.getDouble("ry");
149  //double rz = activeParams.getDouble("rz");
150  //rot_cave->rotate(-angle, G4ThreeVector(rx, ry, rz));
151  new G4PVPlacement(rot_cave, CAVEpos, l_CAVE, "p_CAVE", &topVolume, false, detID);
152 
153  detID++;
154  }
155  }
156  } // cave namespace
158 } // Belle2 namespace
Belle2::cave::CaveCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector cave.
Definition: CaveCreator.h:36
Belle2::cave::CaveFactory
geometry::CreatorFactory< CaveCreator > CaveFactory("CAVECreator")
Creator creates the cave geometry.
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::GearDir
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:41
Belle2::gearbox::Interface::getLength
double getLength(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard length unit.
Definition: Interface.h:261
Belle2::PXD::SensitiveDetector
VXD::SensitiveDetector< PXDSimHit, PXDTrueHit > SensitiveDetector
The PXD Sensitive Detector class.
Definition: SensitiveDetector.h:36
Belle2::geometry::GeometryTypes
GeometryTypes
Flag indiciating the type of geometry to be used.
Definition: GeometryManager.h:39
Belle2::cave::CaveCreator::create
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Function to actually create the geometry, has to be overridden by derived classes.
Definition: CaveCreator.cc:55