Belle II Software  release-08-01-10
GeoHeavyMetalShieldCreator.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 <vxd/geometry/GeoHeavyMetalShieldCreator.h>
10 
11 #include <geometry/Materials.h>
12 #include <geometry/CreatorFactory.h>
13 #include <geometry/utilities.h>
14 #include <framework/gearbox/GearDir.h>
15 #include <framework/gearbox/Unit.h>
16 #include <framework/logging/Logger.h>
17 
18 #include <cmath>
19 #include <boost/format.hpp>
20 #include <boost/algorithm/string.hpp>
21 #include <boost/math/special_functions/sign.hpp>
22 
23 #include <G4LogicalVolume.hh>
24 #include <G4PVPlacement.hh>
25 #include <G4AssemblyVolume.hh>
26 
27 // Shapes
28 #include <G4Box.hh>
29 #include <G4Trd.hh>
30 #include <G4Polycone.hh>
31 #include <G4SubtractionSolid.hh>
32 
33 using namespace std;
34 using namespace boost;
35 
36 namespace Belle2 {
42  using namespace geometry;
43 
45  namespace VXD {
46 
49 
50  HeavyMetalShieldGeometryPar GeoHeavyMetalShieldCreator::createConfiguration(const GearDir& param)
51  {
52  HeavyMetalShieldGeometryPar heavyMetalShieldGeometryPar;
53  //Read the definition of all shields
54  for (const GearDir& shield : param.getNodes("Shield")) {
55  VXDPolyConePar shieldPar(
56  shield.getString("@name"),
57  shield.getString("Material", "Air"),
58  shield.getAngle("minPhi", 0),
59  shield.getAngle("maxPhi", 2 * M_PI),
60  (shield.getNodes("Cutout").size() > 0),
61  //shield.getLength("Cutout/width", 0.),
62  shield.getLength("Cutout/width1", 0.),
63  shield.getLength("Cutout/width2", 0.),
64  shield.getLength("Cutout/height", 0.),
65  shield.getLength("Cutout/depth", 0.)
66  );
67 
68  for (const GearDir& plane : shield.getNodes("Plane")) {
69  VXDPolyConePlanePar planePar(
70  plane.getLength("posZ"),
71  plane.getLength("innerRadius"),
72  plane.getLength("outerRadius")
73  );
74  shieldPar.getPlanes().push_back(planePar);
75  }
76 
77  heavyMetalShieldGeometryPar.getShields().push_back(shieldPar);
78  }
79  return heavyMetalShieldGeometryPar;
80  }
81 
82  void GeoHeavyMetalShieldCreator::createGeometry(const HeavyMetalShieldGeometryPar& parameters, G4LogicalVolume& topVolume,
84  {
85 
86  // Create the shields
87  for (const VXDPolyConePar& shield : parameters.getShields()) {
88 
89  string name = shield.getName();
90  double minZ(0), maxZ(0);
91 
92  // Create a polycone
93  double minPhi = shield.getMinPhi();
94  double dPhi = shield.getMaxPhi() - minPhi;
95  int nPlanes = shield.getPlanes().size();
96  if (nPlanes < 2) {
97  B2ERROR("Polycone needs at least two planes");
98  return ;
99  }
100  std::vector<double> z(nPlanes, 0);
101  std::vector<double> rMin(nPlanes, 0);
102  std::vector<double> rMax(nPlanes, 0);
103  int index(0);
104  minZ = numeric_limits<double>::infinity();
105  maxZ = -numeric_limits<double>::infinity();
106 
107 
108  for (const VXDPolyConePlanePar& plane : shield.getPlanes()) {
109  z[index] = plane.getPosZ() / Unit::mm;
110  minZ = min(minZ, z[index]);
111  maxZ = max(maxZ, z[index]);
112  rMin[index] = plane.getInnerRadius() / Unit::mm;
113  rMax[index] = plane.getOuterRadius() / Unit::mm;
114  ++index;
115  }
116 
117  G4VSolid* geoShield = new G4Polycone(name + " IR Shield", minPhi, dPhi, nPlanes, z.data(), rMin.data(), rMax.data());
118 
119  // Cutouts (if present)
120  if (shield.getDoCutOut()) {
121  //double sizeX = shield.getCutOutWidth() / Unit::mm / 2.;
122  double sizeX1 = shield.getCutOutWidth1() / Unit::mm / 2.;
123  double sizeX2 = shield.getCutOutWidth2() / Unit::mm / 2.;
124  double sizeY = shield.getCutOutHeight() / Unit::mm / 2.;
125  double depth2 = shield.getCutOutDepth() / Unit::mm / 2.;
126  double sizeZ = (maxZ - minZ) / 2.;
127  double sign = math::sign<double>(minZ);
128  double minAbsZ = min(fabs(minZ), fabs(maxZ));
129 
130  G4ThreeVector origin1(0, 0, sign * (minAbsZ + sizeZ));
131  G4ThreeVector origin2(0, 0, sign * (minAbsZ + depth2));
132 
133  //G4Box* box1 = new G4Box("Cutout", sizeX, sizeY, sizeZ);
134  G4Trd* box1 = new G4Trd("Cutout1", sizeX1, sizeX2, sizeY, sizeY, sizeZ);
135  G4Box* box2 = new G4Box("Cutout2", 100 / Unit::mm, sizeY, depth2);
136 
137  geoShield = new G4SubtractionSolid(name + " IR Shield", geoShield, box1, G4Translate3D(origin1));
138  geoShield = new G4SubtractionSolid(name + " IR Shield", geoShield, box2, G4Translate3D(origin2));
139  }
140 
141  string materialName = shield.getMaterial();
142  G4Material* material = Materials::get(materialName);
143  if (!material) B2FATAL("Material '" << materialName << "', required by " << name << " IR Shield could not be found");
144 
145  G4LogicalVolume* volume = new G4LogicalVolume(geoShield, material, name + " IR Shield");
146  setColor(*volume, "#cc0000");
147  //setVisibility(*volume, false);
148  new G4PVPlacement(0, G4ThreeVector(0, 0, 0), volume, name + " IR Shield", &topVolume, false, 0);
149  }
150  }
151  }
153 }
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
The Class for VXD Heavy Metal Shield.
std::vector< VXDPolyConePar > & getShields(void)
Get shields.
The Class for VXD PolyCone, possibly with coutouts.
std::vector< VXDPolyConePlanePar > & getPlanes(void)
Get planes.
The Class for VXD Polycone Plane.
CreatorFactory< GeoHeavyMetalShieldCreator > GeoHeavyMetalShieldFactory("HeavyMetalShieldCreator")
Register the creator.
void setColor(G4LogicalVolume &volume, const std::string &color)
Set the color of a logical volume.
Definition: utilities.cc:100
GeometryTypes
Flag indiciating the type of geometry to be used.
Abstract base class for different kinds of events.
Very simple class to provide an easy way to register creators with the CreatorManager.