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 
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_ = GlobalPoint(fv.translation().x() / cm, fv.translation().y() / cm, fv.translation().z() / cm);
46 
47  // ASSUMPTION: volume names ends with "_NUM" where NUM is the volume number
48  string volName = name;
49  volName.erase(0, volName.rfind('_') + 1);
50  volumeno = std::stoul(volName);
51 
52  for (int i = 0; i < 6; ++i) {
53  isAssigned[i] = false;
54  }
55 
56  if (debug) {
57  cout.precision(7);
58  }
59 
60  referencePlane(fv);
61 
62  if (solid.shape() == DDSolidShape::ddbox) {
63  DDBox box(solid);
64  double halfX = convertUnits(box.halfX());
65  double halfY = convertUnits(box.halfY());
66  double halfZ = convertUnits(box.halfZ());
67  buildBox(halfX, halfY, halfZ);
68  } else if (solid.shape() == DDSolidShape::ddtrap) {
69  DDTrap trap(solid);
70  double x1 = convertUnits(trap.x1());
71  double x2 = convertUnits(trap.x2());
72  double x3 = convertUnits(trap.x3());
73  double x4 = convertUnits(trap.x4());
74  double y1 = convertUnits(trap.y1());
75  double y2 = convertUnits(trap.y2());
76  double theta = trap.theta();
77  double phi = trap.phi();
78  double halfZ = convertUnits(trap.halfZ());
79  double alpha1 = trap.alpha1();
80  double alpha2 = trap.alpha2();
81  buildTrap(x1, x2, x3, x4, y1, y2, theta, phi, halfZ, alpha1, alpha2);
82  } else if (solid.shape() == DDSolidShape::ddcons) {
83  DDCons cons(solid);
84  double zhalf = convertUnits(cons.zhalf());
85  double rInMinusZ = convertUnits(cons.rInMinusZ());
86  double rOutMinusZ = convertUnits(cons.rOutMinusZ());
87  double rInPlusZ = convertUnits(cons.rInPlusZ());
88  double rOutPlusZ = convertUnits(cons.rOutPlusZ());
89  double startPhi = cons.phiFrom();
90  double deltaPhi = cons.deltaPhi();
91  buildCons(zhalf, rInMinusZ, rOutMinusZ, rInPlusZ, rOutPlusZ, startPhi, deltaPhi);
92  } else if (solid.shape() == DDSolidShape::ddtubs) {
93  DDTubs tubs(solid);
94  double zhalf = convertUnits(tubs.zhalf());
95  double rIn = convertUnits(tubs.rIn());
96  double rOut = convertUnits(tubs.rOut());
97  double startPhi = tubs.startPhi();
98  double deltaPhi = tubs.deltaPhi();
99  buildTubs(zhalf, rIn, rOut, startPhi, deltaPhi);
100  } else if (solid.shape() == DDSolidShape::ddpseudotrap) {
101  DDPseudoTrap ptrap(solid);
102  double x1 = convertUnits(ptrap.x1());
103  double x2 = convertUnits(ptrap.x2());
104  double y1 = convertUnits(ptrap.y1());
105  double y2 = convertUnits(ptrap.y2());
106  double halfZ = convertUnits(ptrap.halfZ());
107  double radius = convertUnits(ptrap.radius());
108  bool atMinusZ = ptrap.atMinusZ();
109  buildPseudoTrap(x1, x2, y1, y2, halfZ, radius, atMinusZ);
110  } else if (solid.shape() == DDSolidShape::ddtrunctubs) {
111  DDTruncTubs tubs(solid);
112  double zhalf = convertUnits(tubs.zHalf()); // half of the z-Axis
113  double rIn = convertUnits(tubs.rIn()); // inner radius
114  double rOut = convertUnits(tubs.rOut()); // outer radius
115  double startPhi = tubs.startPhi(); // angular start of the tube-section
116  double deltaPhi = tubs.deltaPhi(); // angular span of the tube-section
117  double cutAtStart = convertUnits(tubs.cutAtStart()); // truncation at begin of the tube-section
118  double cutAtDelta = convertUnits(tubs.cutAtDelta()); // truncation at end of the tube-section
119  bool cutInside = tubs.cutInside(); // true, if truncation is on the inner side of the tube-section
120  buildTruncTubs(zhalf, rIn, rOut, startPhi, deltaPhi, cutAtStart, cutAtDelta, cutInside);
121  } else {
122  cout << "volumeHandle ctor: Unexpected solid: " << DDSolidShapesName::name(solid.shape()) << endl;
123  }
124 
125  // NOTE: Table name and master sector are no longer taken from xml!
126  // DDsvalues_type sv(fv.mergedSpecifics());
127 
128  // { // Extract the name of associated field file.
129  // std::vector<std::string> temp;
130  // std::string pname = "table";
131  // DDValue val(pname);
132  // DDsvalues_type sv(fv.mergedSpecifics());
133  // if (DDfetch(&sv,val)) {
134  // temp = val.strings();
135  // if (temp.size() != 1) {
136  // cout << "*** WARNING: volume has > 1 SpecPar " << pname << endl;
137  // }
138  // magFile = temp[0];
139 
140  // string find="[copyNo]";
141  // std::size_t j;
142  // for ( ; (j = magFile.find(find)) != string::npos ; ) {
143  // stringstream conv;
144  // conv << setfill('0') << setw(2) << copyno;
145  // string repl;
146  // conv >> repl;
147  // magFile.replace(j, find.length(), repl);
148  // }
149 
150  // } else {
151  // cout << "*** WARNING: volume does not have a SpecPar " << pname << endl;
152  // cout << " DDsvalues_type: " << fv.mergedSpecifics() << endl;
153  // }
154  // }
155 
156  // { // Extract the number of the master sector.
157  // std::vector<double> temp;
158  // const std::string pname = "masterSector";
159  // DDValue val(pname);
160  // if (DDfetch(&sv,val)) {
161  // temp = val.doubles();
162  // if (temp.size() != 1) {
163  // cout << "*** WARNING: volume has > 1 SpecPar " << pname << endl;
164  // }
165  // masterSector = int(temp[0]+.5);
166  // } else {
167  // if (MagGeoBuilderFromDDD::debug) {
168  // cout << "Volume does not have a SpecPar " << pname
169  // << " using: " << copyno << endl;
170  // cout << " DDsvalues_type: " << fv.mergedSpecifics() << endl;
171  // }
172  // masterSector = copyno;
173  // }
174  // }
175 
176  // Get material for this volume
177  if (fv.logicalPart().material().name().name() == "Iron")
178  isIronFlag = true;
179 
180  if (debug) {
181  cout << " RMin = " << theRMin << endl;
182  cout << " RMax = " << theRMax << endl;
183 
184  if (theRMin < 0 || theRN < theRMin || theRMax < theRN)
185  cout << "*** WARNING: wrong RMin/RN/RMax , shape: " << DDSolidShapesName::name(shape()) << endl;
186 
187  cout << "Summary: " << name << " " << copyno << " Shape= " << DDSolidShapesName::name(shape()) << " trasl "
188  << center() << " R " << center().perp() << " phi " << center().phi() << " magFile " << magFile
189  << " Material= " << fv.logicalPart().material().name() << " isIron= " << isIronFlag
190  << " masterSector= " << masterSector << std::endl;
191 
192  cout << " Orientation of surfaces:";
193  std::string sideName[3] = {"positiveSide", "negativeSide", "onSurface"};
194  for (int i = 0; i < 6; ++i) {
195  cout << " " << i << ":" << sideName[surfaces[i]->side(center_, 0.3)];
196  }
197  cout << endl;
198  }
199 }
200 
202  // The refPlane is the "main plane" for the solid. It corresponds to the
203  // x,y plane in the DDD local frame, and defines a frame where the local
204  // coordinates are the same as in DDD.
205  // In the geometry version 85l_030919, this plane is normal to the
206  // beam line for all volumes but pseudotraps, so that global R is along Y,
207  // global phi is along -X and global Z along Z:
208  //
209  // Global(for vol at pi/2) Local
210  // +R (+Y) +Y
211  // +phi(-X) -X
212  // +Z +Z
213  //
214  // For pseudotraps the refPlane is parallel to beam line and global R is
215  // along Z, global phi is along +-X and and global Z along Y:
216  //
217  // Global(for vol at pi/2) Local
218  // +R (+Y) +Z
219  // +phi(-X) +X
220  // +Z +Y
221  //
222  // Note that the frame is centered in the DDD volume center, which is
223  // inside the volume for DDD boxes and (pesudo)trapezoids, on the beam line
224  // for tubs, cons and trunctubs.
225 
226  // In geometry version 1103l, trapezoids have X and Z in the opposite direction
227  // than the above. Boxes are either oriented as described above or in some case
228  // have opposite direction for Y and X.
229 
230  // The global position
231  Surface::PositionType &posResult = center_;
232 
233  // The reference plane rotation
234  DD3Vector x, y, z;
235  fv.rotation().GetComponents(x, y, z);
236  if (debug) {
237  if (x.Cross(y).Dot(z) < 0.5) {
238  cout << "*** WARNING: Rotation is not RH " << endl;
239  }
240  }
241 
242  // The global rotation
243  Surface::RotationType rotResult(float(x.X()),
244  float(x.Y()),
245  float(x.Z()),
246  float(y.X()),
247  float(y.Y()),
248  float(y.Z()),
249  float(z.X()),
250  float(z.Y()),
251  float(z.Z()));
252 
253  refPlane = new GloballyPositioned<float>(posResult, rotResult);
254 
255  // Check correct orientation
256  if (debug) {
257  cout << "Refplane pos " << refPlane->position() << endl;
258 
259  // See comments above for the conventions for orientation.
260  LocalVector globalZdir(0., 0., 1.); // Local direction of the axis along global Z
261  if (solid.shape() == DDSolidShape::ddpseudotrap) {
262  globalZdir = LocalVector(0., 1., 0.);
263  }
264  if (refPlane->toGlobal(globalZdir).z() < 0.) {
265  globalZdir = -globalZdir;
266  }
267 
268  float chk = refPlane->toGlobal(globalZdir).dot(GlobalVector(0, 0, 1));
269  if (chk < .999)
270  cout << "*** WARNING RefPlane check failed!***" << chk << endl;
271  }
272 }
273 
274 std::vector<VolumeSide> MagGeoBuilderFromDDD::volumeHandle::sides() const {
275  std::vector<VolumeSide> result;
276  for (int i = 0; i < 6; ++i) {
277  // If this is just a master volume out of wich a 2pi volume
278  // should be built (e.g. central cylinder), skip the phi boundaries.
279  if (expand && (i == phiplus || i == phiminus))
280  continue;
281 
282  // FIXME: Skip null inner degenerate cylindrical surface
283  if (solid.shape() == DDSolidShape::ddtubs && i == SurfaceOrientation::inner && theRMin < 0.001)
284  continue;
285 
286  ReferenceCountingPointer<Surface> s = const_cast<Surface *>(surfaces[i].get());
287  result.push_back(VolumeSide(s, GlobalFace(i), surfaces[i]->side(center_, 0.3)));
288  }
289  return result;
290 }
291 
293 using namespace magneticfield;
294 
295 #include "buildBox.icc"
296 #include "buildTrap.icc"
297 #include "buildTubs.icc"
298 #include "buildCons.icc"
299 #include "buildPseudoTrap.icc"
300 #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:59
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)
Geom::Theta< T > theta() const
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