Belle II Software  release-08-01-10
BgoCreator.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <beast/bgo/geometry/BgoCreator.h>
10 #include <beast/bgo/simulation/SensitiveDetector.h>
11 
12 #include <geometry/CreatorFactory.h>
13 #include <framework/gearbox/GearDir.h>
14 #include <framework/logging/Logger.h>
15 
16 #include <boost/format.hpp>
17 #include <boost/foreach.hpp>
18 #include <boost/algorithm/string.hpp>
19 
20 #include <G4LogicalVolume.hh>
21 #include <G4PVPlacement.hh>
22 
23 //Shapes
24 #include <G4Trap.hh>
25 #include <G4UserLimits.hh>
26 #include "G4NistManager.hh"
27 
28 using namespace std;
29 using namespace boost;
30 
31 namespace Belle2 {
38  namespace bgo {
39 
40  // Register the creator
43 
44  BgoCreator::BgoCreator(): m_sensitive(0)
45  {
46  //m_sensitive = new SensitiveDetector();
47  }
48 
50  {
51  if (m_sensitive) delete m_sensitive;
52  }
53 
54  void BgoCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
55  {
56 
58 
59  // **Materials from the NIST database**
60  G4NistManager* man = G4NistManager::Instance();
61 
62  G4bool isotopes = false;
63 
64  G4Element* O = man->FindOrBuildElement("O", isotopes);
65  G4Element* Bi = man->FindOrBuildElement("Bi", isotopes);
66  G4Element* Ge = man->FindOrBuildElement("Ge", isotopes);
67 
68  G4Material* BGO_BGO = new G4Material("BGO_BGO",//name
69  7.13 * CLHEP::g / CLHEP::cm3, //density
70  3);//number of elements
71  BGO_BGO->AddElement(O, 12);
72  BGO_BGO->AddElement(Bi, 4);
73  BGO_BGO->AddElement(Ge, 3);
74 
75  //lets get the stepsize parameter with a default value of 5 µm
76  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
77 
78  //no get the array. Notice that the default framework unit is cm, so the
79  //values will be automatically converted
80  vector<double> bar = content.getArray("bar");
81  B2INFO("Contents of bar: ");
82  BOOST_FOREACH(double value, bar) {
83  B2INFO("value: " << value);
84  }
85  int detID = 0;
86  //Lets loop over all the Active nodes
87  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
88 
89  //create bgo volume
90  G4Trap* s_BGO = new G4Trap("s_BGO",
91  activeParams.getLength("cDz") / 2.*CLHEP::cm,
92  activeParams.getLength("cDtheta"),
93  activeParams.getLength("cDphi"),
94  activeParams.getLength("cDy1") / 2.*CLHEP::cm,
95  activeParams.getLength("cDx2") / 2.*CLHEP::cm,
96  activeParams.getLength("cDx1") / 2.*CLHEP::cm, 0,
97  activeParams.getLength("cDy2") / 2.*CLHEP::cm,
98  activeParams.getLength("cDx4") / 2.*CLHEP::cm,
99  activeParams.getLength("cDx3") / 2.*CLHEP::cm, 0);
100 
101  //G4LogicalVolume* l_BGO = new G4LogicalVolume(s_BGO, geometry::Materials::get("BGO"), "l_BGO", 0, m_sensitive);
102  G4LogicalVolume* l_BGO = new G4LogicalVolume(s_BGO, BGO_BGO, "l_BGO", 0, m_sensitive);
103 
104  //cout << "BGO volume " << s_BGO->GetCubicVolume() / CLHEP::cm / CLHEP::cm / CLHEP::cm
105  //<< " density " << geometry::Materials::get("BGO")->GetDensity() / CLHEP::g * CLHEP::cm * CLHEP::cm * CLHEP::cm << endl;
106 
107  //Lets limit the Geant4 stepsize inside the volume
108  l_BGO->SetUserLimits(new G4UserLimits(stepSize));
109 
110  //position bgo volume
111  /*
112  G4Transform3D theta_init = G4RotateX3D(- activeParams.getLength("cDtheta"));
113  G4Transform3D phi_init = G4RotateZ3D(activeParams.getLength("k_phi_init"));
114  G4Transform3D tilt_z = G4RotateY3D(activeParams.getLength("k_z_TILTED"));
115  G4Transform3D tilt_phi = G4RotateZ3D(activeParams.getLength("k_phi_TILTED"));
116  G4Transform3D position = G4Translate3D(activeParams.getLength("k_zC") * tan(activeParams.getLength("k_z_TILTED")) * CLHEP::cm, 0,
117  activeParams.getLength("k_zC") * CLHEP::cm);
118  G4Transform3D pos_phi = G4RotateZ3D(activeParams.getLength("k_phiC"));
119  G4Transform3D Tr = pos_phi * position * tilt_phi * tilt_z * phi_init * theta_init;
120  //cout << "rotation " << Tr.getRotation() << " translation " << Tr.getTranslation() << endl;
121  */
122  double px = activeParams.getDouble("px");
123  double py = activeParams.getDouble("py");
124  double pz = activeParams.getDouble("pz");
125  double angle = activeParams.getDouble("angle");
126  double rx = activeParams.getDouble("rx");
127  double ry = activeParams.getDouble("ry");
128  double rz = activeParams.getDouble("rz");
129 
130  G4RotationMatrix* pRot = new G4RotationMatrix();
131  pRot->rotate(-angle, G4ThreeVector(rx, ry, rz));
132  //G4Transform3D transform = G4Translate3D(px, py, pz) * G4Rotate3D(-angle, G4ThreeVector(rx, ry, rz));
133  new G4PVPlacement(pRot, G4ThreeVector(px, py, pz), l_BGO, "p_BGO", &topVolume, false, detID);
134  //new G4PVPlacement(transform, l_BGO, "p_BGO", &topVolume, false, detID);
135  B2INFO("BGO-" << detID << " placed at: (" << px << "," << py << "," << pz << ")" << " mm ");
136  //B2INFO("BGO-" << detID << " placed at: " << transform.getTranslation() << " mm ");
137  B2INFO(" rotation of " << -angle << " degree a long (" << rx << "," << ry << "," << rz << ") axis");
138  detID++;
139  }
140  }
141  } // bgo namespace
143 } // Belle2 namespace
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
virtual ~BgoCreator()
Destructor.
Definition: BgoCreator.cc:49
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Creation of the detector geometry from Gearbox (XML).
Definition: BgoCreator.cc:54
SensitiveDetector * m_sensitive
SensitiveDetector BGO.
Definition: BgoCreator.h:46
Sensitive Detector implementation of the BGO detector.
double getDouble(const std::string &path="") const noexcept(false)
Get the parameter path as a double.
Definition: Interface.cc:41
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:259
geometry::CreatorFactory< BgoCreator > BgoFactory("BGOCreator")
Creator creates the BGO geometry.
GeometryTypes
Flag indiciating the type of geometry to be used.
Abstract base class for different kinds of events.
Very simple class to provide an easy way to register creators with the CreatorManager.