Belle II Software development
BeamabortCreator.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/beamabort/geometry/BeamabortCreator.h>
10#include <beast/beamabort/simulation/SensitiveDetector.h>
11
12#include <geometry/Materials.h>
13#include <geometry/CreatorFactory.h>
14#include <framework/logging/Logger.h>
15
16#include <cmath>
17
18#include <G4LogicalVolume.hh>
19#include <G4PVPlacement.hh>
20
21//Shapes
22#include <G4Box.hh>
23#include "G4SubtractionSolid.hh"
24#include <G4UserLimits.hh>
25
26//Visualization
27#include "G4Colour.hh"
28#include <G4VisAttributes.hh>
29
30using namespace std;
31
32namespace Belle2 {
38 using namespace geometry;
39
41 namespace beamabort {
42
43 // Register the creator
45 geometry::CreatorFactory <BeamabortCreator> BeamabortFactory("BEAMABORTCreator");
46
48 {}
49
51 {}
52
53 void BeamabortCreator::createGeometry(G4LogicalVolume& topVolume, GeometryTypes)
54 {
55 //get parameters from .xml file
56 std::string prep = "Active.";
57
58 B2INFO("BeamabortCreator phase3");
59 //Visualization Attributes
60 G4VisAttributes* orange = new G4VisAttributes(G4Colour(1, 2, 0));
61 orange->SetForceAuxEdgeVisible(true);
62 G4VisAttributes* magenta = new G4VisAttributes(G4Colour(1, 0, 1));
63 magenta->SetForceAuxEdgeVisible(true);
64 //lets get the stepsize parameter with a default value of 5 µm
65 double stepSize = m_config.getParameter("stepSize", 5 * CLHEP::um);
66
67 //Lets loop over all the Active nodes
68 int phase = m_config.getParameter(prep + "phase");
69 //Positioned Diamonds
70 std::vector<double> x_pos;
71 std::vector<double> y_pos;
72 std::vector<double> z_pos;
73 std::vector<double> phi;
74 std::vector<double> thetaX;
75 std::vector<double> thetaY;
76 std::vector<double> thetaZ;
77 std::vector<double> svdAngle;
78 std::vector<double> r;
79 std::vector<double> deltaX;
80 unsigned int dimz = 0;
81
82 for (double z : m_config.getParArray(prep + "z", {0})) {
83 z *= CLHEP::cm;
84 z_pos.push_back(z);
85 dimz++;
86 }
87
88 for (double ThetaZ : m_config.getParArray(prep + "ThetaZ", {0})) {
89 thetaZ.push_back(ThetaZ);
90 }
91
92 if (thetaZ.size() != dimz) { B2ERROR("Diamond data not consistent (i.e. not same number of all position parameters)"); return;}
93
94 if (phase == 2 || phase == 3) {
95
96 for (double Phi : m_config.getParArray(prep + "Phi", {0})) {
97 phi.push_back(Phi);
98 }
99
100 for (double r_dia : m_config.getParArray(prep + "r_dia", {0})) {
101 r_dia *= CLHEP::cm;
102 r.push_back(r_dia);
103 }
104 for (double dX : m_config.getParArray(prep + "deltaX", {0})) {
105 dX *= CLHEP::cm;
106 deltaX.push_back(dX);
107 }
108 for (double addAngle : m_config.getParArray(prep + "addAngle", {0})) {
109 svdAngle.push_back(addAngle);
110 }
111 if (phi.size() != dimz || r.size() != dimz || deltaX.size() != dimz || svdAngle.size() != dimz) { B2ERROR("Diamond data not consistent (i.e. not same number of all position parameters)"); return;}
112 }
113 if (phase == 1) {
114 for (double x : m_config.getParArray(prep + "x", {0})) {
115 x *= CLHEP::cm;
116 x_pos.push_back(x);
117 }
118 for (double y : m_config.getParArray(prep + "y", {0})) {
119 y *= CLHEP::cm;
120 y_pos.push_back(y);
121 }
122
123 for (double ThetaX : m_config.getParArray(prep + "ThetaX", {0})) {
124 thetaX.push_back(ThetaX);
125 }
126 for (double ThetaY : m_config.getParArray(prep + "ThetaY", {0})) {
127 thetaY.push_back(ThetaY);
128 }
129 if (x_pos.size() != dimz || y_pos.size() != dimz || thetaX.size() != dimz || thetaY.size() != dimz) { B2ERROR("Diamond data not consistent (i.e. not same number of all position parameters)"); return;}
130 }
131
132 //create beamabort package
133 G4double dx_opa = 12. / 2. * CLHEP::mm;
134 G4double dy_opa = 18. / 2. * CLHEP::mm;
135 G4double dz_opa = 3.1 / 2. * CLHEP::mm;
136 G4VSolid* s_pa = new G4Box("s_opa", dx_opa, dy_opa, dz_opa);
137 G4double dx_ipa = 6. / 2. * CLHEP::mm;
138 G4double dy_ipa = 6. / 2. * CLHEP::mm;
139 G4double dz_ipa = 5.6 / 2. * CLHEP::mm;
140 G4VSolid* s_ipa = new G4Box("s_ipa", dx_ipa, dy_ipa, dz_ipa);
141 s_pa = new G4SubtractionSolid("s_pa", s_pa, s_ipa, 0, G4ThreeVector(0., 4.0, 0.));
142 G4LogicalVolume* l_pa = new G4LogicalVolume(s_pa, G4Material::GetMaterial("Al6061"), "l_pa");
143 l_pa->SetVisAttributes(magenta);
144 G4Transform3D transform;
145 for (unsigned int i = 0; i < dimz; i++) {
146 B2INFO("DIA-" << i << "RotateZ3D phi: " << phi[i]);
147
148 if (phase == 1) {
149 transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateX3D(thetaX[i]) *
150 G4RotateY3D(thetaY[i]) * G4RotateZ3D(thetaZ[i]);
151 } else if (phase == 2 || phase == 3) {
152 transform = G4Translate3D(deltaX[i], 0, 0) * G4RotateZ3D(phi[i]) *
153 G4Translate3D(r[i], 0, z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ[i]) *
154 G4RotateZ3D(svdAngle[i]) * G4RotateY3D(M_PI / 2);
155 }
156 new G4PVPlacement(transform, l_pa, TString::Format("p_dia_pa_%d", i).Data(), &topVolume, false, 0);
157 B2INFO("DIA-" << i << " placed at: " << transform.getTranslation() << " mm ");
158 }
159
160 //create beamabort volumes
161 G4double dx_ba = 4.5 / 2. * CLHEP::mm;
162 G4double dy_ba = 4.5 / 2. * CLHEP::mm;
163 G4double dz_ba = 0.5 / 2. * CLHEP::mm;
164 G4Box* s_BEAMABORT = new G4Box("s_BEAMABORT", dx_ba, dy_ba, dz_ba);
165 G4LogicalVolume* l_BEAMABORT = new G4LogicalVolume(s_BEAMABORT, geometry::Materials::get("Diamond"), "l_BEAMABORT");
166 l_BEAMABORT->SetVisAttributes(orange);
167
168 //Lets limit the Geant4 stepsize inside the volume
169 l_BEAMABORT->SetUserLimits(new G4UserLimits(stepSize));
170 l_BEAMABORT->SetVisAttributes(orange);
171
172 // make the volume sensitive if the beam background study flag is enabled
173 if (m_config.getParameter("BeamBackgroundStudy")) l_BEAMABORT->SetSensitiveDetector(new SensitiveDetector());
174
175 for (unsigned int i = 0; i < dimz; i++) {
176
177 if (phase == 1)
178 transform = G4Translate3D(x_pos[i], y_pos[i],
179 z_pos[i]) * G4RotateX3D(thetaX[i]) * G4RotateY3D(thetaY[i]) * G4RotateZ3D(thetaZ[i]);
180 if (phase == 2 || phase == 3)
181 transform = G4Translate3D(deltaX[i], 0, 0) * G4RotateZ3D(phi[i]) *
182 G4Translate3D(r[i], 0, z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ[i]) *
183 G4RotateZ3D(svdAngle[i]) * G4RotateY3D(M_PI / 2) * G4Translate3D(0, 4.2, 0);
184
185 new G4PVPlacement(transform, l_BEAMABORT, TString::Format("p_dia_%d", i).Data(), &topVolume, false, i);
186
187 B2INFO("DIA-sensitive volume-" << i << " placed at: " << transform.getTranslation() << " mm "
188 << " at phi angle = " << phi[i] << " at theta angle = "
189 << thetaZ[i]);
190 B2INFO("DIA-sensitive volume-" << i << " G4RotateZ3D of: phi= " << phi[i] << " G4RotateX3D = "
191 << (-M_PI / 2 - thetaZ[i]));
192 }
193 }
194 }
196}
std::vector< double > getParArray(const std::string &name) const
Get parameter value array.
Definition: BeamabortGeo.h:64
double getParameter(const std::string &name) const
Get parameter value.
Definition: BeamabortGeo.h:41
void createGeometry(G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Create detector geometry.
BeamabortGeo m_config
geometry parameters object
Sensitive Detector implementation of the BEAMABORT detector.
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
geometry::CreatorFactory< BeamabortCreator > BeamabortFactory("BEAMABORTCreator")
Creator creates the BEAMABORT geometry.
GeometryTypes
Flag indicating the type of geometry to be used.
Abstract base class for different kinds of events.
STL namespace.