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
19#include <G4LogicalVolume.hh>
20#include <G4PVPlacement.hh>
21
22//Shapes
23#include <G4Box.hh>
24#include <G4UserLimits.hh>
25
26//Visualization Attributes
27#include <G4VisAttributes.hh>
28
29using namespace std;
30
31namespace Belle2 {
38 namespace qcsmonitor {
39
40 // Register the creator
43
45 {
46 //m_sensitive = new SensitiveDetector();
47 }
48
50 {
51 if (m_sensitive) delete m_sensitive;
52 }
53
54 void QcsmonitorCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
55 {
56
58
59 //lets get the stepsize parameter with a default value of 5 µm
60 double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
61
62 G4VisAttributes* red = new G4VisAttributes(G4Colour(1, 0, 0));
63 red->SetForceAuxEdgeVisible(true);
64 G4VisAttributes* green = new G4VisAttributes(G4Colour(0, 1, 0));
65 green->SetForceAuxEdgeVisible(true);
66 G4VisAttributes* gray = new G4VisAttributes(G4Colour(.5, .5, .5));
67 gray->SetForceAuxEdgeVisible(true);
68 G4VisAttributes* coppercolor = new G4VisAttributes(G4Colour(218. / 255., 138. / 255., 103. / 255.));
69 coppercolor->SetForceAuxEdgeVisible(true);
70
71 //Lets loop over all the Active nodes
72 for (const GearDir& activeParams : content.getNodes("Active")) {
73 int phase = activeParams.getInt("phase");
74 G4double dx_scint = activeParams.getLength("dx_scint") / 2.*CLHEP::cm;
75 G4double dy_scint = activeParams.getLength("dy_scint") / 2.*CLHEP::cm;
76 G4double dz_scint = activeParams.getLength("dz_scint") / 2.*CLHEP::cm;
77 double thetaZ = activeParams.getAngle("ThetaZ");
78 G4VSolid* s_scint = new G4Box("s_scint", dx_scint, dy_scint, dz_scint);
79 G4LogicalVolume* l_scint = new G4LogicalVolume(s_scint, geometry::Materials::get("G4_POLYSTYRENE"), "l_scint", 0, m_sensitive);
80 l_scint->SetVisAttributes(green);
81 //Lets limit the Geant4 stepsize inside the volume
82 l_scint->SetUserLimits(new G4UserLimits(stepSize));
83 double x_pos[100];
84 double y_pos[100];
85 double z_pos[100];
86 double r_pos[100];
87 int dim = 0;
88 if (phase == 1) {
89 dim = 0;
90 for (double x : activeParams.getArray("x", {0})) {
91 x *= CLHEP::cm;
92 x_pos[dim] = x;
93 dim++;
94 }
95 dim = 0;
96 for (double y : activeParams.getArray("y", {0})) {
97 y *= CLHEP::cm;
98 y_pos[dim] = y;
99 dim++;
100 }
101 }
102 dim = 0;
103 for (double z : activeParams.getArray("z", {0})) {
104 z *= CLHEP::cm;
105 z_pos[dim] = z;
106 //cout << "QCSS z " << z << " zpos " << z_pos[dim] << endl;
107 if (phase == 1) {
108 r_pos[dim] = sqrt(x_pos[dim] * x_pos[dim] + y_pos[dim] * y_pos[dim]);
109 }
110 dim++;
111 }
112 //if (phase == 1) dim = 1;
113 if (phase == 2) {
114 dim = 0;
115 for (double r : activeParams.getArray("r", {0})) {
116 r *= CLHEP::cm;
117 r_pos[dim] = r + dz_scint;
118 dim++;
119 }
120 for (int i = 0; i < 100; i++) {
121 x_pos[i] = 0;
122 y_pos[i] = 0;
123 //x_off[i] = 0;
124 //y_off[i] = 0;
125 }
126 }
127 int detID = 0;
128 G4Transform3D transform;
129 for (double phi : activeParams.getArray("Phi", {M_PI / 2})) {
130 //phi *= CLHEP::deg;
131 for (int i = 0; i < dim; i++) {
132 if (phase == 2) {
133 transform = G4RotateZ3D(phi - M_PI / 2) * G4Translate3D(0, r_pos[i], z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ);
134 //cout << "phase 2" << endl;
135 }
136 if (phase == 1) {
137 transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateZ3D(phi) * G4RotateX3D(thetaZ);
138 //cout << "phase 1" << endl;
139 }
140 //cout << "QCS r " << r_pos[i] << " width " << dz_scint << " z " << z_pos[i] << " phi " << phi << " x " << x_pos[i] << " y " <<
141 // y_pos[i] << endl;
142 new G4PVPlacement(transform, l_scint, TString::Format("p_scint_%d", detID).Data(), &topVolume, false, detID);
143 B2INFO("QCSS-" << detID << " placed at: " << transform.getTranslation() << " mm ");
144 //cout << " Nb of detector " << detID << endl;
145 detID++;
146 }
147 }
148 }
149 }
150 } // qcsmonitor namespace
152} // 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 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.