CMS 3D CMS Logo

volumeHandle.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author N. Amapane - INFN Torino (original developer)
5  */
6 
7 #include "volumeHandle.h"
13 
18 
19 #include "CLHEP/Units/GlobalSystemOfUnits.h"
20 
21 #include <string>
22 #include <iterator>
23 #include <iomanip>
24 #include <iostream>
25 #include <boost/lexical_cast.hpp>
26 
27 using namespace SurfaceOrientation;
28 using namespace std;
29 
31 
32 namespace {
33  // Old DD returns lengths in mm, but CMS code uses cm
34  template <class NumType>
35  inline constexpr NumType convertUnits(NumType millimeters) // Millimeters -> centimeters
36  {
37  return (geant_units::operators::convertMmToCm(millimeters));
38  }
39 } // namespace
40 
41 MagGeoBuilderFromDDD::volumeHandle::volumeHandle(const DDExpandedView &fv, bool expand2Pi, bool debugVal)
42  : magneticfield::BaseVolumeHandle(expand2Pi, debugVal) {
43  name = fv.logicalPart().name().name();
44  copyno = fv.copyno();
45  solid = fv.logicalPart().solid();
46  center_ = GlobalPoint(fv.translation().x() / cm, fv.translation().y() / cm, fv.translation().z() / cm);
47 
48  // ASSUMPTION: volume names ends with "_NUM" where NUM is the volume number
49  string volName = name;
50  volName.erase(0, volName.rfind('_') + 1);
51  volumeno = boost::lexical_cast<unsigned short>(volName);
52 
53  for (int i = 0; i < 6; ++i) {
54  isAssigned[i] = false;
55  }
56 
57  if (debug) {
58  cout.precision(7);
59  }
60 
61  referencePlane(fv);
62 
63  if (solid.shape() == DDSolidShape::ddbox) {
64  DDBox box(solid);
65  double halfX = convertUnits(box.halfX());
66  double halfY = convertUnits(box.halfY());
67  double halfZ = convertUnits(box.halfZ());
68  buildBox(halfX, halfY, halfZ);
69  } else if (solid.shape() == DDSolidShape::ddtrap) {
70  DDTrap trap(solid);
71  double x1 = convertUnits(trap.x1());
72  double x2 = convertUnits(trap.x2());
73  double x3 = convertUnits(trap.x3());
74  double x4 = convertUnits(trap.x4());
75  double y1 = convertUnits(trap.y1());
76  double y2 = convertUnits(trap.y2());
77  double theta = trap.theta();
78  double phi = trap.phi();
79  double halfZ = convertUnits(trap.halfZ());
80  double alpha1 = trap.alpha1();
81  double alpha2 = trap.alpha2();
82  buildTrap(x1, x2, x3, x4, y1, y2, theta, phi, halfZ, alpha1, alpha2);
83  } else if (solid.shape() == DDSolidShape::ddcons) {
84  DDCons cons(solid);
85  double zhalf = convertUnits(cons.zhalf());
86  double rInMinusZ = convertUnits(cons.rInMinusZ());
87  double rOutMinusZ = convertUnits(cons.rOutMinusZ());
88  double rInPlusZ = convertUnits(cons.rInPlusZ());
89  double rOutPlusZ = convertUnits(cons.rOutPlusZ());
90  double startPhi = cons.phiFrom();
91  double deltaPhi = cons.deltaPhi();
92  buildCons(zhalf, rInMinusZ, rOutMinusZ, rInPlusZ, rOutPlusZ, startPhi, deltaPhi);
93  } else if (solid.shape() == DDSolidShape::ddtubs) {
94  DDTubs tubs(solid);
95  double zhalf = convertUnits(tubs.zhalf());
96  double rIn = convertUnits(tubs.rIn());
97  double rOut = convertUnits(tubs.rOut());
98  double startPhi = tubs.startPhi();
99  double deltaPhi = tubs.deltaPhi();
100  buildTubs(zhalf, rIn, rOut, startPhi, deltaPhi);
101  } else if (solid.shape() == DDSolidShape::ddpseudotrap) {
102  DDPseudoTrap ptrap(solid);
103  double x1 = convertUnits(ptrap.x1());
104  double x2 = convertUnits(ptrap.x2());
105  double y1 = convertUnits(ptrap.y1());
106  double y2 = convertUnits(ptrap.y2());
107  double halfZ = convertUnits(ptrap.halfZ());
108  double radius = convertUnits(ptrap.radius());
109  bool atMinusZ = ptrap.atMinusZ();
110  buildPseudoTrap(x1, x2, y1, y2, halfZ, radius, atMinusZ);
111  } else if (solid.shape() == DDSolidShape::ddtrunctubs) {
112  DDTruncTubs tubs(solid);
113  double zhalf = convertUnits(tubs.zHalf()); // half of the z-Axis
114  double rIn = convertUnits(tubs.rIn()); // inner radius
115  double rOut = convertUnits(tubs.rOut()); // outer radius
116  double startPhi = tubs.startPhi(); // angular start of the tube-section
117  double deltaPhi = tubs.deltaPhi(); // angular span of the tube-section
118  double cutAtStart = convertUnits(tubs.cutAtStart()); // truncation at begin of the tube-section
119  double cutAtDelta = convertUnits(tubs.cutAtDelta()); // truncation at end of the tube-section
120  bool cutInside = tubs.cutInside(); // true, if truncation is on the inner side of the tube-section
121  buildTruncTubs(zhalf, rIn, rOut, startPhi, deltaPhi, cutAtStart, cutAtDelta, cutInside);
122  } else {
123  cout << "volumeHandle ctor: Unexpected solid: " << DDSolidShapesName::name(solid.shape()) << endl;
124  }
125 
126  // NOTE: Table name and master sector are no longer taken from xml!
127  // DDsvalues_type sv(fv.mergedSpecifics());
128 
129  // { // Extract the name of associated field file.
130  // std::vector<std::string> temp;
131  // std::string pname = "table";
132  // DDValue val(pname);
133  // DDsvalues_type sv(fv.mergedSpecifics());
134  // if (DDfetch(&sv,val)) {
135  // temp = val.strings();
136  // if (temp.size() != 1) {
137  // cout << "*** WARNING: volume has > 1 SpecPar " << pname << endl;
138  // }
139  // magFile = temp[0];
140 
141  // string find="[copyNo]";
142  // std::size_t j;
143  // for ( ; (j = magFile.find(find)) != string::npos ; ) {
144  // stringstream conv;
145  // conv << setfill('0') << setw(2) << copyno;
146  // string repl;
147  // conv >> repl;
148  // magFile.replace(j, find.length(), repl);
149  // }
150 
151  // } else {
152  // cout << "*** WARNING: volume does not have a SpecPar " << pname << endl;
153  // cout << " DDsvalues_type: " << fv.mergedSpecifics() << endl;
154  // }
155  // }
156 
157  // { // Extract the number of the master sector.
158  // std::vector<double> temp;
159  // const std::string pname = "masterSector";
160  // DDValue val(pname);
161  // if (DDfetch(&sv,val)) {
162  // temp = val.doubles();
163  // if (temp.size() != 1) {
164  // cout << "*** WARNING: volume has > 1 SpecPar " << pname << endl;
165  // }
166  // masterSector = int(temp[0]+.5);
167  // } else {
168  // if (MagGeoBuilderFromDDD::debug) {
169  // cout << "Volume does not have a SpecPar " << pname
170  // << " using: " << copyno << endl;
171  // cout << " DDsvalues_type: " << fv.mergedSpecifics() << endl;
172  // }
173  // masterSector = copyno;
174  // }
175  // }
176 
177  // Get material for this volume
178  if (fv.logicalPart().material().name().name() == "Iron")
179  isIronFlag = true;
180 
181  if (debug) {
182  cout << " RMin = " << theRMin << endl;
183  cout << " RMax = " << theRMax << endl;
184 
185  if (theRMin < 0 || theRN < theRMin || theRMax < theRN)
186  cout << "*** WARNING: wrong RMin/RN/RMax , shape: " << DDSolidShapesName::name(shape()) << endl;
187 
188  cout << "Summary: " << name << " " << copyno << " Shape= " << DDSolidShapesName::name(shape()) << " trasl "
189  << center() << " R " << center().perp() << " phi " << center().phi() << " magFile " << magFile
190  << " Material= " << fv.logicalPart().material().name() << " isIron= " << isIronFlag
191  << " masterSector= " << masterSector << std::endl;
192 
193  cout << " Orientation of surfaces:";
194  std::string sideName[3] = {"positiveSide", "negativeSide", "onSurface"};
195  for (int i = 0; i < 6; ++i) {
196  cout << " " << i << ":" << sideName[surfaces[i]->side(center_, 0.3)];
197  }
198  cout << endl;
199  }
200 }
201 
203  // The refPlane is the "main plane" for the solid. It corresponds to the
204  // x,y plane in the DDD local frame, and defines a frame where the local
205  // coordinates are the same as in DDD.
206  // In the geometry version 85l_030919, this plane is normal to the
207  // beam line for all volumes but pseudotraps, so that global R is along Y,
208  // global phi is along -X and global Z along Z:
209  //
210  // Global(for vol at pi/2) Local
211  // +R (+Y) +Y
212  // +phi(-X) -X
213  // +Z +Z
214  //
215  // For pseudotraps the refPlane is parallel to beam line and global R is
216  // along Z, global phi is along +-X and and global Z along Y:
217  //
218  // Global(for vol at pi/2) Local
219  // +R (+Y) +Z
220  // +phi(-X) +X
221  // +Z +Y
222  //
223  // Note that the frame is centered in the DDD volume center, which is
224  // inside the volume for DDD boxes and (pesudo)trapezoids, on the beam line
225  // for tubs, cons and trunctubs.
226 
227  // In geometry version 1103l, trapezoids have X and Z in the opposite direction
228  // than the above. Boxes are either oriented as described above or in some case
229  // have opposite direction for Y and X.
230 
231  // The global position
232  Surface::PositionType &posResult = center_;
233 
234  // The reference plane rotation
235  DD3Vector x, y, z;
236  fv.rotation().GetComponents(x, y, z);
237  if (debug) {
238  if (x.Cross(y).Dot(z) < 0.5) {
239  cout << "*** WARNING: Rotation is not RH " << endl;
240  }
241  }
242 
243  // The global rotation
244  Surface::RotationType rotResult(float(x.X()),
245  float(x.Y()),
246  float(x.Z()),
247  float(y.X()),
248  float(y.Y()),
249  float(y.Z()),
250  float(z.X()),
251  float(z.Y()),
252  float(z.Z()));
253 
254  refPlane = new GloballyPositioned<float>(posResult, rotResult);
255 
256  // Check correct orientation
257  if (debug) {
258  cout << "Refplane pos " << refPlane->position() << endl;
259 
260  // See comments above for the conventions for orientation.
261  LocalVector globalZdir(0., 0., 1.); // Local direction of the axis along global Z
262  if (solid.shape() == DDSolidShape::ddpseudotrap) {
263  globalZdir = LocalVector(0., 1., 0.);
264  }
265  if (refPlane->toGlobal(globalZdir).z() < 0.) {
266  globalZdir = -globalZdir;
267  }
268 
269  float chk = refPlane->toGlobal(globalZdir).dot(GlobalVector(0, 0, 1));
270  if (chk < .999)
271  cout << "*** WARNING RefPlane check failed!***" << chk << endl;
272  }
273 }
274 
275 std::vector<VolumeSide> MagGeoBuilderFromDDD::volumeHandle::sides() const {
276  std::vector<VolumeSide> result;
277  for (int i = 0; i < 6; ++i) {
278  // If this is just a master volume out of wich a 2pi volume
279  // should be built (e.g. central cylinder), skip the phi boundaries.
280  if (expand && (i == phiplus || i == phiminus))
281  continue;
282 
283  // FIXME: Skip null inner degenerate cylindrical surface
284  if (solid.shape() == DDSolidShape::ddtubs && i == SurfaceOrientation::inner && theRMin < 0.001)
285  continue;
286 
287  ReferenceCountingPointer<Surface> s = const_cast<Surface *>(surfaces[i].get());
288  result.push_back(VolumeSide(s, GlobalFace(i), surfaces[i]->side(center_, 0.3)));
289  }
290  return result;
291 }
292 
294 using namespace magneticfield;
295 
296 #include "buildBox.icc"
297 #include "buildTrap.icc"
298 #include "buildTubs.icc"
299 #include "buildCons.icc"
300 #include "buildPseudoTrap.icc"
301 #include "buildTruncTubs.icc"
Vector3DBase< float, LocalTag >
GloballyPositioned< float >
magneticfield::BaseVolumeHandle::magFile
std::string magFile
Name of magnetic field table file.
Definition: BaseVolumeHandle.h:66
MagGeoBuilderFromDDD::volumeHandle::shape
DDSolidShape shape() const override
Shape of the solid.
Definition: volumeHandle.h:32
SurfaceOrientation::phiplus
Definition: Surface.h:19
DDTubs::zhalf
double zhalf(void) const
Definition: DDSolid.cc:454
TkRotation< float >
DDPseudoTrap::radius
double radius(void) const
radius of the cut-out (neg.) or rounding (pos.)
Definition: DDSolid.cc:200
DDAxes::y
mps_fire.i
i
Definition: mps_fire.py:428
DDTrap::phi
double phi(void) const
Azimuthal angle of the line joining the centres of the faces at -/+pDz.
Definition: DDSolid.cc:139
Cylinder.h
DDTrap::halfZ
double halfZ(void) const
half of the z-Axis
Definition: DDSolid.cc:135
DDTruncTubs::startPhi
double startPhi(void) const
angular start of the tube-section
Definition: DDSolid.cc:172
DDTrap::x4
double x4(void) const
Half-length along x of the side at y=+pDy2 of the face at +pDz.
Definition: DDSolid.cc:153
DDExpandedView::rotation
const DDRotationMatrix & rotation() const
The absolute rotation of the current node.
Definition: DDExpandedView.cc:48
DDCons
Definition: DDSolid.h:292
DDSolidShape::ddtrap
DDPseudoTrap::y2
double y2(void) const
half length along y on +z
Definition: DDSolid.cc:198
magneticfield::BaseVolumeHandle::debug
const bool debug
Definition: BaseVolumeHandle.h:151
DDCons::deltaPhi
double deltaPhi(void) const
Definition: DDSolid.cc:426
DDTrap::x2
double x2(void) const
Half-length along x of the side at y=+pDy1 of the face at -pDz.
Definition: DDSolid.cc:145
DDPseudoTrap::atMinusZ
bool atMinusZ(void) const
true, if cut-out or rounding is on the -z side
Definition: DDSolid.cc:202
volumeHandle
MagGeoBuilderFromDDD::volumeHandle volumeHandle
Definition: volumeHandle.cc:293
MagGeoBuilderFromDDD::volumeHandle::buildTruncTubs
void buildTruncTubs(double zhalf, double rIn, double rOut, double startPhi, double deltaPhi, double cutAtStart, double cutAtDelta, bool cutInside)
MagGeoBuilderFromDDD::volumeHandle::sides
std::vector< VolumeSide > sides() const override
The surfaces and they orientation, as required to build a MagVolume.
Definition: volumeHandle.cc:275
DDBox::halfZ
double halfZ(void) const
Definition: DDSolid.cc:216
testProducerWithPsetDescEmpty_cfi.x2
x2
Definition: testProducerWithPsetDescEmpty_cfi.py:28
DDSolidShape::ddpseudotrap
magneticfield::BaseVolumeHandle::center
const GlobalPoint & center() const
Return the center of the volume.
Definition: BaseVolumeHandle.cc:42
gather_cfg.cout
cout
Definition: gather_cfg.py:144
magneticfield::volumeHandle
Definition: DD4hep_volumeHandle.h:23
DDLogicalPart::material
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
Definition: DDLogicalPart.cc:118
magneticfield::BaseVolumeHandle::isAssigned
bool isAssigned[6]
Definition: BaseVolumeHandle.h:122
DDPseudoTrap::x2
double x2(void) const
half length along x on +z
Definition: DDSolid.cc:194
magneticfield
Definition: MagFieldConfig.h:22
DDExpandedView.h
CoordinateSets.h
ReferenceCountingPointer< Surface >
DDAxes::x
MagGeoBuilderFromDDD::volumeHandle::solid
DDSolid solid
Definition: volumeHandle.h:78
GlobalVector
Global3DVector GlobalVector
Definition: GlobalVector.h:10
magneticfield::BaseVolumeHandle::volumeno
unsigned short volumeno
volume number
Definition: BaseVolumeHandle.h:69
DDCons::phiFrom
double phiFrom(void) const
Definition: DDSolid.cc:424
DDCons::rInMinusZ
double rInMinusZ(void) const
Definition: DDSolid.cc:416
DDTruncTubs::cutAtStart
double cutAtStart(void) const
truncation at begin of the tube-section
Definition: DDSolid.cc:176
MagGeoBuilderFromDDD::volumeHandle::referencePlane
void referencePlane(const DDExpandedView &fv)
Definition: volumeHandle.cc:202
VolumeSide
Definition: VolumeSide.h:15
DDPseudoTrap
Definition: DDSolid.h:117
magneticfield::BaseVolumeHandle::theRMin
double theRMin
Definition: BaseVolumeHandle.h:132
MagGeoBuilderFromDDD::debug
const bool debug
Definition: MagGeoBuilderFromDDD.h:99
DDCons::rInPlusZ
double rInPlusZ(void) const
Definition: DDSolid.cc:420
MagGeoBuilderFromDDD::volumeHandle::volumeHandle
volumeHandle(const DDExpandedView &fv, bool expand2Pi=false, bool debugVal=false)
Definition: volumeHandle.cc:41
DDTruncTubs::cutInside
bool cutInside(void) const
true, if truncation is on the inner side of the tube-section
Definition: DDSolid.cc:180
DDTruncTubs::zHalf
double zHalf(void) const
half of the z-Axis
Definition: DDSolid.cc:166
Plane.h
DDExpandedView::copyno
int copyno() const
Copy number associated with the current node.
Definition: DDExpandedView.cc:54
DDSolid::shape
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:123
testProducerWithPsetDescEmpty_cfi.x1
x1
Definition: testProducerWithPsetDescEmpty_cfi.py:33
alignCSCRings.s
s
Definition: alignCSCRings.py:92
testProducerWithPsetDescEmpty_cfi.y1
y1
Definition: testProducerWithPsetDescEmpty_cfi.py:29
DDTruncTubs
A truncated tube section.
Definition: DDSolid.h:139
magneticfield::BaseVolumeHandle::theRN
double theRN
Definition: BaseVolumeHandle.h:128
DDCons::rOutMinusZ
double rOutMinusZ(void) const
Definition: DDSolid.cc:418
DDSolidShape::ddtubs
SiPixelRawToDigiRegional_cfi.deltaPhi
deltaPhi
Definition: SiPixelRawToDigiRegional_cfi.py:9
DDBox
Interface to a Box.
Definition: DDSolid.h:167
DDTrap::y1
double y1(void) const
Half-length along y of the face at -pDz.
Definition: DDSolid.cc:141
magneticfield::BaseVolumeHandle::isIronFlag
bool isIronFlag
Definition: BaseVolumeHandle.h:149
SurfaceOrientation::inner
Definition: Surface.h:19
DDTrap::theta
double theta(void) const
Polar angle of the line joining the centres of the faces at -/+pDz.
Definition: DDSolid.cc:137
DDSolidShape::ddtrunctubs
magneticfield::BaseVolumeHandle::theRMax
double theRMax
Definition: BaseVolumeHandle.h:133
DDBase::name
const N & name() const
Definition: DDBase.h:59
DDAxes::z
magneticfield::BaseVolumeHandle::GlobalPoint
Surface::GlobalPoint GlobalPoint
Definition: BaseVolumeHandle.h:24
DDSolid.h
theta
Geom::Theta< T > theta() const
Definition: Basic3DVectorLD.h:150
DDExpandedView
Provides an exploded view of the detector (tree-view)
Definition: DDExpandedView.h:41
Point3DBase< float, GlobalTag >
magneticfield::BaseVolumeHandle::masterSector
int masterSector
The sector for which an interpolator for this class of volumes should be built.
Definition: BaseVolumeHandle.h:111
DDValue.h
DDPseudoTrap::y1
double y1(void) const
half length along y on -z
Definition: DDSolid.cc:196
DDSolidShapesName::name
static const char *const name(DDSolidShape s)
Definition: DDSolidShapes.h:33
DD3Vector
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DD3Vector
Definition: PGeometricDetBuilder.cc:20
DDTubs::rIn
double rIn(void) const
Definition: DDSolid.cc:456
DDTubs::deltaPhi
double deltaPhi(void) const
Definition: DDSolid.cc:462
Cone.h
testProducerWithPsetDescEmpty_cfi.y2
y2
Definition: testProducerWithPsetDescEmpty_cfi.py:30
MagGeoBuilderFromDDD::volumeHandle::buildBox
void buildBox(double halfX, double halfY, double halfZ)
magneticfield::BaseVolumeHandle::center_
GlobalPoint center_
Definition: BaseVolumeHandle.h:142
DDTrap
Interface to a Trapezoid.
Definition: DDSolid.h:88
MagGeoBuilderFromDDD::volumeHandle::buildPseudoTrap
void buildPseudoTrap(double x1, double x2, double y1, double y2, double halfZ, double radius, bool atMinusZ)
MagGeoBuilderFromDDD::volumeHandle::buildCons
void buildCons(double zhalf, double rInMinusZ, double rOutMinusZ, double rInPlusZ, double rOutPlusZ, double startPhi, double deltaPhi)
DDTruncTubs::rIn
double rIn(void) const
inner radius
Definition: DDSolid.cc:168
GeantUnits.h
DDTubs::rOut
double rOut(void) const
Definition: DDSolid.cc:458
DDSolidShape::ddcons
DDMaterial.h
LocalVector
Local3DVector LocalVector
Definition: LocalVector.h:12
DDPseudoTrap::x1
double x1(void) const
half length along x on -z
Definition: DDSolid.cc:192
DDBox::halfX
double halfX(void) const
Definition: DDSolid.cc:212
DDName::name
const std::string & name() const
Returns the name.
Definition: DDName.cc:41
DDTranslation.h
DDTubs::startPhi
double startPhi(void) const
Definition: DDSolid.cc:460
magneticfield::BaseVolumeHandle
Definition: BaseVolumeHandle.h:22
DDBox::halfY
double halfY(void) const
Definition: DDSolid.cc:214
get
#define get
magneticfield::BaseVolumeHandle::copyno
unsigned short copyno
copy number
Definition: BaseVolumeHandle.h:72
Vector3DBase::dot
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
DDAxes::phi
MagGeoBuilderFromDDD::volumeHandle::buildTubs
void buildTubs(double zhalf, double rIn, double rOut, double startPhi, double deltaPhi)
DDPseudoTrap::halfZ
double halfZ(void) const
half of the z-Axis
Definition: DDSolid.cc:190
DDTrap::y2
double y2(void) const
Half-length along y of the face at +pDz.
Definition: DDSolid.cc:149
std
Definition: JetResolutionObject.h:76
DDSolidShape::ddbox
volumeHandle.h
magneticfield::BaseVolumeHandle::name
std::string name
Name of the volume.
Definition: BaseVolumeHandle.h:63
DDTruncTubs::deltaPhi
double deltaPhi(void) const
angular span of the tube-section
Definition: DDSolid.cc:174
DDTruncTubs::cutAtDelta
double cutAtDelta(void) const
truncation at end of the tube-section
Definition: DDSolid.cc:178
MagGeoBuilderFromDDD::volumeHandle::buildTrap
void buildTrap(double x1, double x2, double x3, double x4, double y1, double y2, double theta, double phi, double halfZ, double alpha1, double alpha2)
magneticfield::BaseVolumeHandle::surfaces
RCPS surfaces[6]
Definition: BaseVolumeHandle.h:120
CosmicsPD_Skims.radius
radius
Definition: CosmicsPD_Skims.py:135
DDTruncTubs::rOut
double rOut(void) const
outer radius
Definition: DDSolid.cc:170
SurfaceOrientation::GlobalFace
GlobalFace
Definition: Surface.h:19
DDExpandedView::logicalPart
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the expanded-view.
Definition: DDExpandedView.cc:42
DDTrap::x1
double x1(void) const
Half-length along x of the side at y=-pDy1 of the face at -pDz.
Definition: DDSolid.cc:143
geant_units::operators::convertMmToCm
constexpr NumType convertMmToCm(NumType millimeters)
Definition: GeantUnits.h:63
SurfaceOrientation::phiminus
Definition: Surface.h:19
DDCons::rOutPlusZ
double rOutPlusZ(void) const
Definition: DDSolid.cc:422
DDCons::zhalf
double zhalf(void) const
Definition: DDSolid.cc:414
mps_fire.result
result
Definition: mps_fire.py:311
DDTrap::alpha1
double alpha1(void) const
Angle with respect to the y axis from the centre of the side at y=-pDy1 to the centre at y=+pDy1 of t...
Definition: DDSolid.cc:147
DDExpandedView::translation
const DDTranslation & translation() const
The absolute translation of the current node.
Definition: DDExpandedView.cc:46
DDTrap::alpha2
double alpha2(void) const
Angle with respect to the y axis from the centre of the side at y=-pDy2 to the centre at y=+pDy2 of t...
Definition: DDSolid.cc:155
SurfaceOrientation
Definition: Surface.h:17
DDLogicalPart::solid
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
Definition: DDLogicalPart.cc:120
PV3DBase::perp
T perp() const
Definition: PV3DBase.h:69
PV3DBase::phi
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
DDTrap::x3
double x3(void) const
Half-length along x of the side at y=-pDy2 of the face at +pDz.
Definition: DDSolid.cc:151
DDTubs
Definition: DDSolid.h:266