Belle II Software  release-08-01-10
GeoSTRCreator.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 <structure/geometry/GeoSTRCreator.h>
10 #include "structure/dbobjects/STRGeometryPar.h"
11 
12 #include <geometry/Materials.h>
13 #include <geometry/CreatorFactory.h>
14 #include <geometry/utilities.h>
15 
16 #include <framework/logging/Logger.h>
17 
18 #include <math.h> // PI
19 #include <boost/format.hpp>
20 #include <iostream>
21 #include <fstream>
22 
23 #include <G4NistManager.hh>
24 #include <G4LogicalVolume.hh>
25 #include <G4PVPlacement.hh>
26 #include <G4Polycone.hh>
27 #include <G4VisAttributes.hh>
28 
29 
30 using namespace std;
31 
32 namespace Belle2 {
38  using namespace geometry;
39 
40  namespace structure {
41  //-----------------------------------------------------------------
42  // Register the Creator
43  //-----------------------------------------------------------------
46  //-----------------------------------------------------------------
47  // Implementation
48  //-----------------------------------------------------------------
49  GeoSTRCreator::GeoSTRCreator()
50  {
51 
52  }
53 
54 
55  GeoSTRCreator::~GeoSTRCreator()
56  {
57 
58  }
59 
60 
61 
62  void GeoSTRCreator::createGeometry(const STRGeometryPar& parameters, G4LogicalVolume& topVolume, GeometryTypes)
63  {
64 
65 
66 
67 
68  for (int iShield = 0; iShield < parameters.NECLSHIELDS; iShield++) {
69  float _mass = 0;
70  string side;
71  if (iShield == parameters.FWD_ECLSHIELD) {
72  side = "FWD_Shield";
73  } else if (iShield == parameters.BWD_ECLSHIELD) {
74  side = "BWD_Shield";
75  } else {
76  B2FATAL("Only 2 ECL shields should be defined. Can't retrieve info for shield #" << iShield);
77  }
78 
79  int nLayers = parameters.getNLayers(iShield);
80 
81  for (int iLayer = 0; iLayer < nLayers; iLayer++) {
82  G4Material* LayerMat = Materials::get(parameters.getLayerMaterial(iShield, iLayer));
83  int nPlanes = parameters.getLayerNPlanes(iShield, iLayer);
84 
85  //Thread the strings
86  string shapeName = (boost::format("%1%Layer_%2%") % side % (iLayer + 1)).str();
87  string logiVolName = (boost::format("logi%1%Layer_%2%") % side % (iLayer + 1)).str();
88  string physVolName = (boost::format("phys%1%Layer_%2%") % side % (iLayer + 1)).str();
89 
90  G4Polycone* LayerShape = new G4Polycone(shapeName.c_str(), 0, 2 * M_PI, nPlanes,
91  parameters.getLayerPlaneZ(iShield, iLayer),
92  parameters.getLayerPlaneInnerRadius(iShield, iLayer),
93  parameters.getLayerPlaneOuterRadius(iShield, iLayer));
94 
95  // Create logical volume
96  G4LogicalVolume* logiShieldLayer = new G4LogicalVolume(LayerShape, LayerMat, logiVolName, 0, 0, 0);
97 
98  //Place physical volume
99  new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logiShieldLayer, physVolName, &topVolume, false, 0);
100 
101  B2DEBUG(1, "Mass of " << side << " layer " << iLayer
102  << " = " << logiShieldLayer->GetMass() / CLHEP::kg << ".kg.");
103 
104  _mass += (logiShieldLayer->GetMass() / CLHEP::kg);
105  }
106 
107  B2DEBUG(1, "Total mass of side " << side << " = " << _mass << " kg");
108  }
109 
110  //Place pole pieces
111  for (int iPole = 0; iPole < parameters.NPOLEPIECES; iPole++) {
112 
113  string side;
114  if (iPole == parameters.FWD_POLEPIECE) {
115  side = "PolePieceR";
116  } else if (iPole == parameters.BWD_POLEPIECE) {
117  side = "PolePieceL";
118  } else {
119  B2FATAL("Only 2 pole pieces should be defined. Can't retrieve info for pole #" << iPole);
120  }
121 
122 
123  G4Material* PoleMat = Materials::get(parameters.getPoleMaterial(iPole));
124  int nPlanes = parameters.getPoleNPlanes(iPole);
125 
126  //Thread the strings
127  string shapeName = (boost::format("%1%") % side).str();
128  string logiVolName = (boost::format("logi%1%") % side).str();
129  string physVolName = (boost::format("phys%1%") % side).str();
130 
131  G4Polycone* PoleShape = new G4Polycone(shapeName.c_str(), 0, 2 * M_PI, nPlanes,
132  parameters.getPolePlaneZ(iPole),
133  parameters.getPolePlaneInnerRadius(iPole),
134  parameters.getPolePlaneOuterRadius(iPole));
135 
136  // Create logical volume
137  G4LogicalVolume* logiPole = new G4LogicalVolume(PoleShape, PoleMat, logiVolName, 0, 0, 0);
138 
139  //Place physical volume
140  new G4PVPlacement(0, G4ThreeVector(0, 0, 0), logiPole, physVolName, &topVolume, false, 0);
141 
142  B2DEBUG(1, "Total mass of " << side << " = " << logiPole->GetMass() / CLHEP::kg << " kg");
143  }
144 
145 
146  }
147 
148  STRGeometryPar GeoSTRCreator::createConfiguration(const GearDir& param)
149  {
150  STRGeometryPar strGeometryPar;
151 
152  // Get STR geometry parameters from Gearbox (no calculations here)
153  readShield(param, strGeometryPar, "FWD_Shield");
154  readShield(param, strGeometryPar, "BWD_Shield");
155  readPole(param, strGeometryPar, "PolePieceL");
156  readPole(param, strGeometryPar, "PolePieceR");
157 
158 
159  return strGeometryPar;
160  };
161 
162  void GeoSTRCreator::readPole(const GearDir& content, STRGeometryPar& parameters, std::string side)
163  {
164  // Check if method was called using the correct name for the shields
165  std::size_t foundF = side.find("PolePieceR");
166  std::size_t foundB = side.find("PolePieceL");
167 
168  int iPole = 0;
169  if (foundF != std::string::npos) { iPole = parameters.FWD_POLEPIECE; }
170  else if (foundB != std::string::npos) { iPole = parameters.BWD_POLEPIECE; }
171  else { B2FATAL("No data for the Pole Piece requested " << side << "(not found)");}
172 
173 
174  //Thread the strings
175  std::string polePath = (boost::format("/%1%/") % side).str();
176 
177  // Connect the appropriate Gearbox path
178  GearDir poleContent(content);
179  poleContent.append(polePath);
180 
181  // Retrieve material material
182  parameters.setPoleMaterial(iPole, poleContent.getString("Material", "Air"));
183 
184  // Read the shape parameters
185  const std::vector<GearDir> planes = poleContent.getNodes("Plane");
186  parameters.setPoleNPlanes(iPole, planes.size());
187  B2DEBUG(1, "Number of planes on side " << side << " : " << planes.size());
188 
189  for (unsigned int iPlane = 0; iPlane < planes.size(); iPlane++) {
190  parameters.setPolePlaneZ(iPole, iPlane, planes.at(iPlane).getLength("posZ") / Unit::mm);
191  parameters.setPolePlaneInnerRadius(iPole, iPlane, planes.at(iPlane).getLength("innerRadius") / Unit::mm);
192  parameters.setPolePlaneOuterRadius(iPole, iPlane, planes.at(iPlane).getLength("outerRadius") / Unit::mm);
193 
194  }
195  }
196 
197  void GeoSTRCreator::readShield(const GearDir& content, STRGeometryPar& parameters, std::string side)
198  {
199 
200  // Check if method was called using the correct name for the shields
201  std::size_t foundF = side.find("FWD_Shield");
202  std::size_t foundB = side.find("BWD_Shield");
203 
204  int iShield = 0;
205  if (foundF != std::string::npos) { iShield = parameters.FWD_ECLSHIELD; }
206  else if (foundB != std::string::npos) { iShield = parameters.BWD_ECLSHIELD; }
207  else { B2FATAL("No data for the ECL shield called " << side << "(not found)");}
208 
209 
210  std::string gearPath = (boost::format("%1%/Layers/Layer") % side).str();
211 
212  // Retrieve the number of layers in the shield
213  int nLayers = content.getNumberNodes(gearPath);
214  parameters.setNLayers(iShield, nLayers);
215 
216 
217  for (int iLayer = 0 ; iLayer < nLayers ; ++iLayer) {
218  //Thread the strings
219  std::string layerPath = (boost::format("/%1%[%2%]/") % gearPath % (iLayer + 1)).str();
220 
221  // Connect the appropriate Gearbox path
222  GearDir layerContent(content);
223  layerContent.append(layerPath);
224 
225  // Retrieve material material
226  parameters.setLayerMaterial(iShield, iLayer, layerContent.getString("Material", "Air"));
227 
228  // Read the shape parameters
229  const std::vector<GearDir> planes = layerContent.getNodes("Plane");
230  parameters.setLayerNPlanes(iShield, iLayer, planes.size());
231  B2DEBUG(1, "Number of planes on side " << side << " layer " << iLayer
232  << " : " << planes.size());
233 
234  for (unsigned int iPlane = 0; iPlane < planes.size(); iPlane++) {
235  parameters.setLayerPlaneZ(iShield, iLayer, iPlane, planes.at(iPlane).getLength("posZ") / Unit::mm);
236  parameters.setLayerPlaneInnerRadius(iShield, iLayer, iPlane, planes.at(iPlane).getLength("innerRadius") / Unit::mm);
237  parameters.setLayerPlaneOuterRadius(iShield, iLayer, iPlane, planes.at(iPlane).getLength("outerRadius") / Unit::mm);
238  }
239  }
240  }
241  }
243 }
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
void append(const std::string &path)
Append something to the current path, modifying the GearDir in place.
Definition: GearDir.h:52
virtual std::string getString(const std::string &path="") const noexcept(false) override
Get the parameter path as a string.
Definition: GearDir.h:69
The Class for STR geometry.
std::vector< GearDir > getNodes(const std::string &path="") const
Get vector of GearDirs which point to all the nodes the given path evaluates to.
Definition: Interface.cc:21
GeometryTypes
Flag indiciating the type of geometry to be used.
geometry::CreatorFactory< GeoSTRCreator > GeoSTRFactory("STRCreator")
Create factory instance so that the framework can instantiate the STRCreator.
Abstract base class for different kinds of events.
Very simple class to provide an easy way to register creators with the CreatorManager.