Geant4  10.00.p02
G4TwistedTubs.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4TwistedTubs.cc 81641 2014-06-04 09:11:38Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class source file
32 //
33 //
34 // G4TwistTubsSide.cc
35 //
36 // Author:
37 // 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp)
38 //
39 // History:
40 // 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
41 // from original version in Jupiter-2.5.02 application.
42 // --------------------------------------------------------------------
43 
44 #include "G4TwistedTubs.hh"
45 
46 #include "G4PhysicalConstants.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4GeometryTolerance.hh"
49 #include "G4VoxelLimits.hh"
50 #include "G4AffineTransform.hh"
51 #include "G4SolidExtentList.hh"
52 #include "G4ClippablePolygon.hh"
53 #include "G4VPVParameterisation.hh"
54 #include "meshdefs.hh"
55 
56 #include "G4VGraphicsScene.hh"
57 #include "G4Polyhedron.hh"
58 #include "G4VisExtent.hh"
59 
60 #include "Randomize.hh"
61 
62 //=====================================================================
63 //* constructors ------------------------------------------------------
64 
66  G4double twistedangle,
67  G4double endinnerrad,
68  G4double endouterrad,
69  G4double halfzlen,
70  G4double dphi)
71  : G4VSolid(pname), fDPhi(dphi),
72  fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
73  fFormerTwisted(0), fInnerHype(0), fOuterHype(0),
74  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
75 {
76  if (endinnerrad < DBL_MIN)
77  {
78  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
79  FatalErrorInArgument, "Invalid end-inner-radius!");
80  }
81 
82  G4double sinhalftwist = std::sin(0.5 * twistedangle);
83 
84  G4double endinnerradX = endinnerrad * sinhalftwist;
85  G4double innerrad = std::sqrt( endinnerrad * endinnerrad
86  - endinnerradX * endinnerradX );
87 
88  G4double endouterradX = endouterrad * sinhalftwist;
89  G4double outerrad = std::sqrt( endouterrad * endouterrad
90  - endouterradX * endouterradX );
91 
92  // temporary treatment!!
93  SetFields(twistedangle, innerrad, outerrad, -halfzlen, halfzlen);
95 }
96 
98  G4double twistedangle,
99  G4double endinnerrad,
100  G4double endouterrad,
101  G4double halfzlen,
102  G4int nseg,
103  G4double totphi)
104  : G4VSolid(pname),
105  fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
106  fFormerTwisted(0), fInnerHype(0), fOuterHype(0),
107  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
108 {
109 
110  if (!nseg)
111  {
112  std::ostringstream message;
113  message << "Invalid number of segments." << G4endl
114  << " nseg = " << nseg;
115  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
116  FatalErrorInArgument, message);
117  }
118  if (totphi == DBL_MIN || endinnerrad < DBL_MIN)
119  {
120  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
121  FatalErrorInArgument, "Invalid total-phi or end-inner-radius!");
122  }
123 
124  G4double sinhalftwist = std::sin(0.5 * twistedangle);
125 
126  G4double endinnerradX = endinnerrad * sinhalftwist;
127  G4double innerrad = std::sqrt( endinnerrad * endinnerrad
128  - endinnerradX * endinnerradX );
129 
130  G4double endouterradX = endouterrad * sinhalftwist;
131  G4double outerrad = std::sqrt( endouterrad * endouterrad
132  - endouterradX * endouterradX );
133 
134  // temporary treatment!!
135  fDPhi = totphi / nseg;
136  SetFields(twistedangle, innerrad, outerrad, -halfzlen, halfzlen);
137  CreateSurfaces();
138 }
139 
141  G4double twistedangle,
142  G4double innerrad,
143  G4double outerrad,
144  G4double negativeEndz,
145  G4double positiveEndz,
146  G4double dphi)
147  : G4VSolid(pname), fDPhi(dphi),
148  fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
149  fFormerTwisted(0), fInnerHype(0), fOuterHype(0),
150  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
151 {
152  if (innerrad < DBL_MIN)
153  {
154  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
155  FatalErrorInArgument, "Invalid end-inner-radius!");
156  }
157 
158  SetFields(twistedangle, innerrad, outerrad, negativeEndz, positiveEndz);
159  CreateSurfaces();
160 }
161 
163  G4double twistedangle,
164  G4double innerrad,
165  G4double outerrad,
166  G4double negativeEndz,
167  G4double positiveEndz,
168  G4int nseg,
169  G4double totphi)
170  : G4VSolid(pname),
171  fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
172  fFormerTwisted(0), fInnerHype(0), fOuterHype(0),
173  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
174 {
175  if (!nseg)
176  {
177  std::ostringstream message;
178  message << "Invalid number of segments." << G4endl
179  << " nseg = " << nseg;
180  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
181  FatalErrorInArgument, message);
182  }
183  if (totphi == DBL_MIN || innerrad < DBL_MIN)
184  {
185  G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
186  FatalErrorInArgument, "Invalid total-phi or end-inner-radius!");
187  }
188 
189  fDPhi = totphi / nseg;
190  SetFields(twistedangle, innerrad, outerrad, negativeEndz, positiveEndz);
191  CreateSurfaces();
192 }
193 
194 //=====================================================================
195 //* Fake default constructor ------------------------------------------
196 
198  : G4VSolid(a), fPhiTwist(0.), fInnerRadius(0.), fOuterRadius(0.), fDPhi(0.),
199  fZHalfLength(0.), fInnerStereo(0.), fOuterStereo(0.), fTanInnerStereo(0.),
200  fTanOuterStereo(0.), fKappa(0.), fInnerRadius2(0.), fOuterRadius2(0.),
201  fTanInnerStereo2(0.), fTanOuterStereo2(0.), fLowerEndcap(0), fUpperEndcap(0),
202  fLatterTwisted(0), fFormerTwisted(0), fInnerHype(0), fOuterHype(0),
203  fCubicVolume(0.), fSurfaceArea(0.), fpPolyhedron(0)
204 {
205  fEndZ[0] = 0.; fEndZ[1] = 0.;
206  fEndInnerRadius[0] = 0.; fEndInnerRadius[1] = 0.;
207  fEndOuterRadius[0] = 0.; fEndOuterRadius[1] = 0.;
208  fEndPhi[0] = 0.; fEndPhi[1] = 0.;
209  fEndZ2[0] = 0.; fEndZ2[1] = 0.;
210 }
211 
212 //=====================================================================
213 //* destructor --------------------------------------------------------
214 
216 {
217  if (fLowerEndcap) { delete fLowerEndcap; }
218  if (fUpperEndcap) { delete fUpperEndcap; }
219  if (fLatterTwisted) { delete fLatterTwisted; }
220  if (fFormerTwisted) { delete fFormerTwisted; }
221  if (fInnerHype) { delete fInnerHype; }
222  if (fOuterHype) { delete fOuterHype; }
223  if (fpPolyhedron) { delete fpPolyhedron; }
224 }
225 
226 //=====================================================================
227 //* Copy constructor --------------------------------------------------
228 
230  : G4VSolid(rhs), fPhiTwist(rhs.fPhiTwist),
231  fInnerRadius(rhs.fInnerRadius), fOuterRadius(rhs.fOuterRadius),
232  fDPhi(rhs.fDPhi), fZHalfLength(rhs.fZHalfLength),
233  fInnerStereo(rhs.fInnerStereo), fOuterStereo(rhs.fOuterStereo),
234  fTanInnerStereo(rhs.fTanInnerStereo), fTanOuterStereo(rhs.fTanOuterStereo),
235  fKappa(rhs.fKappa), fInnerRadius2(rhs.fInnerRadius2),
236  fOuterRadius2(rhs.fOuterRadius2), fTanInnerStereo2(rhs.fTanInnerStereo2),
237  fTanOuterStereo2(rhs.fTanOuterStereo2),
238  fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0), fFormerTwisted(0),
239  fInnerHype(0), fOuterHype(0),
240  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
241  fpPolyhedron(0), fLastInside(rhs.fLastInside), fLastNormal(rhs.fLastNormal),
242  fLastDistanceToIn(rhs.fLastDistanceToIn),
243  fLastDistanceToOut(rhs.fLastDistanceToOut),
244  fLastDistanceToInWithV(rhs.fLastDistanceToInWithV),
245  fLastDistanceToOutWithV(rhs.fLastDistanceToOutWithV)
246 {
247  for (size_t i=0; i<2; ++i)
248  {
249  fEndZ[i] = rhs.fEndZ[i];
250  fEndInnerRadius[i] = rhs.fEndInnerRadius[i];
251  fEndOuterRadius[i] = rhs.fEndOuterRadius[i];
252  fEndPhi[i] = rhs.fEndPhi[i];
253  fEndZ2[i] = rhs.fEndZ2[i];
254  }
255  CreateSurfaces();
257 }
258 
259 
260 //=====================================================================
261 //* Assignment operator -----------------------------------------------
262 
264 {
265  // Check assignment to self
266  //
267  if (this == &rhs) { return *this; }
268 
269  // Copy base class data
270  //
271  G4VSolid::operator=(rhs);
272 
273  // Copy data
274  //
275  fPhiTwist= rhs.fPhiTwist;
291 
292  for (size_t i=0; i<2; ++i)
293  {
294  fEndZ[i] = rhs.fEndZ[i];
295  fEndInnerRadius[i] = rhs.fEndInnerRadius[i];
296  fEndOuterRadius[i] = rhs.fEndOuterRadius[i];
297  fEndPhi[i] = rhs.fEndPhi[i];
298  fEndZ2[i] = rhs.fEndZ2[i];
299  }
300 
301  CreateSurfaces();
303 
304  return *this;
305 }
306 
307 //=====================================================================
308 //* ComputeDimensions -------------------------------------------------
309 
311  const G4int /* n */ ,
312  const G4VPhysicalVolume* /* pRep */ )
313 {
314  G4Exception("G4TwistedTubs::ComputeDimensions()",
315  "GeomSolids0001", FatalException,
316  "G4TwistedTubs does not support Parameterisation.");
317 }
318 
319 
320 //=====================================================================
321 //* CalculateExtent ---------------------------------------------------
322 
324  const G4VoxelLimits &voxelLimit,
325  const G4AffineTransform &transform,
326  G4double &min,
327  G4double &max ) const
328 {
329 
330  G4SolidExtentList extentList( axis, voxelLimit );
331  G4double maxEndOuterRad = (fEndOuterRadius[0] > fEndOuterRadius[1] ?
333  G4double maxEndInnerRad = (fEndInnerRadius[0] > fEndInnerRadius[1] ?
335  G4double maxphi = (std::fabs(fEndPhi[0]) > std::fabs(fEndPhi[1]) ?
336  std::fabs(fEndPhi[0]) : std::fabs(fEndPhi[1]));
337 
338  //
339  // Choose phi size of our segment(s) based on constants as
340  // defined in meshdefs.hh
341  //
342  // G4int numPhi = kMaxMeshSections;
343  G4double sigPhi = 2*maxphi + fDPhi;
344  G4double rFudge = 1.0/std::cos(0.5*sigPhi);
345  G4double fudgeEndOuterRad = rFudge * maxEndOuterRad;
346 
347  //
348  // We work around in phi building polygons along the way.
349  // As a reasonable compromise between accuracy and
350  // complexity (=cpu time), the following facets are chosen:
351  //
352  // 1. If fOuterRadius/maxEndOuterRad > 0.95, approximate
353  // the outer surface as a cylinder, and use one
354  // rectangular polygon (0-1) to build its mesh.
355  //
356  // Otherwise, use two trapazoidal polygons that
357  // meet at z = 0 (0-4-1)
358  //
359  // 2. If there is no inner surface, then use one
360  // polygon for each entire endcap. (0) and (1)
361  //
362  // Otherwise, use a trapazoidal polygon for each
363  // phi segment of each endcap. (0-2) and (1-3)
364  //
365  // 3. For the inner surface, if fInnerRadius/maxEndInnerRad > 0.95,
366  // approximate the inner surface as a cylinder of
367  // radius fInnerRadius and use one rectangular polygon
368  // to build each phi segment of its mesh. (2-3)
369  //
370  // Otherwise, use one rectangular polygon centered
371  // at z = 0 (5-6) and two connecting trapazoidal polygons
372  // for each phi segment (2-5) and (3-6).
373  //
374 
375  G4bool splitOuter = (fOuterRadius/maxEndOuterRad < 0.95);
376  G4bool splitInner = (fInnerRadius/maxEndInnerRad < 0.95);
377 
378  //
379  // Vertex assignments (v and w arrays)
380  // [0] and [1] are mandatory
381  // the rest are optional
382  //
383  // + -
384  // [0]------[4]------[1] <--- outer radius
385  // | |
386  // | |
387  // [2]---[5]---[6]---[3] <--- inner radius
388  //
389 
390  G4ClippablePolygon endPoly1, endPoly2;
391 
392  G4double phimax = maxphi + 0.5*fDPhi;
393  if ( phimax > pi/2) { phimax = pi-phimax; }
394  G4double phimin = - phimax;
395 
396  G4ThreeVector v0, v1, v2, v3, v4, v5, v6; // -ve phi verticies for polygon
397  G4ThreeVector w0, w1, w2, w3, w4, w5, w6; // +ve phi verticies for polygon
398 
399  //
400  // decide verticies of -ve phi boundary
401  //
402 
403  G4double cosPhi = std::cos(phimin);
404  G4double sinPhi = std::sin(phimin);
405 
406  // Outer hyperbolic surface
407 
408  v0 = transform.TransformPoint(
409  G4ThreeVector(fudgeEndOuterRad*cosPhi, fudgeEndOuterRad*sinPhi,
410  + fZHalfLength));
411  v1 = transform.TransformPoint(
412  G4ThreeVector(fudgeEndOuterRad*cosPhi, fudgeEndOuterRad*sinPhi,
413  - fZHalfLength));
414  if (splitOuter)
415  {
416  v4 = transform.TransformPoint(
417  G4ThreeVector(fudgeEndOuterRad*cosPhi, fudgeEndOuterRad*sinPhi, 0));
418  }
419 
420  // Inner hyperbolic surface
421 
422  G4double zInnerSplit = 0.;
423  if (splitInner)
424  {
425  v2 = transform.TransformPoint(
426  G4ThreeVector(maxEndInnerRad*cosPhi, maxEndInnerRad*sinPhi,
427  + fZHalfLength));
428  v3 = transform.TransformPoint(
429  G4ThreeVector(maxEndInnerRad*cosPhi, maxEndInnerRad*sinPhi,
430  - fZHalfLength));
431 
432  // Find intersection of tangential line of inner
433  // surface at z = fZHalfLength and line r=fInnerRadius.
435  G4double dz = maxEndInnerRad;
436  zInnerSplit = fZHalfLength + (fInnerRadius - maxEndInnerRad) * dz / dr;
437 
438  // Build associated vertices
439  v5 = transform.TransformPoint(
440  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
441  + zInnerSplit));
442  v6 = transform.TransformPoint(
443  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
444  - zInnerSplit));
445  }
446  else
447  {
448  v2 = transform.TransformPoint(
449  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
450  + fZHalfLength));
451  v3 = transform.TransformPoint(
452  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
453  - fZHalfLength));
454  }
455 
456  //
457  // decide vertices of +ve phi boundary
458  //
459 
460  cosPhi = std::cos(phimax);
461  sinPhi = std::sin(phimax);
462 
463  // Outer hyperbolic surface
464 
465  w0 = transform.TransformPoint(
466  G4ThreeVector(fudgeEndOuterRad*cosPhi, fudgeEndOuterRad*sinPhi,
467  + fZHalfLength));
468  w1 = transform.TransformPoint(
469  G4ThreeVector(fudgeEndOuterRad*cosPhi, fudgeEndOuterRad*sinPhi,
470  - fZHalfLength));
471  if (splitOuter)
472  {
473  G4double r = rFudge*fOuterRadius;
474 
475  w4 = transform.TransformPoint(G4ThreeVector( r*cosPhi, r*sinPhi, 0 ));
476 
477  AddPolyToExtent( v0, v4, w4, w0, voxelLimit, axis, extentList );
478  AddPolyToExtent( v4, v1, w1, w4, voxelLimit, axis, extentList );
479  }
480  else
481  {
482  AddPolyToExtent( v0, v1, w1, w0, voxelLimit, axis, extentList );
483  }
484 
485  // Inner hyperbolic surface
486 
487  if (splitInner)
488  {
489  w2 = transform.TransformPoint(
490  G4ThreeVector(maxEndInnerRad*cosPhi, maxEndInnerRad*sinPhi,
491  + fZHalfLength));
492  w3 = transform.TransformPoint(
493  G4ThreeVector(maxEndInnerRad*cosPhi, maxEndInnerRad*sinPhi,
494  - fZHalfLength));
495 
496  w5 = transform.TransformPoint(
497  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
498  + zInnerSplit));
499  w6 = transform.TransformPoint(
500  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
501  - zInnerSplit));
502 
503  AddPolyToExtent( v3, v6, w6, w3, voxelLimit, axis, extentList );
504  AddPolyToExtent( v6, v5, w5, w6, voxelLimit, axis, extentList );
505  AddPolyToExtent( v5, v2, w2, w5, voxelLimit, axis, extentList );
506 
507  }
508  else
509  {
510  w2 = transform.TransformPoint(
511  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
512  + fZHalfLength));
513  w3 = transform.TransformPoint(
514  G4ThreeVector(fInnerRadius*cosPhi, fInnerRadius*sinPhi,
515  - fZHalfLength));
516 
517  AddPolyToExtent( v3, v2, w2, w3, voxelLimit, axis, extentList );
518  }
519 
520  //
521  // Endplate segments
522  //
523  AddPolyToExtent( v1, v3, w3, w1, voxelLimit, axis, extentList );
524  AddPolyToExtent( v2, v0, w0, w2, voxelLimit, axis, extentList );
525 
526  //
527  // Return min/max value
528  //
529  return extentList.GetExtent( min, max );
530 }
531 
532 
533 //=====================================================================
534 //* AddPolyToExtent ---------------------------------------------------
535 
537  const G4ThreeVector &v1,
538  const G4ThreeVector &w1,
539  const G4ThreeVector &w0,
540  const G4VoxelLimits &voxelLimit,
541  const EAxis axis,
542  G4SolidExtentList &extentList )
543 {
544  // Utility function for CalculateExtent
545  //
546  G4ClippablePolygon phiPoly;
547 
548  phiPoly.AddVertexInOrder( v0 );
549  phiPoly.AddVertexInOrder( v1 );
550  phiPoly.AddVertexInOrder( w1 );
551  phiPoly.AddVertexInOrder( w0 );
552 
553  if (phiPoly.PartialClip( voxelLimit, axis ))
554  {
555  phiPoly.SetNormal( (v1-v0).cross(w0-v0).unit() );
556  extentList.AddSurface( phiPoly );
557  }
558 }
559 
560 
561 //=====================================================================
562 //* Inside ------------------------------------------------------------
563 
565 {
566 
567  const G4double halftol
569  // static G4int timerid = -1;
570  // G4Timer timer(timerid, "G4TwistedTubs", "Inside");
571  // timer.Start();
572 
573  G4ThreeVector *tmpp;
574  EInside *tmpinside;
575  if (fLastInside.p == p)
576  {
577  return fLastInside.inside;
578  }
579  else
580  {
581  tmpp = const_cast<G4ThreeVector*>(&(fLastInside.p));
582  tmpinside = const_cast<EInside*>(&(fLastInside.inside));
583  tmpp->set(p.x(), p.y(), p.z());
584  }
585 
586  EInside outerhypearea = ((G4TwistTubsHypeSide *)fOuterHype)->Inside(p);
587  G4double innerhyperho = ((G4TwistTubsHypeSide *)fInnerHype)->GetRhoAtPZ(p);
588  G4double distanceToOut = p.getRho() - innerhyperho; // +ve: inside
589 
590  if ((outerhypearea == kOutside) || (distanceToOut < -halftol))
591  {
592  *tmpinside = kOutside;
593  }
594  else if (outerhypearea == kSurface)
595  {
596  *tmpinside = kSurface;
597  }
598  else
599  {
600  if (distanceToOut <= halftol)
601  {
602  *tmpinside = kSurface;
603  }
604  else
605  {
606  *tmpinside = kInside;
607  }
608  }
609 
610  return fLastInside.inside;
611 }
612 
613 //=====================================================================
614 //* SurfaceNormal -----------------------------------------------------
615 
617 {
618  //
619  // return the normal unit vector to the Hyperbolical Surface at a point
620  // p on (or nearly on) the surface
621  //
622  // Which of the three or four surfaces are we closest to?
623  //
624 
625  if (fLastNormal.p == p)
626  {
627  return fLastNormal.vec;
628  }
629  G4ThreeVector *tmpp =
630  const_cast<G4ThreeVector*>(&(fLastNormal.p));
631  G4ThreeVector *tmpnormal =
632  const_cast<G4ThreeVector*>(&(fLastNormal.vec));
633  G4VTwistSurface **tmpsurface =
634  const_cast<G4VTwistSurface**>(fLastNormal.surface);
635  tmpp->set(p.x(), p.y(), p.z());
636 
637  G4double distance = kInfinity;
638 
639  G4VTwistSurface *surfaces[6];
640  surfaces[0] = fLatterTwisted;
641  surfaces[1] = fFormerTwisted;
642  surfaces[2] = fInnerHype;
643  surfaces[3] = fOuterHype;
644  surfaces[4] = fLowerEndcap;
645  surfaces[5] = fUpperEndcap;
646 
647  G4ThreeVector xx;
648  G4ThreeVector bestxx;
649  G4int i;
650  G4int besti = -1;
651  for (i=0; i< 6; i++)
652  {
653  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
654  if (tmpdistance < distance)
655  {
656  distance = tmpdistance;
657  bestxx = xx;
658  besti = i;
659  }
660  }
661 
662  tmpsurface[0] = surfaces[besti];
663  *tmpnormal = tmpsurface[0]->GetNormal(bestxx, true);
664 
665  return fLastNormal.vec;
666 }
667 
668 //=====================================================================
669 //* DistanceToIn (p, v) -----------------------------------------------
670 
672  const G4ThreeVector& v ) const
673 {
674 
675  // DistanceToIn (p, v):
676  // Calculate distance to surface of shape from `outside'
677  // along with the v, allowing for tolerance.
678  // The function returns kInfinity if no intersection or
679  // just grazing within tolerance.
680 
681  //
682  // checking last value
683  //
684 
685  G4ThreeVector *tmpp;
686  G4ThreeVector *tmpv;
687  G4double *tmpdist;
688  if ((fLastDistanceToInWithV.p == p) && (fLastDistanceToInWithV.vec == v))
689  {
690  return fLastDistanceToIn.value;
691  }
692  else
693  {
694  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.p));
695  tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.vec));
696  tmpdist = const_cast<G4double*>(&(fLastDistanceToInWithV.value));
697  tmpp->set(p.x(), p.y(), p.z());
698  tmpv->set(v.x(), v.y(), v.z());
699  }
700 
701  //
702  // Calculate DistanceToIn(p,v)
703  //
704 
705  EInside currentside = Inside(p);
706 
707  if (currentside == kInside)
708  {
709  }
710  else
711  {
712  if (currentside == kSurface)
713  {
714  // particle is just on a boundary.
715  // If the particle is entering to the volume, return 0.
716  //
718  if (normal*v < 0)
719  {
720  *tmpdist = 0;
722  }
723  }
724  }
725 
726  // now, we can take smallest positive distance.
727 
728  // Initialize
729  //
730  G4double distance = kInfinity;
731 
732  // find intersections and choose nearest one.
733  //
734  G4VTwistSurface *surfaces[6];
735  surfaces[0] = fLowerEndcap;
736  surfaces[1] = fUpperEndcap;
737  surfaces[2] = fLatterTwisted;
738  surfaces[3] = fFormerTwisted;
739  surfaces[4] = fInnerHype;
740  surfaces[5] = fOuterHype;
741 
742  G4ThreeVector xx;
743  G4ThreeVector bestxx;
744  G4int i;
745  for (i=0; i< 6; i++)
746  {
747  G4double tmpdistance = surfaces[i]->DistanceToIn(p, v, xx);
748  if (tmpdistance < distance)
749  {
750  distance = tmpdistance;
751  bestxx = xx;
752  }
753  }
754  *tmpdist = distance;
755 
757 }
758 
759 //=====================================================================
760 //* DistanceToIn (p) --------------------------------------------------
761 
763 {
764  // DistanceToIn(p):
765  // Calculate distance to surface of shape from `outside',
766  // allowing for tolerance
767 
768  //
769  // checking last value
770  //
771 
772  G4ThreeVector *tmpp;
773  G4double *tmpdist;
774  if (fLastDistanceToIn.p == p)
775  {
776  return fLastDistanceToIn.value;
777  }
778  else
779  {
780  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToIn.p));
781  tmpdist = const_cast<G4double*>(&(fLastDistanceToIn.value));
782  tmpp->set(p.x(), p.y(), p.z());
783  }
784 
785  //
786  // Calculate DistanceToIn(p)
787  //
788 
789  EInside currentside = Inside(p);
790 
791  switch (currentside)
792  {
793  case (kInside) :
794  {}
795  case (kSurface) :
796  {
797  *tmpdist = 0.;
798  return fLastDistanceToIn.value;
799  }
800  case (kOutside) :
801  {
802  // Initialize
803  G4double distance = kInfinity;
804 
805  // find intersections and choose nearest one.
806  G4VTwistSurface *surfaces[6];
807  surfaces[0] = fLowerEndcap;
808  surfaces[1] = fUpperEndcap;
809  surfaces[2] = fLatterTwisted;
810  surfaces[3] = fFormerTwisted;
811  surfaces[4] = fInnerHype;
812  surfaces[5] = fOuterHype;
813 
814  G4int i;
815  G4ThreeVector xx;
816  G4ThreeVector bestxx;
817  for (i=0; i< 6; i++)
818  {
819  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
820  if (tmpdistance < distance)
821  {
822  distance = tmpdistance;
823  bestxx = xx;
824  }
825  }
826  *tmpdist = distance;
827  return fLastDistanceToIn.value;
828  }
829  default :
830  {
831  G4Exception("G4TwistedTubs::DistanceToIn(p)", "GeomSolids0003",
832  FatalException, "Unknown point location!");
833  }
834  } // switch end
835 
836  return kInfinity;
837 }
838 
839 //=====================================================================
840 //* DistanceToOut (p, v) ----------------------------------------------
841 
843  const G4ThreeVector& v,
844  const G4bool calcNorm,
845  G4bool *validNorm,
846  G4ThreeVector *norm ) const
847 {
848  // DistanceToOut (p, v):
849  // Calculate distance to surface of shape from `inside'
850  // along with the v, allowing for tolerance.
851  // The function returns kInfinity if no intersection or
852  // just grazing within tolerance.
853 
854  //
855  // checking last value
856  //
857 
858  G4ThreeVector *tmpp;
859  G4ThreeVector *tmpv;
860  G4double *tmpdist;
862  {
864  }
865  else
866  {
867  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.p));
868  tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.vec));
869  tmpdist = const_cast<G4double*>(&(fLastDistanceToOutWithV.value));
870  tmpp->set(p.x(), p.y(), p.z());
871  tmpv->set(v.x(), v.y(), v.z());
872  }
873 
874  //
875  // Calculate DistanceToOut(p,v)
876  //
877 
878  EInside currentside = Inside(p);
879 
880  if (currentside == kOutside)
881  {
882  }
883  else
884  {
885  if (currentside == kSurface)
886  {
887  // particle is just on a boundary.
888  // If the particle is exiting from the volume, return 0.
889  //
891  G4VTwistSurface *blockedsurface = fLastNormal.surface[0];
892  if (normal*v > 0)
893  {
894  if (calcNorm)
895  {
896  *norm = (blockedsurface->GetNormal(p, true));
897  *validNorm = blockedsurface->IsValidNorm();
898  }
899  *tmpdist = 0.;
901  }
902  }
903  }
904 
905  // now, we can take smallest positive distance.
906 
907  // Initialize
908  //
909  G4double distance = kInfinity;
910 
911  // find intersections and choose nearest one.
912  //
913  G4VTwistSurface *surfaces[6];
914  surfaces[0] = fLatterTwisted;
915  surfaces[1] = fFormerTwisted;
916  surfaces[2] = fInnerHype;
917  surfaces[3] = fOuterHype;
918  surfaces[4] = fLowerEndcap;
919  surfaces[5] = fUpperEndcap;
920 
921  G4int i;
922  G4int besti = -1;
923  G4ThreeVector xx;
924  G4ThreeVector bestxx;
925  for (i=0; i< 6; i++)
926  {
927  G4double tmpdistance = surfaces[i]->DistanceToOut(p, v, xx);
928  if (tmpdistance < distance)
929  {
930  distance = tmpdistance;
931  bestxx = xx;
932  besti = i;
933  }
934  }
935 
936  if (calcNorm)
937  {
938  if (besti != -1)
939  {
940  *norm = (surfaces[besti]->GetNormal(p, true));
941  *validNorm = surfaces[besti]->IsValidNorm();
942  }
943  }
944 
945  *tmpdist = distance;
946 
948 }
949 
950 
951 //=====================================================================
952 //* DistanceToOut (p) ----------------------------------------------
953 
955 {
956  // DistanceToOut(p):
957  // Calculate distance to surface of shape from `inside',
958  // allowing for tolerance
959 
960  //
961  // checking last value
962  //
963 
964  G4ThreeVector *tmpp;
965  G4double *tmpdist;
966  if (fLastDistanceToOut.p == p)
967  {
968  return fLastDistanceToOut.value;
969  }
970  else
971  {
972  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOut.p));
973  tmpdist = const_cast<G4double*>(&(fLastDistanceToOut.value));
974  tmpp->set(p.x(), p.y(), p.z());
975  }
976 
977  //
978  // Calculate DistanceToOut(p)
979  //
980 
981  EInside currentside = Inside(p);
982 
983  switch (currentside)
984  {
985  case (kOutside) :
986  {
987  }
988  case (kSurface) :
989  {
990  *tmpdist = 0.;
991  return fLastDistanceToOut.value;
992  }
993  case (kInside) :
994  {
995  // Initialize
996  G4double distance = kInfinity;
997 
998  // find intersections and choose nearest one.
999  G4VTwistSurface *surfaces[6];
1000  surfaces[0] = fLatterTwisted;
1001  surfaces[1] = fFormerTwisted;
1002  surfaces[2] = fInnerHype;
1003  surfaces[3] = fOuterHype;
1004  surfaces[4] = fLowerEndcap;
1005  surfaces[5] = fUpperEndcap;
1006 
1007  G4int i;
1008  G4ThreeVector xx;
1009  G4ThreeVector bestxx;
1010  for (i=0; i< 6; i++)
1011  {
1012  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
1013  if (tmpdistance < distance)
1014  {
1015  distance = tmpdistance;
1016  bestxx = xx;
1017  }
1018  }
1019  *tmpdist = distance;
1020 
1021  return fLastDistanceToOut.value;
1022  }
1023  default :
1024  {
1025  G4Exception("G4TwistedTubs::DistanceToOut(p)", "GeomSolids0003",
1026  FatalException, "Unknown point location!");
1027  }
1028  } // switch end
1029 
1030  return 0;
1031 }
1032 
1033 //=====================================================================
1034 //* StreamInfo --------------------------------------------------------
1035 
1036 std::ostream& G4TwistedTubs::StreamInfo(std::ostream& os) const
1037 {
1038  //
1039  // Stream object contents to an output stream
1040  //
1041  G4int oldprc = os.precision(16);
1042  os << "-----------------------------------------------------------\n"
1043  << " *** Dump for solid - " << GetName() << " ***\n"
1044  << " ===================================================\n"
1045  << " Solid type: G4TwistedTubs\n"
1046  << " Parameters: \n"
1047  << " -ve end Z : " << fEndZ[0]/mm << " mm \n"
1048  << " +ve end Z : " << fEndZ[1]/mm << " mm \n"
1049  << " inner end radius(-ve z): " << fEndInnerRadius[0]/mm << " mm \n"
1050  << " inner end radius(+ve z): " << fEndInnerRadius[1]/mm << " mm \n"
1051  << " outer end radius(-ve z): " << fEndOuterRadius[0]/mm << " mm \n"
1052  << " outer end radius(+ve z): " << fEndOuterRadius[1]/mm << " mm \n"
1053  << " inner radius (z=0) : " << fInnerRadius/mm << " mm \n"
1054  << " outer radius (z=0) : " << fOuterRadius/mm << " mm \n"
1055  << " twisted angle : " << fPhiTwist/degree << " degrees \n"
1056  << " inner stereo angle : " << fInnerStereo/degree << " degrees \n"
1057  << " outer stereo angle : " << fOuterStereo/degree << " degrees \n"
1058  << " phi-width of a piece : " << fDPhi/degree << " degrees \n"
1059  << "-----------------------------------------------------------\n";
1060  os.precision(oldprc);
1061 
1062  return os;
1063 }
1064 
1065 
1066 //=====================================================================
1067 //* DiscribeYourselfTo ------------------------------------------------
1068 
1070 {
1071  scene.AddSolid (*this);
1072 }
1073 
1074 //=====================================================================
1075 //* GetExtent ---------------------------------------------------------
1076 
1078 {
1079  // Define the sides of the box into which the G4Tubs instance would fit.
1080 
1081  G4double maxEndOuterRad = (fEndOuterRadius[0] > fEndOuterRadius[1] ? 0 : 1);
1082  return G4VisExtent( -maxEndOuterRad, maxEndOuterRad,
1083  -maxEndOuterRad, maxEndOuterRad,
1085 }
1086 
1087 //=====================================================================
1088 //* CreatePolyhedron --------------------------------------------------
1089 
1091 {
1092  // number of meshes
1093  //
1095  const G4int k =
1096  G4int(G4Polyhedron::GetNumberOfRotationSteps() * dA / twopi) + 2;
1097  const G4int n =
1098  G4int(G4Polyhedron::GetNumberOfRotationSteps() * fPhiTwist / twopi) + 2;
1099 
1100  const G4int nnodes = 4*(k-1)*(n-2) + 2*k*k ;
1101  const G4int nfaces = 4*(k-1)*(n-1) + 2*(k-1)*(k-1) ;
1102 
1103  G4Polyhedron *ph=new G4Polyhedron;
1104  typedef G4double G4double3[3];
1105  typedef G4int G4int4[4];
1106  G4double3* xyz = new G4double3[nnodes]; // number of nodes
1107  G4int4* faces = new G4int4[nfaces] ; // number of faces
1108  fLowerEndcap->GetFacets(k,k,xyz,faces,0) ;
1109  fUpperEndcap->GetFacets(k,k,xyz,faces,1) ;
1110  fInnerHype->GetFacets(k,n,xyz,faces,2) ;
1111  fFormerTwisted->GetFacets(k,n,xyz,faces,3) ;
1112  fOuterHype->GetFacets(k,n,xyz,faces,4) ;
1113  fLatterTwisted->GetFacets(k,n,xyz,faces,5) ;
1114 
1115  ph->createPolyhedron(nnodes,nfaces,xyz,faces);
1116 
1117  delete[] xyz;
1118  delete[] faces;
1119 
1120  return ph;
1121 }
1122 
1123 //=====================================================================
1124 //* GetPolyhedron -----------------------------------------------------
1125 
1127 {
1128  if ((!fpPolyhedron) ||
1130  fpPolyhedron->GetNumberOfRotationSteps()))
1131  {
1132  delete fpPolyhedron;
1134  }
1135  return fpPolyhedron;
1136 }
1137 
1138 //=====================================================================
1139 //* CreateSurfaces ----------------------------------------------------
1140 
1142 {
1143  // create 6 surfaces of TwistedTub
1144 
1145  G4ThreeVector x0(0, 0, fEndZ[0]);
1146  G4ThreeVector n (0, 0, -1);
1147 
1148  fLowerEndcap = new G4TwistTubsFlatSide("LowerEndcap",
1150  fDPhi, fEndPhi, fEndZ, -1) ;
1151 
1152  fUpperEndcap = new G4TwistTubsFlatSide("UpperEndcap",
1154  fDPhi, fEndPhi, fEndZ, 1) ;
1155 
1156  G4RotationMatrix rotHalfDPhi;
1157  rotHalfDPhi.rotateZ(0.5*fDPhi);
1158 
1159  fLatterTwisted = new G4TwistTubsSide("LatterTwisted",
1161  fDPhi, fEndPhi, fEndZ,
1163  1 ) ;
1164  fFormerTwisted = new G4TwistTubsSide("FormerTwisted",
1166  fDPhi, fEndPhi, fEndZ,
1168  -1 ) ;
1169 
1170  fInnerHype = new G4TwistTubsHypeSide("InnerHype",
1172  fDPhi, fEndPhi, fEndZ,
1175  fOuterHype = new G4TwistTubsHypeSide("OuterHype",
1177  fDPhi, fEndPhi, fEndZ,
1180 
1181 
1182  // set neighbour surfaces
1183  //
1196 }
1197 
1198 
1199 //=====================================================================
1200 //* GetEntityType -----------------------------------------------------
1201 
1203 {
1204  return G4String("G4TwistedTubs");
1205 }
1206 
1207 //=====================================================================
1208 //* Clone -------------------------------------------------------------
1209 
1211 {
1212  return new G4TwistedTubs(*this);
1213 }
1214 
1215 //=====================================================================
1216 //* GetCubicVolume ----------------------------------------------------
1217 
1219 {
1220  if(fCubicVolume != 0.) {;}
1223  return fCubicVolume;
1224 }
1225 
1226 //=====================================================================
1227 //* GetSurfaceArea ----------------------------------------------------
1228 
1230 {
1231  if(fSurfaceArea != 0.) {;}
1233  return fSurfaceArea;
1234 }
1235 
1236 //=====================================================================
1237 //* GetPointOnSurface -------------------------------------------------
1238 
1240 {
1241 
1243  G4double phi , phimin, phimax ;
1244  G4double x , xmin, xmax ;
1245  G4double r , rmin, rmax ;
1246 
1253 
1254  G4double chose = G4RandFlat::shoot(0.,a1 + a2 + a3 + a4 + a5 + a6) ;
1255 
1256  if(chose < a1)
1257  {
1258 
1259  phimin = fOuterHype->GetBoundaryMin(z) ;
1260  phimax = fOuterHype->GetBoundaryMax(z) ;
1261  phi = G4RandFlat::shoot(phimin,phimax) ;
1262 
1263  return fOuterHype->SurfacePoint(phi,z,true) ;
1264 
1265  }
1266  else if ( (chose >= a1) && (chose < a1 + a2 ) )
1267  {
1268 
1269  phimin = fInnerHype->GetBoundaryMin(z) ;
1270  phimax = fInnerHype->GetBoundaryMax(z) ;
1271  phi = G4RandFlat::shoot(phimin,phimax) ;
1272 
1273  return fInnerHype->SurfacePoint(phi,z,true) ;
1274 
1275  }
1276  else if ( (chose >= a1 + a2 ) && (chose < a1 + a2 + a3 ) )
1277  {
1278 
1279  xmin = fLatterTwisted->GetBoundaryMin(z) ;
1280  xmax = fLatterTwisted->GetBoundaryMax(z) ;
1281  x = G4RandFlat::shoot(xmin,xmax) ;
1282 
1283  return fLatterTwisted->SurfacePoint(x,z,true) ;
1284 
1285  }
1286  else if ( (chose >= a1 + a2 + a3 ) && (chose < a1 + a2 + a3 + a4 ) )
1287  {
1288 
1289  xmin = fFormerTwisted->GetBoundaryMin(z) ;
1290  xmax = fFormerTwisted->GetBoundaryMax(z) ;
1291  x = G4RandFlat::shoot(xmin,xmax) ;
1292 
1293  return fFormerTwisted->SurfacePoint(x,z,true) ;
1294  }
1295  else if( (chose >= a1 + a2 + a3 + a4 )&&(chose < a1 + a2 + a3 + a4 + a5 ) )
1296  {
1297  rmin = GetEndInnerRadius(0) ;
1298  rmax = GetEndOuterRadius(0) ;
1299  r = std::sqrt(G4RandFlat::shoot()*(sqr(rmax)-sqr(rmin))+sqr(rmin));
1300 
1301  phimin = fLowerEndcap->GetBoundaryMin(r) ;
1302  phimax = fLowerEndcap->GetBoundaryMax(r) ;
1303  phi = G4RandFlat::shoot(phimin,phimax) ;
1304 
1305  return fLowerEndcap->SurfacePoint(phi,r,true) ;
1306  }
1307  else
1308  {
1309  rmin = GetEndInnerRadius(1) ;
1310  rmax = GetEndOuterRadius(1) ;
1311  r = rmin + (rmax-rmin)*std::sqrt(G4RandFlat::shoot());
1312 
1313  phimin = fUpperEndcap->GetBoundaryMin(r) ;
1314  phimax = fUpperEndcap->GetBoundaryMax(r) ;
1315  phi = G4RandFlat::shoot(phimin,phimax) ;
1316 
1317  return fUpperEndcap->SurfacePoint(phi,r,true) ;
1318  }
1319 }
G4double fZHalfLength
G4String GetName() const
ThreeVector shoot(const G4int Ap, const G4int Af)
G4double GetEndInnerRadius() const
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double fPhiTwist
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetPointOnSurface() const
CLHEP::HepRotation G4RotationMatrix
virtual G4ThreeVector GetNormal(const G4ThreeVector &xx, G4bool isGlobal)=0
G4double z
Definition: TRTMaterials.hh:39
static void AddPolyToExtent(const G4ThreeVector &v0, const G4ThreeVector &v1, const G4ThreeVector &w1, const G4ThreeVector &w0, const G4VoxelLimits &voxellimit, const EAxis axis, G4SolidExtentList &extentlist)
static const G4double a1
std::ostream & StreamInfo(std::ostream &os) const
G4bool GetExtent(G4double &min, G4double &max) const
const G4double pi
G4Polyhedron * CreatePolyhedron() const
G4VTwistSurface * fFormerTwisted
void SetNormal(const G4ThreeVector &newNormal)
G4double fEndOuterRadius[2]
LastValueWithDoubleVector fLastDistanceToInWithV
static const G4double a4
virtual G4bool PartialClip(const G4VoxelLimits &voxelLimit, const EAxis IgnoreMe)
virtual void AddVertexInOrder(const G4ThreeVector vertex)
LastValue fLastDistanceToIn
G4double a
Definition: TRTMaterials.hh:39
EInside Inside(const G4ThreeVector &p) const
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4double fEndInnerRadius[2]
G4bool IsValidNorm() const
G4GeometryType GetEntityType() const
G4VTwistSurface * fUpperEndcap
G4double fEndPhi[2]
virtual G4double DistanceToIn(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
G4double fInnerRadius
G4double fEndZ[2]
G4double GetEndOuterRadius() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4VTwistSurface * fOuterHype
G4double fSurfaceArea
virtual G4double GetSurfaceArea()=0
G4double fInnerRadius2
G4double fOuterRadius2
void SetNeighbours(G4VTwistSurface *axis0min, G4VTwistSurface *axis1min, G4VTwistSurface *axis0max, G4VTwistSurface *axis1max)
virtual G4double DistanceTo(const G4ThreeVector &gp, G4ThreeVector &gxx)
void DescribeYourselfTo(G4VGraphicsScene &scene) const
virtual G4double DistanceToOut(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
virtual G4ThreeVector SurfacePoint(G4double, G4double, G4bool isGlobal=false)=0
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
G4VTwistSurface ** surface
void SetFields(G4double phitwist, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz)
void AddSurface(const G4ClippablePolygon &surface)
LastValueWithDoubleVector fLastDistanceToOutWithV
G4double fEndZ2[2]
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=G4bool(false), G4bool *validnorm=0, G4ThreeVector *n=0) const
const G4int n
G4VSolid * Clone() const
virtual G4double GetBoundaryMax(G4double)=0
static const G4double a3
G4double GetSurfaceArea()
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
LastValue fLastDistanceToOut
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)=0
virtual G4double GetBoundaryMin(G4double)=0
LastState fLastInside
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4VTwistSurface * fLatterTwisted
G4double fOuterRadius
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual ~G4TwistedTubs()
T max(const T t1, const T t2)
brief Return the largest of the two arguments
EInside
Definition: geomdefs.hh:58
G4VTwistSurface * fInnerHype
EAxis
Definition: geomdefs.hh:54
#define DBL_MIN
Definition: templates.hh:75
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double fTanOuterStereo
static const double degree
Definition: G4SIunits.hh:125
G4double fTanInnerStereo2
G4Polyhedron * fpPolyhedron
G4TwistedTubs(const G4String &pname, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
#define G4endl
Definition: G4ios.hh:61
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
LastVector fLastNormal
G4double fTanOuterStereo2
G4double fTanInnerStereo
T sqr(const T &x)
Definition: templates.hh:145
static const G4double a5
double G4double
Definition: G4Types.hh:76
G4double fInnerStereo
G4double GetCubicVolume()
G4double fCubicVolume
G4VTwistSurface * fLowerEndcap
G4VisExtent GetExtent() const
static const double mm
Definition: G4SIunits.hh:102
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
G4double fOuterStereo
static const G4double a2
static G4GeometryTolerance * GetInstance()
G4Polyhedron * GetPolyhedron() const
G4bool CalculateExtent(const EAxis paxis, const G4VoxelLimits &pvoxellimit, const G4AffineTransform &ptransform, G4double &pmin, G4double &pmax) const
G4TwistedTubs & operator=(const G4TwistedTubs &rhs)