Belle II Software  release-05-01-25
BgoCreator.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, FaHui Lin *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <beast/bgo/geometry/BgoCreator.h>
12 #include <beast/bgo/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 <G4Trap.hh>
27 #include <G4UserLimits.hh>
28 #include "G4NistManager.hh"
29 
30 using namespace std;
31 using namespace boost;
32 
33 namespace Belle2 {
40  namespace bgo {
41 
42  // Register the creator
44  geometry::CreatorFactory<BgoCreator> BgoFactory("BGOCreator");
45 
46  BgoCreator::BgoCreator(): m_sensitive(0)
47  {
49  }
50 
51  BgoCreator::~BgoCreator()
52  {
53  if (m_sensitive) delete m_sensitive;
54  }
55 
56  void BgoCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
57  {
58  // **Materials from the NIST database**
59  G4NistManager* man = G4NistManager::Instance();
60 
61  G4bool isotopes = false;
62 
63  G4Element* O = man->FindOrBuildElement("O" , isotopes);
64  G4Element* Bi = man->FindOrBuildElement("Bi", isotopes);
65  G4Element* Ge = man->FindOrBuildElement("Ge", isotopes);
66 
67  G4Material* BGO_BGO = new G4Material("BGO_BGO",//name
68  7.13 * CLHEP::g / CLHEP::cm3, //density
69  3);//number of elements
70  BGO_BGO->AddElement(O, 12);
71  BGO_BGO->AddElement(Bi, 4);
72  BGO_BGO->AddElement(Ge, 3);
73 
74  //lets get the stepsize parameter with a default value of 5 µm
75  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
76 
77  //no get the array. Notice that the default framework unit is cm, so the
78  //values will be automatically converted
79  vector<double> bar = content.getArray("bar");
80  B2INFO("Contents of bar: ");
81  BOOST_FOREACH(double value, bar) {
82  B2INFO("value: " << value);
83  }
84  int detID = 0;
85  //Lets loop over all the Active nodes
86  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
87 
88  //create bgo volume
89  G4Trap* s_BGO = new G4Trap("s_BGO",
90  activeParams.getLength("cDz") / 2.*CLHEP::cm ,
91  activeParams.getLength("cDtheta") ,
92  activeParams.getLength("cDphi") ,
93  activeParams.getLength("cDy1") / 2.*CLHEP::cm ,
94  activeParams.getLength("cDx2") / 2.*CLHEP::cm ,
95  activeParams.getLength("cDx1") / 2.*CLHEP::cm , 0,
96  activeParams.getLength("cDy2") / 2.*CLHEP::cm ,
97  activeParams.getLength("cDx4") / 2.*CLHEP::cm ,
98  activeParams.getLength("cDx3") / 2.*CLHEP::cm , 0);
99 
100  //G4LogicalVolume* l_BGO = new G4LogicalVolume(s_BGO, geometry::Materials::get("BGO"), "l_BGO", 0, m_sensitive);
101  G4LogicalVolume* l_BGO = new G4LogicalVolume(s_BGO, BGO_BGO, "l_BGO", 0, m_sensitive);
102 
103  //cout << "BGO volume " << s_BGO->GetCubicVolume() / CLHEP::cm / CLHEP::cm / CLHEP::cm
104  //<< " density " << geometry::Materials::get("BGO")->GetDensity() / CLHEP::g * CLHEP::cm * CLHEP::cm * CLHEP::cm << endl;
105 
106  //Lets limit the Geant4 stepsize inside the volume
107  l_BGO->SetUserLimits(new G4UserLimits(stepSize));
108 
109  //position bgo volume
110  /*
111  G4Transform3D theta_init = G4RotateX3D(- activeParams.getLength("cDtheta"));
112  G4Transform3D phi_init = G4RotateZ3D(activeParams.getLength("k_phi_init"));
113  G4Transform3D tilt_z = G4RotateY3D(activeParams.getLength("k_z_TILTED"));
114  G4Transform3D tilt_phi = G4RotateZ3D(activeParams.getLength("k_phi_TILTED"));
115  G4Transform3D position = G4Translate3D(activeParams.getLength("k_zC") * tan(activeParams.getLength("k_z_TILTED")) * CLHEP::cm, 0,
116  activeParams.getLength("k_zC") * CLHEP::cm);
117  G4Transform3D pos_phi = G4RotateZ3D(activeParams.getLength("k_phiC"));
118  G4Transform3D Tr = pos_phi * position * tilt_phi * tilt_z * phi_init * theta_init;
119  //cout << "rotation " << Tr.getRotation() << " translation " << Tr.getTranslation() << endl;
120  */
121  double px = activeParams.getDouble("px");
122  double py = activeParams.getDouble("py");
123  double pz = activeParams.getDouble("pz");
124  double angle = activeParams.getDouble("angle");
125  double rx = activeParams.getDouble("rx");
126  double ry = activeParams.getDouble("ry");
127  double rz = activeParams.getDouble("rz");
128 
129  G4RotationMatrix* pRot = new G4RotationMatrix();
130  pRot->rotate(-angle, G4ThreeVector(rx, ry, rz));
131  //G4Transform3D transform = G4Translate3D(px, py, pz) * G4Rotate3D(-angle, G4ThreeVector(rx, ry, rz));
132  new G4PVPlacement(pRot, G4ThreeVector(px, py, pz), l_BGO, "p_BGO", &topVolume, false, detID);
133  //new G4PVPlacement(transform, l_BGO, "p_BGO", &topVolume, false, detID);
134  B2INFO("BGO-" << detID << " placed at: (" << px << "," << py << "," << pz << ")" << " mm ");
135  //B2INFO("BGO-" << detID << " placed at: " << transform.getTranslation() << " mm ");
136  B2INFO(" rotation of " << -angle << " degree a long (" << rx << "," << ry << "," << rz << ") axis");
137  detID++;
138  }
139  }
140  } // bgo namespace
142 } // Belle2 namespace
Belle2::bgo::BgoCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector BGO.
Definition: BgoCreator.h:36
Belle2::gearbox::Interface::getDouble
double getDouble(const std::string &path="") const noexcept(false)
Get the parameter path as a double.
Definition: Interface.cc:51
Belle2::bgo::BgoFactory
geometry::CreatorFactory< BgoCreator > BgoFactory("BGOCreator")
Creator creates the BGO 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::bgo::BgoCreator::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: BgoCreator.cc:56
Belle2::geometry::GeometryTypes
GeometryTypes
Flag indiciating the type of geometry to be used.
Definition: GeometryManager.h:39