Belle II Software  release-05-01-25
QcsmonitorCreator.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter, Igal Jaegle *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <beast/qcsmonitor/geometry/QcsmonitorCreator.h>
12 #include <beast/qcsmonitor/simulation/SensitiveDetector.h>
13 
14 #include <geometry/Materials.h>
15 #include <geometry/CreatorFactory.h>
16 #include <framework/gearbox/GearDir.h>
17 #include <framework/logging/Logger.h>
18 
19 #include <cmath>
20 #include <boost/format.hpp>
21 #include <boost/foreach.hpp>
22 #include <boost/algorithm/string.hpp>
23 
24 #include <G4LogicalVolume.hh>
25 #include <G4PVPlacement.hh>
26 
27 //Shapes
28 #include <G4Box.hh>
29 #include <G4UserLimits.hh>
30 
31 //Visualization Attributes
32 #include <G4VisAttributes.hh>
33 
34 using namespace std;
35 using namespace boost;
36 
37 namespace Belle2 {
44  namespace qcsmonitor {
45 
46  // Register the creator
48  geometry::CreatorFactory<QcsmonitorCreator> QcsmonitorFactory("QCSMONITORCreator");
49 
50  QcsmonitorCreator::QcsmonitorCreator(): m_sensitive(0)
51  {
53  }
54 
55  QcsmonitorCreator::~QcsmonitorCreator()
56  {
57  if (m_sensitive) delete m_sensitive;
58  }
59 
60  void QcsmonitorCreator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
61  {
62  //lets get the stepsize parameter with a default value of 5 µm
63  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
64 
65  G4VisAttributes* red = new G4VisAttributes(G4Colour(1, 0, 0));
66  red->SetForceAuxEdgeVisible(true);
67  G4VisAttributes* green = new G4VisAttributes(G4Colour(0, 1, 0));
68  green->SetForceAuxEdgeVisible(true);
69  G4VisAttributes* gray = new G4VisAttributes(G4Colour(.5, .5, .5));
70  gray->SetForceAuxEdgeVisible(true);
71  G4VisAttributes* coppercolor = new G4VisAttributes(G4Colour(218. / 255., 138. / 255., 103. / 255.));
72  coppercolor->SetForceAuxEdgeVisible(true);
73 
74  //Lets loop over all the Active nodes
75  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
76  int phase = activeParams.getInt("phase");
77  G4double dx_scint = activeParams.getLength("dx_scint") / 2.*CLHEP::cm;
78  G4double dy_scint = activeParams.getLength("dy_scint") / 2.*CLHEP::cm;
79  G4double dz_scint = activeParams.getLength("dz_scint") / 2.*CLHEP::cm;
80  double thetaZ = activeParams.getAngle("ThetaZ");
81  G4VSolid* s_scint = new G4Box("s_scint", dx_scint, dy_scint, dz_scint);
82  G4LogicalVolume* l_scint = new G4LogicalVolume(s_scint, geometry::Materials::get("G4_POLYSTYRENE"), "l_scint", 0, m_sensitive);
83  l_scint->SetVisAttributes(green);
84  //Lets limit the Geant4 stepsize inside the volume
85  l_scint->SetUserLimits(new G4UserLimits(stepSize));
86  double x_pos[100];
87  double y_pos[100];
88  double z_pos[100];
89  double r_pos[100];
90  int dim = 0;
91  if (phase == 1) {
92  dim = 0;
93  for (double x : activeParams.getArray("x", {0})) {
94  x *= CLHEP::cm;
95  x_pos[dim] = x;
96  dim++;
97  }
98  dim = 0;
99  for (double y : activeParams.getArray("y", {0})) {
100  y *= CLHEP::cm;
101  y_pos[dim] = y;
102  dim++;
103  }
104  }
105  dim = 0;
106  for (double z : activeParams.getArray("z", {0})) {
107  z *= CLHEP::cm;
108  z_pos[dim] = z;
109  //cout << "QCSS z " << z << " zpos " << z_pos[dim] << endl;
110  if (phase == 1) {
111  r_pos[dim] = sqrt(x_pos[dim] * x_pos[dim] + y_pos[dim] * y_pos[dim]);
112  }
113  dim++;
114  }
115  //if (phase == 1) dim = 1;
116  if (phase == 2) {
117  dim = 0;
118  for (double r : activeParams.getArray("r", {0})) {
119  r *= CLHEP::cm;
120  r_pos[dim] = r + dz_scint;
121  dim++;
122  }
123  for (int i = 0; i < 100; i++) {
124  x_pos[i] = 0;
125  y_pos[i] = 0;
126  //x_off[i] = 0;
127  //y_off[i] = 0;
128  }
129  }
130  int detID = 0;
131  G4Transform3D transform;
132  for (double phi : activeParams.getArray("Phi", {M_PI / 2})) {
133  //phi *= CLHEP::deg;
134  for (int i = 0; i < dim; i++) {
135  if (phase == 2) {
136  transform = G4RotateZ3D(phi - M_PI / 2) * G4Translate3D(0, r_pos[i], z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ);
137  //cout << "phase 2" << endl;
138  }
139  if (phase == 1) {
140  transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateZ3D(phi) * G4RotateX3D(thetaZ);
141  //cout << "phase 1" << endl;
142  }
143  //cout << "QCS r " << r_pos[i] << " width " << dz_scint << " z " << z_pos[i] << " phi " << phi << " x " << x_pos[i] << " y " <<
144  // y_pos[i] << endl;
145  new G4PVPlacement(transform, l_scint, TString::Format("p_scint_%d", detID).Data() , &topVolume, false, detID);
146  B2INFO("QCSS-" << detID << " placed at: " << transform.getTranslation() << " mm ");
147  //cout << " Nb of detector " << detID << endl;
148  detID++;
149  }
150  }
151  }
152  }
153  } // qcsmonitor namespace
155 } // Belle2 namespace
Belle2::gearbox::Interface::getAngle
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:301
Belle2::gearbox::Interface::getInt
int getInt(const std::string &path="") const noexcept(false)
Get the parameter path as a int.
Definition: Interface.cc:70
Belle2::qcsmonitor::QcsmonitorFactory
geometry::CreatorFactory< QcsmonitorCreator > QcsmonitorFactory("QCSMONITORCreator")
Creator creates the QCSMONITOR geometry.
Belle2::geometry::Materials::get
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:65
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::GearDir
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:41
Belle2::qcsmonitor::QcsmonitorCreator::create
virtual void create(const GearDir &content, G4LogicalVolume &topVolume, geometry::GeometryTypes type)
Function to actually create the geometry, has to be overridden by derived classes.
Definition: QcsmonitorCreator.cc:60
Belle2::qcsmonitor::QcsmonitorCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector QCSMONITOR.
Definition: QcsmonitorCreator.h:36
Belle2::gearbox::Interface::getLength
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:261
Belle2::gearbox::Interface::getArray
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:133
Belle2::PXD::SensitiveDetector
VXD::SensitiveDetector< PXDSimHit, PXDTrueHit > SensitiveDetector
The PXD Sensitive Detector class.
Definition: SensitiveDetector.h:36
Belle2::geometry::GeometryTypes
GeometryTypes
Flag indiciating the type of geometry to be used.
Definition: GeometryManager.h:39