Belle II Software  release-05-01-25
Fei4Creator.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/srsensor/geometry/Fei4Creator.h>
12 #include <beast/srsensor/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 <boost/format.hpp>
20 #include <boost/foreach.hpp>
21 #include <boost/algorithm/string.hpp>
22 
23 #include <G4LogicalVolume.hh>
24 #include <G4PVPlacement.hh>
25 
26 //Shapes
27 #include <G4Box.hh>
28 #include <G4UserLimits.hh>
29 
30 using namespace std;
31 using namespace boost;
32 
33 namespace Belle2 {
40  namespace srsensor {
41 
42  // Register the creator
44  geometry::CreatorFactory<Fei4Creator> Fei4Factory("FEI4Creator");
45 
46  Fei4Creator::Fei4Creator(): m_sensitive(0)
47  {
49  }
50 
51  Fei4Creator::~Fei4Creator()
52  {
53  if (m_sensitive) delete m_sensitive;
54  }
55 
56  void Fei4Creator::create(const GearDir& content, G4LogicalVolume& topVolume, geometry::GeometryTypes /* type */)
57  {
58  //lets get the stepsize parameter with a default value of 5 µm
59  double stepSize = content.getLength("stepSize", 5 * CLHEP::um);
60 
61  //no get the array. Notice that the default framework unit is cm, so the
62  //values will be automatically converted
63  vector<double> bar_fei4 = content.getArray("bar_fei4");
64  B2INFO("Contents of bar_fei4: ");
65  BOOST_FOREACH(double value, bar_fei4) {
66  B2INFO("value: " << value);
67  }
68  int fei4Nb = 100;
69  //Lets loop over all the Active nodes
70  BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
71 
72  //create fei4 volume
73  G4Box* s_FEI4 = new G4Box("s_FEI4",
74  activeParams.getLength("fei4_dx")*CLHEP::cm,
75  activeParams.getLength("fei4_dy")*CLHEP::cm,
76  activeParams.getLength("fei4_dz")*CLHEP::cm);
77 
78  string matFEI4 = activeParams.getString("MaterialFEI4");
79  G4LogicalVolume* l_FEI4 = new G4LogicalVolume(s_FEI4, geometry::Materials::get(matFEI4), "l_FEI4", 0, m_sensitive);
80 
81  //Lets limit the Geant4 stepsize inside the volume
82  l_FEI4->SetUserLimits(new G4UserLimits(stepSize));
83 
84  //position fei4 volume
85  G4ThreeVector FEI4pos = G4ThreeVector(
86  activeParams.getLength("x_fei4") * CLHEP::cm,
87  activeParams.getLength("y_fei4") * CLHEP::cm,
88  activeParams.getLength("z_fei4") * CLHEP::cm
89  );
90 
91  G4RotationMatrix* rot_fei4 = new G4RotationMatrix();
92  rot_fei4->rotateX(activeParams.getAngle("AngleX"));
93  rot_fei4->rotateY(activeParams.getAngle("AngleY"));
94  rot_fei4->rotateZ(activeParams.getAngle("AngleZ"));
95  //geometry::setColor(*l_FEI4, "#006699");
96 
97  new G4PVPlacement(rot_fei4, FEI4pos, l_FEI4, "p_FEI4", &topVolume, false, fei4Nb);
98 
99  fei4Nb++;
100  }
101 
102  }
103  } // srsensor namespace
105 } // Belle2 namespace
Belle2::srsensor::Fei4Factory
geometry::CreatorFactory< Fei4Creator > Fei4Factory("FEI4Creator")
Creator creates the FEI4 geometry.
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::geometry::Materials::get
static G4Material * get(const std::string &name)
Find given material.
Definition: Materials.h:65
Belle2::srsensor::Fei4Creator::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: Fei4Creator.cc:56
Belle2::srsensor::DiamondCreator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector DIAMOND.
Definition: DiamondCreator.h:36
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::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::GearDir::getString
virtual std::string getString(const std::string &path="") const noexcept(false) override
Get the parameter path as a string.
Definition: GearDir.h:79
Belle2::srsensor::Fei4Creator::m_sensitive
SensitiveDetector * m_sensitive
SensitiveDetector FEI4.
Definition: Fei4Creator.h:36
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