Belle II Software  release-05-01-25
BeamabortCreator.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/beamabort/geometry/BeamabortCreator.h>
12 #include <beast/beamabort/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 "G4SubtractionSolid.hh"
30 #include <G4UserLimits.hh>
31 
32 //Visualization
33 #include "G4Colour.hh"
34 #include <G4VisAttributes.hh>
35 
36 using namespace std;
37 using namespace boost;
38 
39 namespace Belle2 {
46  namespace beamabort {
47 
48  // Register the creator
50  geometry::CreatorFactory <BeamabortCreator> BeamabortFactory("BEAMABORTCreator");
51 
52  BeamabortCreator::BeamabortCreator() : m_sensitive(0)
53  {
55  }
56 
57  BeamabortCreator::~BeamabortCreator()
58  {
59  if (m_sensitive) delete m_sensitive;
60  }
61 
62  void BeamabortCreator::create(const GearDir& content, G4LogicalVolume& topVolume,
63  geometry::GeometryTypes /* type */)
64  {
65 
66  B2INFO("BeamabortCreator phase2");
67  //Visualization Attributes
68  G4VisAttributes* orange = new G4VisAttributes(G4Colour(1, 2, 0));
69  orange->SetForceAuxEdgeVisible(true);
70  G4VisAttributes* magenta = new G4VisAttributes(G4Colour(1, 0, 1));
71  magenta->SetForceAuxEdgeVisible(true);
72  //lets get the stepsize parameter with a default value of 5 µm
73  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
74 
75  //Lets loop over all the Active nodes
76  BOOST_FOREACH(
77  const GearDir & activeParams, content.getNodes("Active")) {
78 
79  int phase = activeParams.getInt("phase");
80  //Positioned PIN diodes
81  std::vector<double> x_pos;
82  std::vector<double> y_pos;
83  std::vector<double> z_pos;
84  std::vector<double> phi;
85  std::vector<double> thetaX;
86  std::vector<double> thetaY;
87  std::vector<double> thetaZ;
88  std::vector<double> svdAngle;
89  std::vector<double> r;
90  std::vector<double> deltaX;
91  unsigned int dimz = 0;
92 
93  for (double z : activeParams.getArray("z", {0})) {
94  z *= CLHEP::cm;
95  z_pos.push_back(z);
96  dimz++;
97  }
98 
99  for (double ThetaZ : activeParams.getArray("ThetaZ", {0})) {
100  thetaZ.push_back(ThetaZ);
101  }
102 
103  if (thetaZ.size() != dimz) { B2ERROR("Diamond data not consistent (i.e. not same number of all position parmeters)"); return;}
104 
105  if (phase == 2 || phase == 3) {
106 
107  for (double Phi : activeParams.getArray("Phi", {0})) {
108  phi.push_back(Phi);
109  }
110 
111  for (double r_dia : activeParams.getArray("r_dia", {0})) {
112  r_dia *= CLHEP::cm;
113  r.push_back(r_dia);
114  }
115  for (double dX : activeParams.getArray("deltaX", {0})) {
116  dX *= CLHEP::cm;
117  deltaX.push_back(dX);
118  }
119  for (double addAngle : activeParams.getArray("addAngle", {0})) {
120  svdAngle.push_back(addAngle);
121  }
122  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 parmeters)"); return;}
123  }
124  if (phase == 1) {
125  for (double x : activeParams.getArray("x", {0})) {
126  x *= CLHEP::cm;
127  x_pos.push_back(x);
128  }
129  for (double y : activeParams.getArray("y", {0})) {
130  y *= CLHEP::cm;
131  y_pos.push_back(y);
132  }
133 
134  for (double ThetaX : activeParams.getArray("ThetaX", {0})) {
135  thetaX.push_back(ThetaX);
136  }
137  for (double ThetaY : activeParams.getArray("ThetaY", {0})) {
138  thetaY.push_back(ThetaY);
139  }
140  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 parmeters)"); return;}
141  }
142 
143 
144  //create beamabort package
145  G4double dx_opa = 12. / 2. * CLHEP::mm;
146  G4double dy_opa = 18. / 2. * CLHEP::mm;
147  G4double dz_opa = 3.1 / 2. * CLHEP::mm;
148  G4VSolid* s_pa = new G4Box("s_opa", dx_opa, dy_opa, dz_opa);
149  G4double dx_ipa = 6. / 2. * CLHEP::mm;
150  G4double dy_ipa = 6. / 2. * CLHEP::mm;
151  G4double dz_ipa = 5.6 / 2. * CLHEP::mm;
152  G4VSolid* s_ipa = new G4Box("s_ipa", dx_ipa, dy_ipa, dz_ipa);
153  s_pa = new G4SubtractionSolid("s_pa", s_pa, s_ipa, 0, G4ThreeVector(0., 4.0, 0.));
154  G4LogicalVolume* l_pa = new G4LogicalVolume(s_pa, G4Material::GetMaterial("Al6061"), "l_pa");
155  l_pa->SetVisAttributes(magenta);
156  G4Transform3D transform;
157  for (unsigned int i = 0; i < dimz; i++) {
158  B2INFO("DIA-" << i << "RotateZ3D phi: " << phi[i]);
159 
160  if (phase == 1) {
161  transform = G4Translate3D(x_pos[i], y_pos[i], z_pos[i]) * G4RotateX3D(thetaX[i]) *
162  G4RotateY3D(thetaY[i]) * G4RotateZ3D(thetaZ[i]);
163  } else if (phase == 2 || phase == 3) {
164  transform = G4Translate3D(deltaX[i], 0, 0) * G4RotateZ3D(phi[i]) *
165  G4Translate3D(r[i], 0, z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ[i]) *
166  G4RotateZ3D(svdAngle[i]) * G4RotateY3D(M_PI / 2);
167  }
168  new G4PVPlacement(transform, l_pa, TString::Format("p_dia_pa_%d", i).Data(), &topVolume, false, 0);
169  B2INFO("DIA-" << i << " placed at: " << transform.getTranslation() << " mm ");
170  }
171 
172  //create beamabort volumes
173  G4double dx_ba = 4.5 / 2. * CLHEP::mm;
174  G4double dy_ba = 4.5 / 2. * CLHEP::mm;
175  G4double dz_ba = 0.5 / 2. * CLHEP::mm;
176  G4Box* s_BEAMABORT = new G4Box("s_BEAMABORT", dx_ba, dy_ba, dz_ba);
177  G4LogicalVolume* l_BEAMABORT = new G4LogicalVolume(s_BEAMABORT, geometry::Materials::get("Diamond"),
178  "l_BEAMABORT", 0, m_sensitive);
179  l_BEAMABORT->SetVisAttributes(orange);
180 
181  //Lets limit the Geant4 stepsize inside the volume
182  l_BEAMABORT->SetUserLimits(new G4UserLimits(stepSize));
183  l_BEAMABORT->SetVisAttributes(orange);
184 
185  for (unsigned int i = 0; i < dimz; i++) {
186 
187  if (phase == 1)
188  transform = G4Translate3D(x_pos[i], y_pos[i],
189  z_pos[i]) * G4RotateX3D(thetaX[i]) * G4RotateY3D(thetaY[i]) * G4RotateZ3D(thetaZ[i]);
190  if (phase == 2 || phase == 3)
191  transform = G4Translate3D(deltaX[i], 0, 0) * G4RotateZ3D(phi[i]) *
192  G4Translate3D(r[i], 0, z_pos[i]) * G4RotateX3D(-M_PI / 2 - thetaZ[i]) *
193  G4RotateZ3D(svdAngle[i]) * G4RotateY3D(M_PI / 2) * G4Translate3D(0, 4.2, 0);
194 
195  new G4PVPlacement(transform, l_BEAMABORT, TString::Format("p_dia_%d", i).Data(), &topVolume, false, i);
196 
197  B2INFO("DIA-sensitive volume-" << i << " placed at: " << transform.getTranslation() << " mm "
198  << " at phi angle = " << phi[i] << " at theta angle = "
199  << thetaZ[i]);
200  B2INFO("DIA-sensitive volume-" << i << " G4RotateZ3D of: phi= " << phi[i] << " G4RotateX3D = "
201  << (-M_PI / 2 - thetaZ[i]));
202  }
203  }
204  }
205  }
207 }
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::geometry::Materials::get
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:65
Belle2::beamabort::BeamabortFactory
geometry::CreatorFactory< BeamabortCreator > BeamabortFactory("BEAMABORTCreator")
Creator creates the BEAMABORT geometry.
Belle2::beamabort::BeamabortCreator::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: BeamabortCreator.cc:62
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::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::beamabort::BeamabortCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector BEAMABORT.
Definition: BeamabortCreator.h:36
Belle2::geometry::GeometryTypes
GeometryTypes
Flag indiciating the type of geometry to be used.
Definition: GeometryManager.h:39