Belle II Software development
SddCreator.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/srsensor/geometry/SddCreator.h>
10#include <beast/srsensor/simulation/SensitiveDetector.h>
11
12#include <geometry/Materials.h>
13#include <geometry/CreatorFactory.h>
14#include <framework/gearbox/GearDir.h>
15#include <framework/logging/Logger.h>
16
17#include <cmath>
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 <G4UserLimits.hh>
27#include "G4Tubs.hh"
28
29using namespace std;
30using namespace boost;
31
32namespace Belle2 {
39 namespace srsensor {
40
41 // Register the creator
44
45 SddCreator::SddCreator(): m_sensitive(0)
46 {
47 //m_sensitive = new SensitiveDetector();
48 }
49
51 {
52 if (m_sensitive) delete m_sensitive;
53 }
54
55 void SddCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
56 {
57
59
60 //lets get the stepsize parameter with a default value of 5 µm
61 double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
62
63 //no get the array. Notice that the default framework unit is cm, so the
64 //values will be automatically converted
65 vector<double> bar_sdd = content.getArray("bar_sdd");
66 B2INFO("Contents of bar_sdd: ");
67 BOOST_FOREACH(double value, bar_sdd) {
68 B2INFO("value: " << value);
69 }
70 int sddNb = 0;
71 //Lets loop over all the Active nodes
72 BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
73
74 //create sdd volume
75 G4double startAngle = 0.*CLHEP::deg;
76 G4double spanningAngle = 360.*CLHEP::deg;
77 G4Tubs* s_SDD = new G4Tubs("s_SDD",
78 activeParams.getLength("sdd_innerRadius")*CLHEP::cm,
79 activeParams.getLength("sdd_outerRadius")*CLHEP::cm,
80 activeParams.getLength("sdd_hz")*CLHEP::cm,
81 startAngle, spanningAngle);
82
83 string matSDD = activeParams.getString("MaterialSDD");
84 G4LogicalVolume* l_SDD = new G4LogicalVolume(s_SDD, geometry::Materials::get(matSDD), "l_SDD", 0, m_sensitive);
85
86 //Lets limit the Geant4 stepsize inside the volume
87 l_SDD->SetUserLimits(new G4UserLimits(stepSize));
88
89 //position sdd volume
90 G4ThreeVector SDDpos = G4ThreeVector(
91 activeParams.getLength("x_sdd") * CLHEP::cm,
92 activeParams.getLength("y_sdd") * CLHEP::cm,
93 activeParams.getLength("z_sdd") * CLHEP::cm
94 );
95
96 G4RotationMatrix* rot_sdd = new G4RotationMatrix();
97 rot_sdd->rotateX(activeParams.getAngle("AngleX"));
98 rot_sdd->rotateY(activeParams.getAngle("AngleY"));
99 rot_sdd->rotateZ(activeParams.getAngle("AngleZ"));
100 //geometry::setColor(*l_SDD, "#006699");
101
102 new G4PVPlacement(rot_sdd, SDDpos, l_SDD, "p_SDD", &topVolume, false, sddNb);
103
104 sddNb++;
105 }
106 }
107 } // srsensor namespace
109} // Belle2 namespace
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
virtual std::string getString(const std::string &path="") const noexcept(false) override
Get the parameter path as a string.
Definition: GearDir.h:69
double getAngle(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard angle unit.
Definition: Interface.h:299
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
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
virtual ~SddCreator()
Destructor.
Definition: SddCreator.cc:50
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Creation of the detector geometry from Gearbox (XML).
Definition: SddCreator.cc:55
SensitiveDetector * m_sensitive
SensitiveDetector SDD.
Definition: SddCreator.h:46
Sensitive Detector implementation of the SRSENSOR detector.
GeometryTypes
Flag indiciating the type of geometry to be used.
geometry::CreatorFactory< SddCreator > SddFactory("SDDCreator")
Creator creates the SDD geometry.
Abstract base class for different kinds of events.
STL namespace.
Very simple class to provide an easy way to register creators with the CreatorManager.