Belle II Software  release-08-01-10
GeoSVDCreator.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 <svd/geometry/GeoSVDCreator.h>
10 #include <vxd/geometry/GeoCache.h>
11 #include <svd/geometry/SensorInfo.h>
12 #include <svd/simulation/SensitiveDetector.h>
13 
14 #include <simulation/background/BkgSensitiveDetector.h>
15 
16 #include <geometry/Materials.h>
17 #include <geometry/CreatorFactory.h>
18 #include <geometry/utilities.h>
19 #include <framework/gearbox/GearDir.h>
20 #include <framework/gearbox/Unit.h>
21 #include <framework/logging/Logger.h>
22 
23 #include <cmath>
24 #include <boost/format.hpp>
25 
26 #include <G4LogicalVolume.hh>
27 #include <G4PVPlacement.hh>
28 
29 //Shapes
30 #include <G4Box.hh>
31 #include <G4Tubs.hh>
32 #include <G4Torus.hh>
33 #include <G4Polycone.hh>
34 #include <G4UnionSolid.hh>
35 #include <G4SubtractionSolid.hh>
36 #include <G4Region.hh>
37 
38 using namespace std;
39 using namespace boost;
40 
41 namespace Belle2 {
46  using namespace geometry;
48  namespace SVD {
49 
52 
53  GeoSVDCreator::~GeoSVDCreator()
54  {
55  for (SensorInfo* sensorInfo : m_SensorInfo) delete sensorInfo;
56  m_SensorInfo.clear();
57  }
58 
59  VXD::SensorInfoBase* GeoSVDCreator::createSensorInfo(const VXDGeoSensorPar& sensor)
60  {
61  const SVDSensorInfoPar& infoPar = dynamic_cast<const SVDSensorInfoPar&>(*sensor.getSensorInfo());
62 
63  SensorInfo* info = new SensorInfo(
64  VxdID(0, 0, 0),
65  infoPar.getWidth(),
66  infoPar.getLength(),
67  infoPar.getThickness(),
68  infoPar.getUCells(),
69  infoPar.getVCells(),
70  infoPar.getWidth2()
71  );
72  info->setSensorParams(
73  infoPar.getStripEdgeU(),
74  infoPar.getStripEdgeV(),
75  infoPar.getDepletionVoltage(),
76  infoPar.getBiasVoltage(),
77  infoPar.getBackplaneCapacitanceU(),
78  infoPar.getInterstripCapacitanceU(),
79  infoPar.getCouplingCapacitanceU(),
80  infoPar.getBackplaneCapacitanceV(),
81  infoPar.getInterstripCapacitanceV(),
82  infoPar.getCouplingCapacitanceV(),
83  infoPar.getAduEquivalentU(),
84  infoPar.getAduEquivalentV(),
85  infoPar.getElectronicNoiseU(),
86  infoPar.getElectronicNoiseV(),
87  infoPar.getAduEquivalentSbwU(),
88  infoPar.getAduEquivalentSbwV(),
89  infoPar.getElectronicNoiseSbwU(),
90  infoPar.getElectronicNoiseSbwV()
91  );
92  m_SensorInfo.push_back(info);
93  return info;
94  }
95 
96  SVDSensorInfoPar* GeoSVDCreator::readSensorInfo(const GearDir& sensor)
97  {
98 
99  const double unit_pFcm = 1;
100  // This was 1000 * Unit::fC / Unit::V / Unit::cm; // pF/cm.
101  // We only use ratios of capacities, and this gives nicer numbers.
103  VxdID(0, 0, 0),
104  sensor.getLength("width"),
105  sensor.getLength("length"),
106  sensor.getLength("height"),
107  sensor.getInt("stripsU"),
108  sensor.getInt("stripsV"),
109  sensor.getLength("width2", 0)
110  );
111 
112  info->setSensorParams(
113  sensor.getWithUnit("stripEdgeU"),
114  sensor.getWithUnit("stripEdgeV"),
115  sensor.getWithUnit("DepletionVoltage"),
116  sensor.getWithUnit("BiasVoltage"),
117  sensor.getDouble("BackplaneCapacitanceU") * unit_pFcm,
118  sensor.getDouble("InterstripCapacitanceU") * unit_pFcm,
119  sensor.getDouble("CouplingCapacitanceU") * unit_pFcm,
120  sensor.getDouble("BackplaneCapacitanceV") * unit_pFcm,
121  sensor.getDouble("InterstripCapacitanceV") * unit_pFcm,
122  sensor.getDouble("CouplingCapacitanceV") * unit_pFcm,
123  sensor.getWithUnit("ADUEquivalentU"),
124  sensor.getWithUnit("ADUEquivalentV"),
125  sensor.getWithUnit("ElectronicNoiseU"),
126  sensor.getWithUnit("ElectronicNoiseV"),
127  sensor.getWithUnit("ADUEquivalentSbwU", 0),
128  sensor.getWithUnit("ADUEquivalentSbwV", 0),
129  sensor.getWithUnit("ElectronicNoiseSbwU", 0),
130  sensor.getWithUnit("ElectronicNoiseSbwV", 0)
131  );
132 
133  return info;
134  }
135 
136  VXD::SensitiveDetectorBase* GeoSVDCreator::createSensitiveDetector(
137  VxdID sensorID, const VXDGeoSensor& sensor, const VXDGeoSensorPlacement&)
138  {
139  SensorInfo* sensorInfo = new SensorInfo(dynamic_cast<const SensorInfo&>(*sensor.getSensorInfo()));
140  sensorInfo->setID(sensorID);
141  SensitiveDetector* sensitive = new SensitiveDetector(sensorInfo);
142  return sensitive;
143  }
144 
145  SVDGeometryPar GeoSVDCreator::createConfiguration(const GearDir& content)
146  {
147  SVDGeometryPar svdGeometryPar;
148 
149  //Read prefix ('SVD' or 'PXD')
150  svdGeometryPar.setPrefix(m_prefix);
151 
152  //Read some global parameters
153  VXDGlobalPar globals((float)content.getDouble("ElectronTolerance", 100),
154  (float)content.getDouble("MinimumElectrons", 10),
155  content.getLength("ActiveStepSize", 0.005),
156  content.getBool("ActiveChips", false),
157  content.getBool("SeeNeutrons", false),
158  content.getBool("OnlyPrimaryTrueHits", false),
159  content.getBool("OnlyActiveMaterial", false),
160  (float)content.getLength("DistanceTolerance", 0.005),
161  content.getString("DefaultMaterial", "Air")
162  );
163  svdGeometryPar.setGlobalParams(globals);
164 
165  //Read envelope parameters
166  GearDir envelopeParams(content, "Envelope/");
167  VXDEnvelopePar envelope(envelopeParams.getString("Name", ""),
168  envelopeParams.getString("Material", "Air"),
169  envelopeParams.getString("Color", ""),
170  envelopeParams.getAngle("minPhi", 0),
171  envelopeParams.getAngle("maxPhi", 2 * M_PI),
172  (envelopeParams.getNodes("InnerPoints/point").size() > 0)
173  );
174 
175  for (const GearDir& point : envelopeParams.getNodes("InnerPoints/point")) {
176  pair<double, double> ZXPoint(point.getLength("z"), point.getLength("x"));
177  envelope.getInnerPoints().push_back(ZXPoint);
178  }
179  for (const GearDir& point : envelopeParams.getNodes("OuterPoints/point")) {
180  pair<double, double> ZXPoint(point.getLength("z"), point.getLength("x"));
181  envelope.getOuterPoints().push_back(ZXPoint);
182  }
183  svdGeometryPar.setEnvelope(envelope);
184 
185  // Read alignment for detector m_prefix ('PXD' or 'SVD')
186  string pathAlign = (boost::format("Align[@component='%1%']/") % m_prefix).str();
187  GearDir paramsAlign(GearDir(content, "Alignment/"), pathAlign);
188  if (!paramsAlign) {
189  B2WARNING("Could not find alignment parameters for component " << m_prefix);
190  return svdGeometryPar;
191  }
192  svdGeometryPar.getAlignmentMap()[m_prefix] = VXDAlignmentPar(paramsAlign.getLength("du"),
193  paramsAlign.getLength("dv"),
194  paramsAlign.getLength("dw"),
195  paramsAlign.getAngle("alpha"),
196  paramsAlign.getAngle("beta"),
197  paramsAlign.getAngle("gamma")
198  );
199 
200  //Read the definition of all sensor types
201  GearDir components(content, "Components/");
202  for (const GearDir& paramsSensor : components.getNodes("Sensor")) {
203  string sensorTypeID = paramsSensor.getString("@type");
204 
205  VXDGeoSensorPar sensor(paramsSensor.getString("Material"),
206  paramsSensor.getString("Color", ""),
207  paramsSensor.getLength("width"),
208  paramsSensor.getLength("width2", 0),
209  paramsSensor.getLength("length"),
210  paramsSensor.getLength("height"),
211  paramsSensor.getAngle("angle", 0),
212  paramsSensor.getBool("@slanted", false)
213  );
214  sensor.setActive(VXDGeoComponentPar(
215  paramsSensor.getString("Material"),
216  paramsSensor.getString("Active/Color", "#f00"),
217  paramsSensor.getLength("Active/width"),
218  paramsSensor.getLength("Active/width2", 0),
219  paramsSensor.getLength("Active/length"),
220  paramsSensor.getLength("Active/height")
222  "Active",
223  paramsSensor.getLength("Active/u"),
224  paramsSensor.getLength("Active/v"),
225  paramsSensor.getString("Active/w", "center"),
226  paramsSensor.getLength("Active/woffset", 0)
227  ));
228 
229  SVDSensorInfoPar* svdInfo = readSensorInfo(GearDir(paramsSensor, "Active"));
230  sensor.setSensorInfo(svdInfo);
231  sensor.setComponents(getSubComponents(paramsSensor));
232  svdGeometryPar.getSensorMap()[sensorTypeID] = sensor;
233  svdGeometryPar.getSensorInfos().push_back(svdInfo);
234  }
235 
236  //Build all ladders including Sensors
237  GearDir support(content, "Support/");
238  readHalfShellSupport(support, svdGeometryPar);
239 
240  for (const GearDir& shell : content.getNodes("HalfShell")) {
241 
242  string shellName = m_prefix + "." + shell.getString("@name");
243  string pathShell = (boost::format("Align[@component='%1%']/") % shellName).str();
244  GearDir paramsShell(GearDir(content, "Alignment/"), pathShell);
245  if (!paramsShell) {
246  B2WARNING("Could not find alignment parameters for component " << shellName);
247  return svdGeometryPar;
248  }
249  svdGeometryPar.getAlignmentMap()[shellName] = VXDAlignmentPar(paramsShell.getLength("du"),
250  paramsShell.getLength("dv"),
251  paramsShell.getLength("dw"),
252  paramsShell.getAngle("alpha"),
253  paramsShell.getAngle("beta"),
254  paramsShell.getAngle("gamma")
255  );
256 
257  VXDHalfShellPar halfShell(shell.getString("@name"), shell.getAngle("shellAngle", 0));
258 
259  for (const GearDir& layer : shell.getNodes("Layer")) {
260  int layerID = layer.getInt("@id");
261 
262  readLadder(layerID, components, svdGeometryPar);
263  readLayerSupport(layerID, support, svdGeometryPar);
264  readLadderSupport(layerID, support, svdGeometryPar);
265 
266  //Loop over defined ladders
267  for (const GearDir& ladder : layer.getNodes("Ladder")) {
268  int ladderID = ladder.getInt("@id");
269  double phi = ladder.getAngle("phi", 0);
270  readLadderComponents(layerID, ladderID, content, svdGeometryPar);
271  halfShell.addLadder(layerID, ladderID, phi);
272  }
273  }
274  svdGeometryPar.getHalfShells().push_back(halfShell);
275  }
276 
277  //Create diamond radiation sensors if defined and in background mode
278  GearDir radiationDir(content, "RadiationSensors");
279  if (svdGeometryPar.getGlobalParams().getActiveChips() && radiationDir) {
280  VXDGeoRadiationSensorsPar radiationSensors(
281  m_prefix,
282  radiationDir.getBool("insideEnvelope"),
283  radiationDir.getLength("width"),
284  radiationDir.getLength("length"),
285  radiationDir.getLength("height"),
286  radiationDir.getString("material")
287  );
288 
289  //Add radiation sensor positions
290  for (GearDir& position : radiationDir.getNodes("position")) {
291  VXDGeoRadiationSensorsPositionPar diamonds(position.getLength("z"),
292  position.getLength("radius"),
293  position.getAngle("theta")
294  );
295 
296  //Loop over all phi positions
297  for (GearDir& sensor : position.getNodes("phi")) {
298  //Add sensor with angle and id
299  diamonds.addSensor(sensor.getInt("@id"), sensor.getAngle());
300  }
301  radiationSensors.addPosition(diamonds);
302  }
303  svdGeometryPar.setRadiationSensors(radiationSensors);
304  }
305  return svdGeometryPar;
306  }
307 
308  void GeoSVDCreator::createGeometry(const SVDGeometryPar& parameters, G4LogicalVolume& topVolume, geometry::GeometryTypes)
309  {
310 
311  m_activeStepSize = parameters.getGlobalParams().getActiveStepSize() / Unit::mm;
312  m_activeChips = parameters.getGlobalParams().getActiveChips();
313  m_seeNeutrons = parameters.getGlobalParams().getSeeNeutrons();
314  m_onlyPrimaryTrueHits = parameters.getGlobalParams().getOnlyPrimaryTrueHits();
315  m_distanceTolerance = parameters.getGlobalParams().getDistanceTolerance();
316  m_electronTolerance = parameters.getGlobalParams().getElectronTolerance();
317  m_minimumElectrons = parameters.getGlobalParams().getMinimumElectrons();
318  m_onlyActiveMaterial = parameters.getGlobalParams().getOnlyActiveMaterial();
319  m_defaultMaterial = parameters.getGlobalParams().getDefaultMaterial();
320 
321  G4Material* material = Materials::get(m_defaultMaterial);
322  if (!material) B2FATAL("Default Material of VXD, '" << m_defaultMaterial << "', could not be found");
323 
324 
325  //Build envelope
326  G4LogicalVolume* envelope(0);
327  G4VPhysicalVolume* physEnvelope{nullptr};
328  if (!parameters.getEnvelope().getExists()) {
329  B2INFO("Could not find definition for " + m_prefix + " Envelope, placing directly in top volume");
330  envelope = &topVolume;
331  } else {
332  double minZ(0), maxZ(0);
333  G4Polycone* envelopeCone = geometry::createRotationSolid("Envelope",
334  parameters.getEnvelope().getInnerPoints(),
335  parameters.getEnvelope().getOuterPoints(),
336  parameters.getEnvelope().getMinPhi(),
337  parameters.getEnvelope().getMaxPhi(),
338  minZ, maxZ
339  );
340  envelope = new G4LogicalVolume(envelopeCone, material, m_prefix + ".Envelope");
341  setVisibility(*envelope, false);
342  physEnvelope = new G4PVPlacement(getAlignment(parameters.getAlignment(m_prefix)), envelope, m_prefix + ".Envelope",
343  &topVolume, false, 1);
344 
345  // Set up region for production cuts
346  G4Region* aRegion = new G4Region("SVDEnvelope");
347  envelope->SetRegion(aRegion);
348  aRegion->AddRootLogicalVolume(envelope);
349  }
350 
351  //Read the definition of all sensor types
352  for (const pair<const string, VXDGeoSensorPar>& typeAndSensor : parameters.getSensorMap()) {
353  const string& sensorTypeID = typeAndSensor.first;
354  const VXDGeoSensorPar& paramsSensor = typeAndSensor.second;
355  VXDGeoSensor sensor(
356  paramsSensor.getMaterial(),
357  paramsSensor.getColor(),
358  paramsSensor.getWidth() / Unit::mm,
359  paramsSensor.getWidth2() / Unit::mm,
360  paramsSensor.getLength() / Unit::mm,
361  paramsSensor.getHeight() / Unit::mm,
362  paramsSensor.getSlanted()
363  );
364  sensor.setActive(VXDGeoComponent(
365  paramsSensor.getMaterial(),
366  paramsSensor.getActiveArea().getColor(),
367  paramsSensor.getActiveArea().getWidth() / Unit::mm,
368  paramsSensor.getActiveArea().getWidth2() / Unit::mm,
369  paramsSensor.getActiveArea().getLength() / Unit::mm,
370  paramsSensor.getActiveArea().getHeight() / Unit::mm
371  ), VXDGeoPlacement(
372  "Active",
373  paramsSensor.getActivePlacement().getU() / Unit::mm,
374  paramsSensor.getActivePlacement().getV() / Unit::mm,
375  paramsSensor.getActivePlacement().getW(),
376  paramsSensor.getActivePlacement().getWOffset() / Unit::mm
377  ));
378  sensor.setSensorInfo(createSensorInfo(paramsSensor));
379 
380  vector<VXDGeoPlacement> subcomponents;
381  for (const VXDGeoPlacementPar& component : paramsSensor.getComponents()) {
382  subcomponents.push_back(VXDGeoPlacement(
383  component.getName(),
384  component.getU() / Unit::mm,
385  component.getV() / Unit::mm,
386  component.getW(),
387  component.getWOffset() / Unit::mm
388  ));
389  }
390  sensor.setComponents(subcomponents);
391  m_sensorMap[sensorTypeID] = sensor;
392  }
393 
394  //Read the component cache from DB
395  for (const string& name : parameters.getComponentInsertOder()) {
396  if (m_componentCache.find(name) != m_componentCache.end()) {
397  // already created due to being a sub component of a previous
398  // component. Seems fishy since the information of this component
399  // is in the db at least twice so we could run into
400  // inconsistencies.
401  B2WARNING("Component " << name << " already created from previous subcomponents, should not be here");
402  continue;
403  }
404  const VXDGeoComponentPar& paramsComponent = parameters.getComponent(name);
405  VXDGeoComponent c(
406  paramsComponent.getMaterial(),
407  paramsComponent.getColor(),
408  paramsComponent.getWidth() / Unit::mm,
409  paramsComponent.getWidth2() / Unit::mm,
410  paramsComponent.getLength() / Unit::mm,
411  paramsComponent.getHeight() / Unit::mm
412  );
413  double angle = paramsComponent.getAngle();
414 
415 
416  if (c.getWidth() <= 0 || c.getLength() <= 0 || c.getHeight() <= 0) {
417  B2DEBUG(100, "One dimension empty, using auto resize for component");
418  } else {
419  G4VSolid* solid = createTrapezoidal(m_prefix + "." + name, c.getWidth(), c.getWidth2(), c.getLength(), c.getHeight(), angle);
420  c.setVolume(new G4LogicalVolume(solid, Materials::get(c.getMaterial()), m_prefix + "." + name));
421  }
422 
423  vector<VXDGeoPlacement> subComponents;
424  for (const VXDGeoPlacementPar& paramsSubComponent : paramsComponent.getSubComponents()) {
425  subComponents.push_back(VXDGeoPlacement(
426  paramsSubComponent.getName(),
427  paramsSubComponent.getU() / Unit::mm,
428  paramsSubComponent.getV() / Unit::mm,
429  paramsSubComponent.getW(),
430  paramsSubComponent.getWOffset() / Unit::mm
431  ));
432 
433  }
434 
435  createSubComponents(m_prefix + "." + name, c, subComponents);
436  if (m_activeChips && parameters.getSensitiveChipID(name) >= 0) {
437  int chipID = parameters.getSensitiveChipID(name);
438  B2DEBUG(50, "Creating BkgSensitiveDetector for component " << name << " with chipID " << chipID);
439  BkgSensitiveDetector* sensitive = new BkgSensitiveDetector(m_prefix.c_str(), chipID);
440  c.getVolume()->SetSensitiveDetector(sensitive);
441  m_sensitive.push_back(sensitive);
442  }
443 
444  m_componentCache[name] = c;
445  }
446 
447  //Build all ladders including Sensors
448  VXD::GeoVXDAssembly shellSupport = createHalfShellSupport(parameters);
449 
450  //const std::vector<VXDHalfShellPar>& HalfShells = parameters.getHalfShells();
451  for (const VXDHalfShellPar& shell : parameters.getHalfShells()) {
452  string shellName = shell.getName();
453  m_currentHalfShell = m_prefix + "." + shellName;
454  G4Transform3D shellAlignment = getAlignment(parameters.getAlignment(m_currentHalfShell));
455 
456  // Remember shell coordinate system (into which ladders are inserted)
457  VXD::GeoCache::getInstance().addHalfShellPlacement(m_halfShellVxdIDs[m_currentHalfShell], shellAlignment);
458 
459  //Place shell support
460  double shellAngle = shell.getShellAngle();
461  if (!m_onlyActiveMaterial) shellSupport.place(envelope, shellAlignment * G4RotateZ3D(shellAngle));
462 
463  //const std::map< int, std::vector<std::pair<int, double>> >& Layers = shell.getLayers();
464  for (const std::pair<const int, std::vector<std::pair<int, double>> >& layer : shell.getLayers()) {
465  int layerID = layer.first;
466  const std::vector<std::pair<int, double>>& Ladders = layer.second;
467 
468 
469  setCurrentLayer(layerID, parameters);
470 
471  //Place Layer support
472  VXD::GeoVXDAssembly layerSupport = createLayerSupport(layerID, parameters);
473  if (!m_onlyActiveMaterial) layerSupport.place(envelope, shellAlignment * G4RotateZ3D(shellAngle));
474  VXD::GeoVXDAssembly ladderSupport = createLadderSupport(layerID, parameters);
475 
476  //Loop over defined ladders
477  for (const std::pair<int, double>& ladder : Ladders) {
478  int ladderID = ladder.first;
479  double phi = ladder.second;
480 
481  G4Transform3D ladderPlacement = placeLadder(ladderID, phi, envelope, shellAlignment, parameters);
482  if (!m_onlyActiveMaterial) ladderSupport.place(envelope, ladderPlacement);
483  }
484 
485  }
486  }
487 
488  //Now build cache with all transformations
489  if (physEnvelope) {
490  VXD::GeoCache::getInstance().findVolumes(physEnvelope);
491  } else {
492  //create a temporary placement of the top volume.
493  G4PVPlacement topPlacement(nullptr, G4ThreeVector(0, 0, 0), &topVolume,
494  "temp_Top", nullptr, false, 1, false);
495  //and search for all VXD sensitive sensors within
496  VXD::GeoCache::getInstance().findVolumes(&topPlacement);
497  }
498 
499  //Create diamond radiation sensors if defined and in background mode
500  if (m_activeChips) {
501  if (parameters.getRadiationSensors().getSubDetector() == "") {
502  B2DEBUG(10, "Apparently no radiation sensors defined, skipping");
503  } else {
504  createDiamonds(parameters.getRadiationSensors(), topVolume, *envelope);
505  }
506  }
507  }
508 
509  void GeoSVDCreator::readHalfShellSupport(const GearDir& support, SVDGeometryPar& svdGeometryPar)
510  {
511  if (!support) return;
512 
513  for (const GearDir& params : support.getNodes("HalfShell/RotationSolid")) {
514 
515  VXDRotationSolidPar rotationSolidPar(params.getString("Name", ""),
516  params.getString("Material", "Air"),
517  params.getString("Color", ""),
518  params.getAngle("minPhi", 0),
519  params.getAngle("maxPhi", 2 * M_PI),
520  (params.getNodes("InnerPoints/point").size() > 0)
521  );
522 
523  for (const GearDir& point : params.getNodes("InnerPoints/point")) {
524  pair<double, double> ZXPoint(point.getLength("z"), point.getLength("x"));
525  rotationSolidPar.getInnerPoints().push_back(ZXPoint);
526  }
527  for (const GearDir& point : params.getNodes("OuterPoints/point")) {
528  pair<double, double> ZXPoint(point.getLength("z"), point.getLength("x"));
529  rotationSolidPar.getOuterPoints().push_back(ZXPoint);
530  }
531  svdGeometryPar.getRotationSolids().push_back(rotationSolidPar);
532  }
533  return;
534  }
535 
536  void GeoSVDCreator::readLayerSupport(int layer, const GearDir& support, SVDGeometryPar& svdGeometryPar)
537  {
538  if (!support) return;
539 
540  //Check if there are any endrings defined for this layer. If not we don't create any
541  GearDir endrings(support, (boost::format("Endrings/Layer[@id='%1%']") % layer).str());
542  if (endrings) {
543  svdGeometryPar.getEndrings()[layer] = SVDEndringsPar(support.getString("Endrings/Material"),
544  support.getLength("Endrings/length"),
545  support.getLength("Endrings/gapWidth"),
546  support.getLength("Endrings/baseThickness")
547  );
548 
549  //Create the endrings
550  for (const GearDir& endring : endrings.getNodes("Endring")) {
551  SVDEndringsTypePar endringPar(endring.getString("@name"),
552  endring.getLength("z"),
553  endring.getLength("baseRadius"),
554  endring.getLength("innerRadius"),
555  endring.getLength("outerRadius"),
556  endring.getLength("horizontalBar"),
557  endring.getLength("verticalBar")
558  );
559  svdGeometryPar.getEndrings()[layer].getTypes().push_back(endringPar);
560  }
561  }
562 
563  // Now let's add the cooling pipes to the Support
564  GearDir pipes(support, (boost::format("CoolingPipes/Layer[@id='%1%']") % layer).str());
565  if (pipes) {
566  svdGeometryPar.getCoolingPipes()[layer] = SVDCoolingPipesPar(support.getString("CoolingPipes/Material"),
567  support.getLength("CoolingPipes/outerDiameter"),
568  support.getLength("CoolingPipes/wallThickness"),
569  pipes.getInt("nPipes"),
570  pipes.getAngle("startPhi"),
571  pipes.getAngle("deltaPhi"),
572  pipes.getLength("radius"),
573  pipes.getLength("zstart"),
574  pipes.getLength("zend")
575  );
576  if (pipes.exists("deltaL")) svdGeometryPar.getCoolingPipes()[layer].setDeltaL(pipes.getLength("deltaL"));
577  }
578  return;
579  }
580 
581  void GeoSVDCreator::readLadderSupport(int layer, const GearDir& support, SVDGeometryPar& svdGeometryPar)
582  {
583  if (!support) return;
584 
585  // Check if there are any support ribs defined for this layer. If not return empty assembly
586  GearDir params(support, (boost::format("SupportRibs/Layer[@id='%1%']") % layer).str());
587  if (params) {
588  svdGeometryPar.getSupportRibs()[layer] = SVDSupportRibsPar(support.getLength("SupportRibs/spacing"),
589  support.getLength("SupportRibs/height"),
590  support.getLength("SupportRibs/inner/width"),
591  support.getLength("SupportRibs/outer/width"),
592  support.getLength("SupportRibs/inner/tabLength"),
593  support.getString("SupportRibs/outer/Material"),
594  support.getString("SupportRibs/inner/Material"),
595  support.getString("SupportRibs/outer/Color"),
596  support.getString("SupportRibs/inner/Color"),
597  support.getString("SupportRibs/endmount/Material")
598  );
599 
600  // Get values for the layer if available
601  if (params.exists("spacing")) svdGeometryPar.getSupportRibs()[layer].setSpacing(params.getLength("spacing"));
602  if (params.exists("height")) svdGeometryPar.getSupportRibs()[layer].setHeight(params.getLength("height"));
603 
604  for (const GearDir& box : params.getNodes("box")) {
605  SVDSupportBoxPar boxPar(box.getAngle("theta"),
606  box.getLength("z"),
607  box.getLength("r"),
608  box.getLength("length")
609  );
610  svdGeometryPar.getSupportRibs()[layer].getBoxes().push_back(boxPar);
611  }
612 
613  for (const GearDir& tab : params.getNodes("tab")) {
614  SVDSupportTabPar tabPar(tab.getAngle("theta"),
615  tab.getLength("z"),
616  tab.getLength("r")
617  );
618  svdGeometryPar.getSupportRibs()[layer].getTabs().push_back(tabPar);
619  }
620 
621  for (const GearDir& endmount : params.getNodes("Endmount")) {
622  SVDEndmountPar mountPar(endmount.getString("@name"),
623  endmount.getLength("height"),
624  endmount.getLength("width"),
625  endmount.getLength("length"),
626  endmount.getLength("z"),
627  endmount.getLength("r")
628  );
629  svdGeometryPar.getSupportRibs()[layer].getEndmounts().push_back(mountPar);
630  }
631  }
632  return;
633  }
634 
635  VXD::GeoVXDAssembly GeoSVDCreator::createHalfShellSupport(const SVDGeometryPar& parameters)
636  {
637  VXD::GeoVXDAssembly supportAssembly;
638 
639  //Half shell support is easy as we just add all the defined RotationSolids
640  double minZ(0), maxZ(0);
641 
642  const std::vector<VXDRotationSolidPar>& RotationSolids = parameters.getRotationSolids();
643 
644  for (const VXDRotationSolidPar& component : RotationSolids) {
645 
646  string name = component.getName();
647  string material = component.getMaterial();
648 
649  G4Polycone* solid = geometry::createRotationSolid(name,
650  component.getInnerPoints(),
651  component.getOuterPoints(),
652  component.getMinPhi(),
653  component.getMaxPhi(),
654  minZ, maxZ
655  );
656 
657  G4LogicalVolume* volume = new G4LogicalVolume(
658  solid, geometry::Materials::get(material), m_prefix + ". " + name);
659  geometry::setColor(*volume, component.getColor());
660  supportAssembly.add(volume);
661  }
662  return supportAssembly;
663  }
664 
665 
666 
667  VXD::GeoVXDAssembly GeoSVDCreator::createLayerSupport(int layer, const SVDGeometryPar& parameters)
668  {
669  VXD::GeoVXDAssembly supportAssembly;
670 
671  //Check if there are any endrings defined for this layer. If not we don't create any
672  if (parameters.getEndringsExist(layer)) {
673  const SVDEndringsPar& support = parameters.getEndrings(layer);
674 
675  string material = support.getMaterial();
676  double length = support.getLength() / Unit::mm / 2.0;
677  double gapWidth = support.getGapWidth() / Unit::mm;
678  double baseThickness = support.getBaseThickness() / Unit::mm / 2.0;
679 
680  //Create the endrings
681  const std::vector<SVDEndringsTypePar>& Endrings = support.getTypes();
682  for (const SVDEndringsTypePar& endring : Endrings) {
683  double z = endring.getZ() / Unit::mm;
684  double baseRadius = endring.getBaseRadius() / Unit::mm;
685  double innerRadius = endring.getInnerRadius() / Unit::mm;
686  double outerRadius = endring.getOuterRadius() / Unit::mm;
687  double horiBarWidth = endring.getHorizontalBarWidth() / Unit::mm / 2.0;
688  double vertBarWidth = endring.getVerticalBarWidth() / Unit::mm / 2.0;
689 
690  double angle = asin(gapWidth / innerRadius);
691  G4VSolid* endringSolid = new G4Tubs("OuterEndring", innerRadius, outerRadius, length, -M_PI / 2 + angle, M_PI - 2 * angle);
692  angle = asin(gapWidth / baseRadius);
693  G4VSolid* endringBase = new G4Tubs("InnerEndring", baseRadius, baseRadius + baseThickness, length, -M_PI / 2 + angle,
694  M_PI - 2 * angle);
695  endringSolid = new G4UnionSolid("Endring", endringSolid, endringBase);
696 
697  //Now we need the bars which connect the two rings
698  double height = (innerRadius - baseRadius) / 2.0;
699  double x = vertBarWidth + gapWidth;
700  G4Box* verticalBar = new G4Box("VerticalBar", vertBarWidth, height, length);
701  G4Box* horizontalBar = new G4Box("HorizontalBar", height, horiBarWidth, length);
702  endringSolid = new G4UnionSolid("Endring", endringSolid, verticalBar, G4Translate3D(x, baseRadius + height, 0));
703  endringSolid = new G4UnionSolid("Endring", endringSolid, verticalBar, G4Translate3D(x, -(baseRadius + height), 0));
704  endringSolid = new G4UnionSolid("Endring", endringSolid, horizontalBar, G4Translate3D((baseRadius + height), 0, 0));
705 
706  //Finally create the volume and add it to the assembly at the correct z position
707  G4LogicalVolume* endringVolume = new G4LogicalVolume(
708  endringSolid, geometry::Materials::get(material),
709  (boost::format("%1%.Layer%2%.%3%") % m_prefix % layer % endring.getName()).str());
710  supportAssembly.add(endringVolume, G4TranslateZ3D(z));
711  }
712  }
713 
714  //Check if there are any coling pipes defined for this layer. If not we don't create any
715  if (parameters.getCoolingPipesExist(layer)) {
716  const SVDCoolingPipesPar& pipes = parameters.getCoolingPipes(layer);
717 
718  string material = pipes.getMaterial();
719  double outerRadius = pipes.getOuterDiameter() / Unit::mm / 2.0;
720  double innerRadius = outerRadius - pipes.getWallThickness() / Unit::mm;
721  int nPipes = pipes.getNPipes();
722  double startPhi = pipes.getStartPhi();
723  double deltaPhi = pipes.getDeltaPhi();
724  double radius = pipes.getRadius() / Unit::mm;
725  double zstart = pipes.getZStart() / Unit::mm;
726  double zend = pipes.getZEnd() / Unit::mm;
727  double zlength = (zend - zstart) / 2.0;
728 
729 
730  // There are two parts: the straight pipes and the bendings. So we only need two different volumes
731  // which we place multiple times
732  G4Tubs* pipeSolid = new G4Tubs("CoolingPipe", innerRadius, outerRadius, zlength, 0, 2 * M_PI);
733  G4LogicalVolume* pipeVolume = new G4LogicalVolume(
734  pipeSolid, geometry::Materials::get(material),
735  (boost::format("%1%.Layer%2%.CoolingPipe") % m_prefix % layer).str());
736  geometry::setColor(*pipeVolume, "#ccc");
737 
738  G4Torus* bendSolid = new G4Torus("CoolingBend", innerRadius, outerRadius, sin(deltaPhi / 2.0)*radius, -M_PI / 2, M_PI);
739  G4LogicalVolume* bendVolume = new G4LogicalVolume(
740  bendSolid, geometry::Materials::get(material),
741  (boost::format("%1%.Layer%2%.CoolingBend") % m_prefix % layer).str());
742 
743  // Last pipe may be closer, thus we need additional bending
744  if (pipes.getDeltaL() > 0) {
745  double deltaL = pipes.getDeltaL() / Unit::mm;
746  G4Torus* bendSolidLast = new G4Torus("CoolingBendLast", innerRadius, outerRadius, sin(deltaPhi / 2.0) * radius - deltaL / 2.0,
747  -M_PI / 2, M_PI);
748  G4LogicalVolume* bendVolumeLast = new G4LogicalVolume(bendSolidLast, geometry::Materials::get(material),
749  (boost::format("%1%.Layer%2%.CoolingBendLast") % m_prefix % layer).str());
750  --nPipes;
751 
752  // Place the last straight pipe
753  G4Transform3D placement_pipe = G4RotateZ3D(startPhi + (nPipes - 0.5) * deltaPhi) * G4Translate3D(cos(deltaPhi / 2.0) * radius,
754  sin(deltaPhi / 2.0) * radius - deltaL, zstart + zlength);
755  supportAssembly.add(pipeVolume, placement_pipe);
756 
757  // Place forward or backward bend
758  double zpos = nPipes % 2 > 0 ? zend : zstart;
759  // Calculate transformation
760  G4Transform3D placement = G4RotateZ3D(startPhi + (nPipes - 0.5) * deltaPhi) * G4Translate3D(cos(deltaPhi / 2.0) * radius,
761  -deltaL / 2.0, zpos) * G4RotateY3D(M_PI / 2);
762  // If we are at the forward side we rotate the bend by 180 degree
763  if (nPipes % 2 > 0) {
764  placement = placement * G4RotateZ3D(M_PI);
765  }
766  // And place the bend
767  supportAssembly.add(bendVolumeLast, placement);
768  }
769 
770  for (int i = 0; i < nPipes; ++i) {
771  // Place the straight pipes
772  G4Transform3D placement_pipe = G4RotateZ3D(startPhi + i * deltaPhi) * G4Translate3D(radius, 0, zstart + zlength);
773  supportAssembly.add(pipeVolume, placement_pipe);
774 
775  // This was the easy part, now lets add the connection between the pipes. We only need n-1 bendings
776  if (i > 0) {
777  // Place forward or backward bend
778  double zpos = i % 2 > 0 ? zend : zstart;
779  // Calculate transformation
780  G4Transform3D placement = G4RotateZ3D(startPhi + (i - 0.5) * deltaPhi) * G4Translate3D(cos(deltaPhi / 2.0) * radius, 0,
781  zpos) * G4RotateY3D(M_PI / 2);
782  // If we are at the forward side we rotate the bend by 180 degree
783  if (i % 2 > 0) {
784  placement = placement * G4RotateZ3D(M_PI);
785  }
786  // And place the bend
787  supportAssembly.add(bendVolume, placement);
788  }
789  }
790  }
791 
792  return supportAssembly;
793  }
794 
795 
796 
797  VXD::GeoVXDAssembly GeoSVDCreator::createLadderSupport(int layer, const SVDGeometryPar& parameters)
798  {
799  VXD::GeoVXDAssembly supportAssembly;
800 
801  if (!parameters.getSupportRibsExist(layer)) return supportAssembly;
802  const SVDSupportRibsPar& support = parameters.getSupportRibs(layer);
803 
804  // Get the common values for all layers
805  double spacing = support.getSpacing() / Unit::mm / 2.0;
806  double height = support.getHeight() / Unit::mm / 2.0;
807  double innerWidth = support.getInnerWidth() / Unit::mm / 2.0;
808  double outerWidth = support.getOuterWidth() / Unit::mm / 2.0;
809  double tabLength = support.getTabLength() / Unit::mm / 2.0;
810  G4VSolid* inner(0);
811  G4VSolid* outer(0);
812  G4Transform3D placement;
813 
814 
815  // Now lets create the ribs by adding all boxes to form one union solid
816  const std::vector<SVDSupportBoxPar>& Boxes = support.getBoxes();
817  for (const SVDSupportBoxPar& box : Boxes) {
818  double theta = box.getTheta();
819  double zpos = box.getZ() / Unit::mm;
820  double rpos = box.getR() / Unit::mm;
821  double length = box.getLength() / Unit::mm / 2.0;
822  G4Box* innerBox = new G4Box("innerBox", height, innerWidth, length);
823  G4Box* outerBox = new G4Box("outerBox", height, outerWidth, length);
824  if (!inner) {
825  inner = innerBox;
826  outer = outerBox;
827  placement = G4Translate3D(rpos, 0, zpos) * G4RotateY3D(theta);
828  } else {
829  G4Transform3D relative = placement.inverse() * G4Translate3D(rpos, 0, zpos) * G4RotateY3D(theta);
830  inner = new G4UnionSolid("innerBox", inner, innerBox, relative);
831  outer = new G4UnionSolid("outerBox", outer, outerBox, relative);
832  }
833  }
834 
835  // Now lets add the tabs
836  const std::vector<SVDSupportTabPar>& Tabs = support.getTabs();
837  for (const SVDSupportTabPar& tab : Tabs) {
838  double theta = tab.getTheta();
839  double zpos = tab.getZ() / Unit::mm;
840  double rpos = tab.getR() / Unit::mm;
841  G4Box* innerBox = new G4Box("innerBox", height, innerWidth, tabLength);
842  if (!inner) {
843  inner = innerBox;
844  placement = G4Translate3D(rpos, 0, zpos) * G4RotateY3D(theta);
845  } else {
846  G4Transform3D relative = placement.inverse() * G4Translate3D(rpos, 0, zpos) * G4RotateY3D(theta);
847  inner = new G4UnionSolid("innerBox", inner, innerBox, relative);
848  }
849  }
850 
851  // Now lets create forward and backward endmounts for the ribs
852  const std::vector<SVDEndmountPar>& Endmounts = support.getEndmounts();
853  for (const SVDEndmountPar& endmount : Endmounts) {
854  double endMountHeight = endmount.getHeight() / Unit::mm / 2.0;
855  double endMountWidth = endmount.getWidth() / Unit::mm / 2.0;
856  double endMountLength = endmount.getLength() / Unit::mm / 2.0;
857  double zpos = endmount.getZ() / Unit::mm;
858  double rpos = endmount.getR() / Unit::mm;
859  G4VSolid* endmountBox = new G4Box("endmountBox", endMountHeight, endMountWidth, endMountLength);
860  if (outer) { // holes for the ribs
861  endmountBox = new G4SubtractionSolid("endmountBox", endmountBox, outer, G4TranslateY3D(-spacing)*placement * G4Translate3D(-rpos, 0,
862  -zpos));
863  endmountBox = new G4SubtractionSolid("endmountBox", endmountBox, outer, G4TranslateY3D(spacing)*placement * G4Translate3D(-rpos, 0,
864  -zpos));
865  }
866  G4LogicalVolume* endmountVolume = new G4LogicalVolume(
867  endmountBox, geometry::Materials::get(support.getEndmountMaterial()),
868  (boost::format("%1%.Layer%2%.%3%Endmount") % m_prefix % layer % endmount.getName()).str());
869  supportAssembly.add(endmountVolume, G4Translate3D(rpos, 0, zpos));
870  }
871 
872  // If there has been at least one Box, create the volumes and add them to the assembly
873  if (inner) {
874  outer = new G4SubtractionSolid("outerBox", outer, inner);
875  G4LogicalVolume* outerVolume = new G4LogicalVolume(
876  outer, geometry::Materials::get(support.getOuterMaterial()),
877  (boost::format("%1%.Layer%2%.SupportRib") % m_prefix % layer).str());
878  G4LogicalVolume* innerVolume = new G4LogicalVolume(
879  inner, geometry::Materials::get(support.getInnerMaterial()),
880  (boost::format("%1%.Layer%2%.SupportRib.Airex") % m_prefix % layer).str());
881  geometry::setColor(*outerVolume, support.getOuterColor());
882  geometry::setColor(*innerVolume, support.getInnerColor());
883  supportAssembly.add(innerVolume, G4TranslateY3D(-spacing)*placement);
884  supportAssembly.add(innerVolume, G4TranslateY3D(spacing)*placement);
885  supportAssembly.add(outerVolume, G4TranslateY3D(-spacing)*placement);
886  supportAssembly.add(outerVolume, G4TranslateY3D(spacing)*placement);
887  }
888 
889  // Done, return the finished assembly
890  return supportAssembly;
891  }
892 
893 
894 
895  }
897 }
The Class for BeamBackground Sensitive Detector.
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
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 SVD Cooling Pipes.
double getStartPhi() const
get start phi
double getOuterDiameter() const
get outer diameter
double getDeltaPhi() const
get delta phi
int getNPipes() const
get nPipes
double getDeltaL() const
get deltal
double getZEnd() const
get zend
double getZStart() const
get zstart
double getWallThickness() const
get wall thickness
const std::string & getMaterial() const
get material
double getRadius() const
get radius
The Class for SVD Support Rib Endmounts.
The Class for SVD Endring.
const std::vector< SVDEndringsTypePar > & getTypes() const
get types (forward/backward)
double getGapWidth() const
get gap width
double getBaseThickness() const
get base thickness
const std::string & getMaterial() const
get material
double getLength() const
get length
The Class for SVD Endring Type.
The Class for VXD geometry.
const std::map< int, SVDSupportRibsPar > & getSupportRibs() const
get support ribs
const std::map< int, SVDCoolingPipesPar > & getCoolingPipes() const
get cooling pipes
std::vector< SVDSensorInfoPar * > & getSensorInfos()
get sensorInfos
const std::vector< VXDRotationSolidPar > & getRotationSolids() const
get SVD halfshell Rotation Solids
const std::map< int, SVDEndringsPar > & getEndrings() const
get endrings
The Class for VXD geometry.
double getBiasVoltage() const
Return the bias voltage on the sensor.
double getAduEquivalentSbwU() const
Return ADU equivalent for U strips in Sbw barrel sensor.
double getElectronicNoiseSbwU() const
Return electronic noise in e- for u strips in bw barrel sensors.
double getBackplaneCapacitanceV() const
Return the backplane capacitance/cm for V-side strips.
double getAduEquivalentU() const
Return ADU equivalent for U strips.
double getElectronicNoiseV() const
Return electronic noise in e- for v strips.
double getStripEdgeV() const
Return the distance between end of strip and edge of active area.
double getAduEquivalentSbwV() const
Return ADU equivalent for V strips in Sbw barrel sensor.
double getBackplaneCapacitanceU() const
Return the backplane capacitance/cm for U-side strips.
double getInterstripCapacitanceU() const
Return the interstrip capacitance/cm for U-side strips.
double getElectronicNoiseU() const
Return electronic noise in e- for u strips.
double getStripEdgeU() const
Return the distance between end of strip and edge of active area.
double getCouplingCapacitanceU() const
Return the coupling capacitance/cm for U-side strips.
double getCouplingCapacitanceV() const
Return the coupling capacitance/cm for V-side strips.
double getAduEquivalentV() const
Return ADU equivalent for V strips.
double getDepletionVoltage() const
Return the depletion voltage of the sensor.
double getInterstripCapacitanceV() const
Return the interstrip capacitance/cm for V-side strips.
double getElectronicNoiseSbwV() const
Return electronic noise in e- for v strips in bw barrel sensors.
The Class for SVD Support Box.
The Class for SVD Support Ribs (one layer)
double getHeight() const
get height
double getInnerWidth() const
get inner width
double getTabLength() const
get tabLength
double getOuterWidth() const
get outer width
const std::vector< SVDSupportBoxPar > & getBoxes() const
get boxes
const std::string & getEndmountMaterial() const
get the name of endmount material
const std::vector< SVDEndmountPar > & getEndmounts() const
get endmounts
const std::vector< SVDSupportTabPar > & getTabs() const
get tabs
const std::string & getOuterMaterial() const
get the name of outer material
double getSpacing() const
get spacing
const std::string & getInnerMaterial() const
get the name of inner material
const std::string & getOuterColor() const
get the name of outer color
const std::string & getInnerColor() const
get the name of inner color
The Class for SVD Support Rib Tab.
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
Definition: SensorInfo.h:25
void setID(VxdID id)
Change the SensorID.
Definition: SensorInfo.h:85
The Class for VXD Alignment payload.
The Class for VXD Envelope parameters.
The Class for VXD geometry component.
double getWidth() const
get the width of the component
double getHeight() const
get the height of the component
double getAngle() const
get the angle of the component
const std::vector< VXDGeoPlacementPar > & getSubComponents() const
get sub components
double getWidth2() const
get the forward width of the component, 0 for rectangular
const std::string & getColor() const
get the name of the color for the component
const std::string & getMaterial() const
get the name of the Material for the component
double getLength() const
get the length of the component
Class holding all parameters for an VXD geometry component.
The Class for VXD placement payload.
double getU() const
get local u coordinate where to place the component
const std::string & getW() const
get local w position where to place the component
double getV() const
get local v coordinate where to place the component
double getWOffset() const
get offset to local w position where to place the component
Class holding all parameters to place a VXD geometry subcomponent.
The Class for VXD Radiation Sensor parameters.
void addPosition(const VXDGeoRadiationSensorsPositionPar &position)
add radiation sensor position
The Class for VXD Radiation Sensor Position parameters.
void addSensor(int id, double phi)
add sensor with individual id
The Class for VXD Sensor payload.
const VXDGeoPlacementPar & getActivePlacement() const
get the placement description for the active area
const std::vector< VXDGeoPlacementPar > & getComponents() const
get the list of sub components
bool getSlanted() const
return wether or not the sensor is slanted (usually only the first sensor in layers 4-6)
const VXDGeoComponentPar & getActiveArea() const
get the component description for the active area
Struct holding the information where a sensor should be placed inside the ladder.
Struct holding all parameters for a completeVXD Sensor.
const std::map< std::string, VXDGeoSensorPar > & getSensorMap() const
get sensor map
void setGlobalParams(const VXDGlobalPar &globals)
set global parameters
void setRadiationSensors(const VXDGeoRadiationSensorsPar &diamonds)
set radiation sensor parameters
const VXDGlobalPar & getGlobalParams() const
get global parameters
const std::vector< VXDHalfShellPar > & getHalfShells() const
get half-shell
void setEnvelope(const VXDEnvelopePar &envelope)
set envelope parameters
void setPrefix(const std::string &prefix)
set prefix
std::map< std::string, VXDAlignmentPar > & getAlignmentMap()
get alignmant map
The Class for VXD global paramter payload.
Definition: VXDGlobalPar.h:24
bool getActiveChips() const
Get whether chips are sensitive
Definition: VXDGlobalPar.h:48
The Class for VXD half shell payload.
void addLadder(int layerID, int ladderID, double phi)
add ladder
The Class for VXD Envelope parameters.
const std::list< std::pair< double, double > > & getInnerPoints() const
get inner XZ points
const std::list< std::pair< double, double > > & getOuterPoints() const
get outer XZ points
double getWidth() const
Return the (backward) width of the sensor.
int getVCells() const
Return number of pixel/strips in v direction.
int getUCells() const
Return number of pixel/strips in u direction.
double getThickness() const
Return the thickness of the sensor.
double getWidth2() const
Return forward width for a slanted sensor.
double getLength() const
Return the length of the sensor.
Class to group some Geant4 volumes and place them all at once with a given transformation matrix.
void place(G4LogicalVolume *mother, const G4Transform3D &transform)
Place all the volumes already added to the assembly in the given mother.
void add(G4LogicalVolume *volume, const G4Transform3D &transform=G4Transform3D())
Add a volume to the assembly.
Base class for Sensitive Detector implementation of PXD and SVD.
Sensitive Detector implementation of PXD and SVD.
Base class to provide Sensor Information for PXD and SVD.
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
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:299
bool exists(const std::string &path="") const
Check if a given parameter path exists.
Definition: Interface.h:54
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:259
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
bool getBool(const std::string &path="") const noexcept(false)
Get the parameter path as a bool.
Definition: Interface.cc:80
int getInt(const std::string &path="") const noexcept(false)
Get the parameter path as a int.
Definition: Interface.cc:60
geometry::CreatorFactory< GeoSVDCreator > GeoSVDFactory("SVDCreator")
Register the creator.
void setVisibility(G4LogicalVolume &volume, bool visible)
Helper function to quickly set the visibility of a given volume.
Definition: utilities.cc:108
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.