Belle II Software development
QcsmonitorCreator.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/qcsmonitor/geometry/QcsmonitorCreator.h>
10#include <beast/qcsmonitor/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 <G4Box.hh>
27#include <G4UserLimits.hh>
28
29//Visualization Attributes
30#include <G4VisAttributes.hh>
31
32using namespace std;
33using namespace boost;
34
35namespace Belle2 {
42 namespace qcsmonitor {
43
44 // Register the creator
47
49 {
50 //m_sensitive = new SensitiveDetector();
51 }
52
54 {
55 if (m_sensitive) delete m_sensitive;
56 }
57
58 void QcsmonitorCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
59 {
60
62
63 //lets get the stepsize parameter with a default value of 5 µm
64 double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
65
66 G4VisAttributes* red = new G4VisAttributes(G4Colour(1, 0, 0));
67 red->SetForceAuxEdgeVisible(true);
68 G4VisAttributes* green = new G4VisAttributes(G4Colour(0, 1, 0));
69 green->SetForceAuxEdgeVisible(true);
70 G4VisAttributes* gray = new G4VisAttributes(G4Colour(.5, .5, .5));
71 gray->SetForceAuxEdgeVisible(true);
72 G4VisAttributes* coppercolor = new G4VisAttributes(G4Colour(218. / 255., 138. / 255., 103. / 255.));
73 coppercolor->SetForceAuxEdgeVisible(true);
74
75 //Lets loop over all the Active nodes
76 BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
77 int phase = activeParams.getInt("phase");
78 G4double dx_scint = activeParams.getLength("dx_scint") / 2.*CLHEP::cm;
79 G4double dy_scint = activeParams.getLength("dy_scint") / 2.*CLHEP::cm;
80 G4double dz_scint = activeParams.getLength("dz_scint") / 2.*CLHEP::cm;
81 double thetaZ = activeParams.getAngle("ThetaZ");
82 G4VSolid* s_scint = new G4Box("s_scint", dx_scint, dy_scint, dz_scint);
83 G4LogicalVolume* l_scint = new G4LogicalVolume(s_scint, geometry::Materials::get("G4_POLYSTYRENE"), "l_scint", 0, m_sensitive);
84 l_scint->SetVisAttributes(green);
85 //Lets limit the Geant4 stepsize inside the volume
86 l_scint->SetUserLimits(new G4UserLimits(stepSize));
87 double x_pos[100];
88 double y_pos[100];
89 double z_pos[100];
90 double r_pos[100];
91 int dim = 0;
92 if (phase == 1) {
93 dim = 0;
94 for (double x : activeParams.getArray("x", {0})) {
95 x *= CLHEP::cm;
96 x_pos[dim] = x;
97 dim++;
98 }
99 dim = 0;
100 for (double y : activeParams.getArray("y", {0})) {
101 y *= CLHEP::cm;
102 y_pos[dim] = y;
103 dim++;
104 }
105 }
106 dim = 0;
107 for (double z : activeParams.getArray("z", {0})) {
108 z *= CLHEP::cm;
109 z_pos[dim] = z;
110 //cout << "QCSS z " << z << " zpos " << z_pos[dim] << endl;
111 if (phase == 1) {
112 r_pos[dim] = sqrt(x_pos[dim] * x_pos[dim] + y_pos[dim] * y_pos[dim]);
113 }
114 dim++;
115 }
116 //if (phase == 1) dim = 1;
117 if (phase == 2) {
118 dim = 0;
119 for (double r : activeParams.getArray("r", {0})) {
120 r *= CLHEP::cm;
121 r_pos[dim] = r + dz_scint;
122 dim++;
123 }
124 for (int i = 0; i < 100; i++) {
125 x_pos[i] = 0;
126 y_pos[i] = 0;
127 //x_off[i] = 0;
128 //y_off[i] = 0;
129 }
130 }
131 int detID = 0;
132 G4Transform3D transform;
133 for (double phi : activeParams.getArray("Phi", {M_PI / 2})) {
134 //phi *= CLHEP::deg;
135 for (int i = 0; i < dim; i++) {
136 if (phase == 2) {
137 transform = G4RotateZ3D(phi - M_PI / 2) * G4Translate3D(0, r_pos[i], z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ);
138 //cout << "phase 2" << endl;
139 }
140 if (phase == 1) {
141 transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateZ3D(phi) * G4RotateX3D(thetaZ);
142 //cout << "phase 1" << endl;
143 }
144 //cout << "QCS r " << r_pos[i] << " width " << dz_scint << " z " << z_pos[i] << " phi " << phi << " x " << x_pos[i] << " y " <<
145 // y_pos[i] << endl;
146 new G4PVPlacement(transform, l_scint, TString::Format("p_scint_%d", detID).Data(), &topVolume, false, detID);
147 B2INFO("QCSS-" << detID << " placed at: " << transform.getTranslation() << " mm ");
148 //cout << " Nb of detector " << detID << endl;
149 detID++;
150 }
151 }
152 }
153 }
154 } // qcsmonitor namespace
156} // Belle2 namespace
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
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
std::vector< double > getArray(const std::string &path) const noexcept(false)
Get the parameter path as a list of double values converted to the standard unit.
Definition: Interface.cc:123
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
int getInt(const std::string &path="") const noexcept(false)
Get the parameter path as a int.
Definition: Interface.cc:60
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:63
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Creation of the detector geometry from Gearbox (XML).
SensitiveDetector * m_sensitive
SensitiveDetector QCSMONITOR.
Sensitive Detector implementation of the QCSMONITOR detector.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
GeometryTypes
Flag indiciating the type of geometry to be used.
geometry::CreatorFactory< QcsmonitorCreator > QcsmonitorFactory("QCSMONITORCreator")
Creator creates the QCSMONITOR 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.