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/SystemOfUnits.h>
20 
21 #include <string>
22 #include <iterator>
23 #include <iomanip>
24 #include <iostream>
25 
26 using namespace SurfaceOrientation;
27 using namespace std;
28 
30 
31 namespace {
32  // Old DD returns lengths in mm, but CMS code uses cm
33  template <class NumType>
34  inline constexpr NumType convertUnits(NumType millimeters) // Millimeters -> centimeters
35  {
36  return (geant_units::operators::convertMmToCm(millimeters));
37  }
38 } // namespace
39 
40 MagGeoBuilderFromDDD::volumeHandle::volumeHandle(const DDExpandedView &fv, bool expand2Pi, bool debugVal)
41  : magneticfield::BaseVolumeHandle(expand2Pi, debugVal) {
42  name = fv.logicalPart().name().name();
43  copyno = fv.copyno();
44  solid = fv.logicalPart().solid();
45  center_ =
46  GlobalPoint(fv.translation().x() / CLHEP::cm, fv.translation().y() / CLHEP::cm, fv.translation().z() / CLHEP::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 = std::stoul(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"
double radius(void) const
radius of the cut-out (neg.) or rounding (pos.)
Definition: DDSolid.cc:200
MagGeoBuilderFromDDD::volumeHandle volumeHandle
Local3DVector LocalVector
Definition: LocalVector.h:12
T perp() const
Definition: PV3DBase.h:69
A truncated tube section.
Definition: DDSolid.h:139
double halfZ(void) const
Definition: DDSolid.cc:216
const DDTranslation & translation() const
The absolute translation of the current node.
void buildTrap(double x1, double x2, double x3, double x4, double y1, double y2, double theta, double phi, double halfZ, double alpha1, double alpha2)
double startPhi(void) const
angular start of the tube-section
Definition: DDSolid.cc:172
double rInPlusZ(void) const
Definition: DDSolid.cc:420
double rIn(void) const
Definition: DDSolid.cc:456
const DDLogicalPart & logicalPart() const
The logical-part of the current node in the expanded-view.
PreciseFloatType< T, U >::Type dot(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:99
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
volumeHandle(const DDExpandedView &fv, bool expand2Pi=false, bool debugVal=false)
Definition: volumeHandle.cc:40
std::vector< VolumeSide > sides() const override
The surfaces and they orientation, as required to build a MagVolume.
double rOutMinusZ(void) const
Definition: DDSolid.cc:418
double rOut(void) const
Definition: DDSolid.cc:458
double cutAtStart(void) const
truncation at begin of the tube-section
Definition: DDSolid.cc:176
double phi(void) const
Azimuthal angle of the line joining the centres of the faces at -/+pDz.
Definition: DDSolid.cc:139
double halfZ(void) const
half of the z-Axis
Definition: DDSolid.cc:135
bool cutInside(void) const
true, if truncation is on the inner side of the tube-section
Definition: DDSolid.cc:180
const GlobalPoint & center() const
Return the center of the volume.
unsigned short volumeno
volume number
unsigned short copyno
copy number
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > DD3Vector
double halfZ(void) const
half of the z-Axis
Definition: DDSolid.cc:190
double x3(void) const
Half-length along x of the side at y=-pDy2 of the face at +pDz.
Definition: DDSolid.cc:151
static const char *const name(DDSolidShape s)
Definition: DDSolidShapes.h:33
double x2(void) const
half length along x on +z
Definition: DDSolid.cc:194
const std::string & name() const
Returns the name.
Definition: DDName.cc:41
double phiFrom(void) const
Definition: DDSolid.cc:424
double y2(void) const
half length along y on +z
Definition: DDSolid.cc:198
void buildBox(double halfX, double halfY, double halfZ)
Surface::GlobalPoint GlobalPoint
double zHalf(void) const
half of the z-Axis
Definition: DDSolid.cc:166
Interface to a Trapezoid.
Definition: DDSolid.h:88
double deltaPhi(void) const
Definition: DDSolid.cc:426
int masterSector
The sector for which an interpolator for this class of volumes should be built.
void buildCons(double zhalf, double rInMinusZ, double rOutMinusZ, double rInPlusZ, double rOutPlusZ, double startPhi, double deltaPhi)
double y2(void) const
Half-length along y of the face at +pDz.
Definition: DDSolid.cc:149
void referencePlane(const DDExpandedView &fv)
double rInMinusZ(void) const
Definition: DDSolid.cc:416
double x4(void) const
Half-length along x of the side at y=+pDy2 of the face at +pDz.
Definition: DDSolid.cc:153
DDSolidShape shape(void) const
The type of the solid.
Definition: DDSolid.cc:123
const DDMaterial & material(void) const
Returns a reference object of the material this LogicalPart is made of.
double rIn(void) const
inner radius
Definition: DDSolid.cc:168
const N & name() const
Definition: DDBase.h:58
int copyno() const
Copy number associated with the current node.
double rOutPlusZ(void) const
Definition: DDSolid.cc:422
double deltaPhi(void) const
Definition: DDSolid.cc:462
constexpr NumType convertMmToCm(NumType millimeters)
Definition: angle_units.h:44
double x2(void) const
Half-length along x of the side at y=+pDy1 of the face at -pDz.
Definition: DDSolid.cc:145
double rOut(void) const
outer radius
Definition: DDSolid.cc:170
void buildTruncTubs(double zhalf, double rIn, double rOut, double startPhi, double deltaPhi, double cutAtStart, double cutAtDelta, bool cutInside)
Interface to a Box.
Definition: DDSolid.h:167
void buildPseudoTrap(double x1, double x2, double y1, double y2, double halfZ, double radius, bool atMinusZ)
double x1(void) const
half length along x on -z
Definition: DDSolid.cc:192
std::string name
Name of the volume.
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
const DDRotationMatrix & rotation() const
The absolute rotation of the current node.
DDSolidShape shape() const override
Shape of the solid.
Definition: volumeHandle.h:32
const DDSolid & solid(void) const
Returns a reference object of the solid being the shape of this LogicalPart.
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
double halfX(void) const
Definition: DDSolid.cc:212
double zhalf(void) const
Definition: DDSolid.cc:454
double startPhi(void) const
Definition: DDSolid.cc:460
double halfY(void) const
Definition: DDSolid.cc:214
double deltaPhi(void) const
angular span of the tube-section
Definition: DDSolid.cc:174
double y1(void) const
Half-length along y of the face at -pDz.
Definition: DDSolid.cc:141
bool atMinusZ(void) const
true, if cut-out or rounding is on the -z side
Definition: DDSolid.cc:202
Provides an exploded view of the detector (tree-view)
double y1(void) const
half length along y on -z
Definition: DDSolid.cc:196
double theta(void) const
Polar angle of the line joining the centres of the faces at -/+pDz.
Definition: DDSolid.cc:137
double cutAtDelta(void) const
truncation at end of the tube-section
Definition: DDSolid.cc:178
void buildTubs(double zhalf, double rIn, double rOut, double startPhi, double deltaPhi)
std::string magFile
Name of magnetic field table file.
Global3DVector GlobalVector
Definition: GlobalVector.h:10
double zhalf(void) const
Definition: DDSolid.cc:414
double x1(void) const
Half-length along x of the side at y=-pDy1 of the face at -pDz.
Definition: DDSolid.cc:143