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 <G4LogicalVolume.hh>
18#include <G4PVPlacement.hh>
19
20//Shapes
21#include <G4UserLimits.hh>
22#include "G4Tubs.hh"
23
24using namespace std;
25
26namespace Belle2 {
33 namespace srsensor {
34
35 // Register the creator
38
39 SddCreator::SddCreator(): m_sensitive(0)
40 {
41 //m_sensitive = new SensitiveDetector();
42 }
43
45 {
46 if (m_sensitive) delete m_sensitive;
47 }
48
49 void SddCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
50 {
51
53
54 //lets get the stepsize parameter with a default value of 5 µm
55 double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
56
57 //no get the array. Notice that the default framework unit is cm, so the
58 //values will be automatically converted
59 vector<double> bar_sdd = content.getArray("bar_sdd");
60 B2INFO("Contents of bar_sdd: ");
61 for (double value : bar_sdd) {
62 B2INFO("value: " << value);
63 }
64 int sddNb = 0;
65 //Lets loop over all the Active nodes
66 for (const GearDir& activeParams : content.getNodes("Active")) {
67
68 //create sdd volume
69 G4double startAngle = 0.*CLHEP::deg;
70 G4double spanningAngle = 360.*CLHEP::deg;
71 G4Tubs* s_SDD = new G4Tubs("s_SDD",
72 activeParams.getLength("sdd_innerRadius")*CLHEP::cm,
73 activeParams.getLength("sdd_outerRadius")*CLHEP::cm,
74 activeParams.getLength("sdd_hz")*CLHEP::cm,
75 startAngle, spanningAngle);
76
77 string matSDD = activeParams.getString("MaterialSDD");
78 G4LogicalVolume* l_SDD = new G4LogicalVolume(s_SDD, geometry::Materials::get(matSDD), "l_SDD", 0, m_sensitive);
79
80 //Lets limit the Geant4 stepsize inside the volume
81 l_SDD->SetUserLimits(new G4UserLimits(stepSize));
82
83 //position sdd volume
84 G4ThreeVector SDDpos = G4ThreeVector(
85 activeParams.getLength("x_sdd") * CLHEP::cm,
86 activeParams.getLength("y_sdd") * CLHEP::cm,
87 activeParams.getLength("z_sdd") * CLHEP::cm
88 );
89
90 G4RotationMatrix* rot_sdd = new G4RotationMatrix();
91 rot_sdd->rotateX(activeParams.getAngle("AngleX"));
92 rot_sdd->rotateY(activeParams.getAngle("AngleY"));
93 rot_sdd->rotateZ(activeParams.getAngle("AngleZ"));
94 //geometry::setColor(*l_SDD, "#006699");
95
96 new G4PVPlacement(rot_sdd, SDDpos, l_SDD, "p_SDD", &topVolume, false, sddNb);
97
98 sddNb++;
99 }
100 }
101 } // srsensor namespace
103} // Belle2 namespace
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
virtual ~SddCreator()
Destructor.
Definition: SddCreator.cc:44
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Creation of the detector geometry from Gearbox (XML).
Definition: SddCreator.cc:49
SensitiveDetector * m_sensitive
SensitiveDetector SDD.
Definition: SddCreator.h:46
Sensitive Detector implementation of the SRSENSOR detector.
GeometryTypes
Flag indicating 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.