9#include <beast/qcsmonitor/geometry/QcsmonitorCreator.h>
10#include <beast/qcsmonitor/simulation/SensitiveDetector.h>
12#include <geometry/Materials.h>
13#include <geometry/CreatorFactory.h>
14#include <framework/gearbox/GearDir.h>
15#include <framework/logging/Logger.h>
18#include <boost/format.hpp>
19#include <boost/foreach.hpp>
20#include <boost/algorithm/string.hpp>
22#include <G4LogicalVolume.hh>
23#include <G4PVPlacement.hh>
27#include <G4UserLimits.hh>
30#include <G4VisAttributes.hh>
42 namespace qcsmonitor {
64 double stepSize = content.getLength(
"stepSize", 5 * CLHEP::um);
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);
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);
84 l_scint->SetVisAttributes(green);
86 l_scint->SetUserLimits(
new G4UserLimits(stepSize));
94 for (
double x : activeParams.
getArray(
"x", {0})) {
100 for (
double y : activeParams.
getArray(
"y", {0})) {
107 for (
double z : activeParams.
getArray(
"z", {0})) {
112 r_pos[dim] =
sqrt(x_pos[dim] * x_pos[dim] + y_pos[dim] * y_pos[dim]);
119 for (
double r : activeParams.
getArray(
"r", {0})) {
121 r_pos[dim] = r + dz_scint;
124 for (
int i = 0; i < 100; i++) {
132 G4Transform3D transform;
133 for (
double phi : activeParams.
getArray(
"Phi", {M_PI / 2})) {
135 for (
int i = 0; i < dim; i++) {
137 transform = G4RotateZ3D(phi - M_PI / 2) * G4Translate3D(0, r_pos[i], z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ);
141 transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateZ3D(phi) * G4RotateX3D(thetaZ);
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 ");
GearDir is the basic class used for accessing the parameter store.
double getAngle(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard angle unit.
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.
double getLength(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard length unit.
int getInt(const std::string &path="") const noexcept(false)
Get the parameter path as a int.
static G4Material * get(const std::string &name)
Find given material.
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Creation of the detector geometry from Gearbox (XML).
virtual ~QcsmonitorCreator()
Destructor.
QcsmonitorCreator()
Constructor.
SensitiveDetector * m_sensitive
SensitiveDetector QCSMONITOR.
Sensitive Detector implementation of the QCSMONITOR detector.
double sqrt(double a)
sqrt for double
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.
Very simple class to provide an easy way to register creators with the CreatorManager.