Geant4  10.02.p02
G4VTwistedFaceted.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 // $Id: G4VTwistedFaceted.cc 83572 2014-09-01 15:23:27Z gcosmo $
27 //
28 //
29 // --------------------------------------------------------------------
30 // GEANT 4 class source file
31 //
32 //
33 // G4VTwistedFaceted.cc
34 //
35 // Author:
36 //
37 // 04-Nov-2004 - O.Link (Oliver.Link@cern.ch)
38 //
39 // --------------------------------------------------------------------
40 
41 #include "G4VTwistedFaceted.hh"
42 
43 #include "G4PhysicalConstants.hh"
44 #include "G4SystemOfUnits.hh"
45 #include "G4VoxelLimits.hh"
46 #include "G4AffineTransform.hh"
47 #include "G4SolidExtentList.hh"
48 #include "G4ClippablePolygon.hh"
49 #include "G4VPVParameterisation.hh"
50 #include "G4GeometryTolerance.hh"
51 #include "meshdefs.hh"
52 
53 #include "G4VGraphicsScene.hh"
54 #include "G4Polyhedron.hh"
55 #include "G4VisExtent.hh"
56 
57 #include "Randomize.hh"
58 
59 #include "G4AutoLock.hh"
60 
61 namespace
62 {
63  G4Mutex polyhedronMutex = G4MUTEX_INITIALIZER;
64 }
65 
66 //=====================================================================
67 //* constructors ------------------------------------------------------
68 
70 G4VTwistedFaceted( const G4String &pname, // Name of instance
71  G4double PhiTwist, // twist angle
72  G4double pDz, // half z length
73  G4double pTheta, // direction between end planes
74  G4double pPhi, // defined by polar and azim. angles
75  G4double pDy1, // half y length at -pDz
76  G4double pDx1, // half x length at -pDz,-pDy
77  G4double pDx2, // half x length at -pDz,+pDy
78  G4double pDy2, // half y length at +pDz
79  G4double pDx3, // half x length at +pDz,-pDy
80  G4double pDx4, // half x length at +pDz,+pDy
81  G4double pAlph // tilt angle
82  )
83  : G4VSolid(pname), fRebuildPolyhedron(false), fpPolyhedron(0),
84  fLowerEndcap(0), fUpperEndcap(0), fSide0(0),
85  fSide90(0), fSide180(0), fSide270(0),
86  fSurfaceArea(0.)
87 {
88 
89  G4double pDytmp ;
90  G4double fDxUp ;
91  G4double fDxDown ;
92 
93  fDx1 = pDx1 ;
94  fDx2 = pDx2 ;
95  fDx3 = pDx3 ;
96  fDx4 = pDx4 ;
97  fDy1 = pDy1 ;
98  fDy2 = pDy2 ;
99  fDz = pDz ;
100 
101  G4double kAngTolerance
103 
104  // maximum values
105  //
106  fDxDown = ( fDx1 > fDx2 ? fDx1 : fDx2 ) ;
107  fDxUp = ( fDx3 > fDx4 ? fDx3 : fDx4 ) ;
108  fDx = ( fDxUp > fDxDown ? fDxUp : fDxDown ) ;
109  fDy = ( fDy1 > fDy2 ? fDy1 : fDy2 ) ;
110 
111  // planarity check
112  //
113  if ( fDx1 != fDx2 && fDx3 != fDx4 )
114  {
115  pDytmp = fDy1 * ( fDx3 - fDx4 ) / ( fDx1 - fDx2 ) ;
116  if ( std::fabs(pDytmp - fDy2) > kCarTolerance )
117  {
118  std::ostringstream message;
119  message << "Not planar surface in untwisted Trapezoid: "
120  << GetName() << G4endl
121  << "fDy2 is " << fDy2 << " but should be "
122  << pDytmp << ".";
123  G4Exception("G4VTwistedFaceted::G4VTwistedFaceted()", "GeomSolids0002",
124  FatalErrorInArgument, message);
125  }
126  }
127 
128 #ifdef G4TWISTDEBUG
129  if ( fDx1 == fDx2 && fDx3 == fDx4 )
130  {
131  G4cout << "Trapezoid is a box" << G4endl ;
132  }
133 
134 #endif
135 
136  if ( ( fDx1 == fDx2 && fDx3 != fDx4 ) || ( fDx1 != fDx2 && fDx3 == fDx4 ) )
137  {
138  std::ostringstream message;
139  message << "Not planar surface in untwisted Trapezoid: "
140  << GetName() << G4endl
141  << "One endcap is rectangular, the other is a trapezoid." << G4endl
142  << "For planarity reasons they have to be rectangles or trapezoids "
143  << "on both sides.";
144  G4Exception("G4VTwistedFaceted::G4VTwistedFaceted()", "GeomSolids0002",
145  FatalErrorInArgument, message);
146  }
147 
148  // twist angle
149  //
150  fPhiTwist = PhiTwist ;
151 
152  // tilt angle
153  //
154  fAlph = pAlph ;
155  fTAlph = std::tan(fAlph) ;
156 
157  fTheta = pTheta ;
158  fPhi = pPhi ;
159 
160  // dx in surface equation
161  //
162  fdeltaX = 2 * fDz * std::tan(fTheta) * std::cos(fPhi) ;
163 
164  // dy in surface equation
165  //
166  fdeltaY = 2 * fDz * std::tan(fTheta) * std::sin(fPhi) ;
167 
168  if ( ! ( ( fDx1 > 2*kCarTolerance)
169  && ( fDx2 > 2*kCarTolerance)
170  && ( fDx3 > 2*kCarTolerance)
171  && ( fDx4 > 2*kCarTolerance)
172  && ( fDy1 > 2*kCarTolerance)
173  && ( fDy2 > 2*kCarTolerance)
174  && ( fDz > 2*kCarTolerance)
175  && ( std::fabs(fPhiTwist) > 2*kAngTolerance )
176  && ( std::fabs(fPhiTwist) < pi/2 )
177  && ( std::fabs(fAlph) < pi/2 )
178  && ( fTheta < pi/2 && fTheta >= 0 ) )
179  )
180  {
181  std::ostringstream message;
182  message << "Invalid dimensions. Too small, or twist angle too big: "
183  << GetName() << G4endl
184  << "fDx 1-4 = " << fDx1/cm << ", " << fDx2/cm << ", "
185  << fDx3/cm << ", " << fDx4/cm << " cm" << G4endl
186  << "fDy 1-2 = " << fDy1/cm << ", " << fDy2/cm << ", "
187  << " cm" << G4endl
188  << "fDz = " << fDz/cm << " cm" << G4endl
189  << " twistangle " << fPhiTwist/deg << " deg" << G4endl
190  << " phi,theta = " << fPhi/deg << ", " << fTheta/deg << " deg";
191  G4Exception("G4TwistedTrap::G4VTwistedFaceted()",
192  "GeomSolids0002", FatalErrorInArgument, message);
193  }
194  CreateSurfaces();
195  fCubicVolume = 2 * fDz * ( ( fDx1 + fDx2 ) * fDy1 + ( fDx3 + fDx4 ) * fDy2 );
196 }
197 
198 
199 //=====================================================================
200 //* Fake default constructor ------------------------------------------
201 
203  : G4VSolid(a), fRebuildPolyhedron(false), fpPolyhedron(0),
204  fTheta(0.), fPhi(0.), fDy1(0.),
205  fDx1(0.), fDx2(0.), fDy2(0.), fDx3(0.), fDx4(0.),
206  fDz(0.), fDx(0.), fDy(0.), fAlph(0.),
207  fTAlph(0.), fdeltaX(0.), fdeltaY(0.), fPhiTwist(0.),
208  fLowerEndcap(0), fUpperEndcap(0), fSide0(0), fSide90(0), fSide180(0),
209  fSide270(0), fCubicVolume(0.), fSurfaceArea(0.)
210 {
211 }
212 
213 //=====================================================================
214 //* destructor --------------------------------------------------------
215 
217 {
218  if (fLowerEndcap) { delete fLowerEndcap ; }
219  if (fUpperEndcap) { delete fUpperEndcap ; }
220 
221  if (fSide0) { delete fSide0 ; }
222  if (fSide90) { delete fSide90 ; }
223  if (fSide180) { delete fSide180 ; }
224  if (fSide270) { delete fSide270 ; }
225  if (fpPolyhedron) { delete fpPolyhedron; fpPolyhedron = 0; }
226 }
227 
228 //=====================================================================
229 //* Copy constructor --------------------------------------------------
230 
232  : G4VSolid(rhs), fRebuildPolyhedron(false), fpPolyhedron(0),
233  fTheta(rhs.fTheta), fPhi(rhs.fPhi),
234  fDy1(rhs.fDy1), fDx1(rhs.fDx1), fDx2(rhs.fDx2), fDy2(rhs.fDy2),
235  fDx3(rhs.fDx3), fDx4(rhs.fDx4), fDz(rhs.fDz), fDx(rhs.fDx), fDy(rhs.fDy),
236  fAlph(rhs.fAlph), fTAlph(rhs.fTAlph), fdeltaX(rhs.fdeltaX),
237  fdeltaY(rhs.fdeltaY), fPhiTwist(rhs.fPhiTwist), fLowerEndcap(0),
238  fUpperEndcap(0), fSide0(0), fSide90(0), fSide180(0), fSide270(0),
239  fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
240  fLastInside(rhs.fLastInside), fLastNormal(rhs.fLastNormal),
241  fLastDistanceToIn(rhs.fLastDistanceToIn),
242  fLastDistanceToOut(rhs.fLastDistanceToOut),
243  fLastDistanceToInWithV(rhs.fLastDistanceToInWithV),
244  fLastDistanceToOutWithV(rhs.fLastDistanceToOutWithV)
245 {
246  CreateSurfaces();
247 }
248 
249 
250 //=====================================================================
251 //* Assignment operator -----------------------------------------------
252 
254 {
255  // Check assignment to self
256  //
257  if (this == &rhs) { return *this; }
258 
259  // Copy base class data
260  //
261  G4VSolid::operator=(rhs);
262 
263  // Copy data
264  //
265  fTheta = rhs.fTheta; fPhi = rhs.fPhi;
266  fDy1= rhs.fDy1; fDx1= rhs.fDx1; fDx2= rhs.fDx2; fDy2= rhs.fDy2;
267  fDx3= rhs.fDx3; fDx4= rhs.fDx4; fDz= rhs.fDz; fDx= rhs.fDx; fDy= rhs.fDy;
268  fAlph= rhs.fAlph; fTAlph= rhs.fTAlph; fdeltaX= rhs.fdeltaX;
270  fUpperEndcap= 0; fSide0= 0; fSide90= 0; fSide180= 0; fSide270= 0;
272  fRebuildPolyhedron = false;
273  delete fpPolyhedron; fpPolyhedron= 0;
279 
280  CreateSurfaces();
281 
282  return *this;
283 }
284 
285 //=====================================================================
286 //* ComputeDimensions -------------------------------------------------
287 
289  const G4int ,
290  const G4VPhysicalVolume* )
291 {
292  G4Exception("G4VTwistedFaceted::ComputeDimensions()",
293  "GeomSolids0001", FatalException,
294  "G4VTwistedFaceted does not support Parameterisation.");
295 }
296 
297 //=====================================================================
298 //* CalculateExtent ---------------------------------------------------
299 
300 G4bool
302  const G4VoxelLimits &pVoxelLimit,
303  const G4AffineTransform &pTransform,
304  G4double &pMin,
305  G4double &pMax ) const
306 {
307  G4double maxRad = std::sqrt( fDx*fDx + fDy*fDy);
308 
309  if (!pTransform.IsRotated())
310  {
311  // Special case handling for unrotated boxes
312  // Compute x/y/z mins and maxs respecting limits, with early returns
313  // if outside limits. Then switch() on pAxis
314 
315  G4double xoffset,xMin,xMax;
316  G4double yoffset,yMin,yMax;
317  G4double zoffset,zMin,zMax;
318 
319  xoffset = pTransform.NetTranslation().x() ;
320  xMin = xoffset - maxRad ;
321  xMax = xoffset + maxRad ;
322 
323  if (pVoxelLimit.IsXLimited())
324  {
325  if ( xMin > pVoxelLimit.GetMaxXExtent()+kCarTolerance ||
326  xMax < pVoxelLimit.GetMinXExtent()-kCarTolerance ) return false;
327  else
328  {
329  if (xMin < pVoxelLimit.GetMinXExtent())
330  {
331  xMin = pVoxelLimit.GetMinXExtent() ;
332  }
333  if (xMax > pVoxelLimit.GetMaxXExtent())
334  {
335  xMax = pVoxelLimit.GetMaxXExtent() ;
336  }
337  }
338  }
339  yoffset = pTransform.NetTranslation().y() ;
340  yMin = yoffset - maxRad ;
341  yMax = yoffset + maxRad ;
342 
343  if (pVoxelLimit.IsYLimited())
344  {
345  if ( yMin > pVoxelLimit.GetMaxYExtent()+kCarTolerance ||
346  yMax < pVoxelLimit.GetMinYExtent()-kCarTolerance ) return false;
347  else
348  {
349  if (yMin < pVoxelLimit.GetMinYExtent())
350  {
351  yMin = pVoxelLimit.GetMinYExtent() ;
352  }
353  if (yMax > pVoxelLimit.GetMaxYExtent())
354  {
355  yMax = pVoxelLimit.GetMaxYExtent() ;
356  }
357  }
358  }
359  zoffset = pTransform.NetTranslation().z() ;
360  zMin = zoffset - fDz ;
361  zMax = zoffset + fDz ;
362 
363  if (pVoxelLimit.IsZLimited())
364  {
365  if ( zMin > pVoxelLimit.GetMaxZExtent()+kCarTolerance ||
366  zMax < pVoxelLimit.GetMinZExtent()-kCarTolerance ) return false;
367  else
368  {
369  if (zMin < pVoxelLimit.GetMinZExtent())
370  {
371  zMin = pVoxelLimit.GetMinZExtent() ;
372  }
373  if (zMax > pVoxelLimit.GetMaxZExtent())
374  {
375  zMax = pVoxelLimit.GetMaxZExtent() ;
376  }
377  }
378  }
379  switch (pAxis)
380  {
381  case kXAxis:
382  pMin = xMin ;
383  pMax = xMax ;
384  break ;
385  case kYAxis:
386  pMin=yMin;
387  pMax=yMax;
388  break;
389  case kZAxis:
390  pMin=zMin;
391  pMax=zMax;
392  break;
393  default:
394  break;
395  }
396  pMin -= kCarTolerance ;
397  pMax += kCarTolerance ;
398 
399  return true;
400  }
401  else // General rotated case - create and clip mesh to boundaries
402  {
403  G4bool existsAfterClip = false ;
404  G4ThreeVectorList* vertices ;
405 
406  pMin = +kInfinity ;
407  pMax = -kInfinity ;
408 
409  // Calculate rotated vertex coordinates
410 
411  vertices = CreateRotatedVertices(pTransform) ;
412  ClipCrossSection(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
413  ClipCrossSection(vertices,4,pVoxelLimit,pAxis,pMin,pMax) ;
414  ClipBetweenSections(vertices,0,pVoxelLimit,pAxis,pMin,pMax) ;
415 
416  if (pVoxelLimit.IsLimited(pAxis) == false)
417  {
418  if ( pMin != kInfinity || pMax != -kInfinity )
419  {
420  existsAfterClip = true ;
421 
422  // Add 2*tolerance to avoid precision troubles
423 
424  pMin -= kCarTolerance;
425  pMax += kCarTolerance;
426  }
427  }
428  else
429  {
430  G4ThreeVector clipCentre(
431  ( pVoxelLimit.GetMinXExtent()+pVoxelLimit.GetMaxXExtent())*0.5,
432  ( pVoxelLimit.GetMinYExtent()+pVoxelLimit.GetMaxYExtent())*0.5,
433  ( pVoxelLimit.GetMinZExtent()+pVoxelLimit.GetMaxZExtent())*0.5 );
434 
435  if ( pMin != kInfinity || pMax != -kInfinity )
436  {
437  existsAfterClip = true ;
438 
439  // Check to see if endpoints are in the solid
440 
441  clipCentre(pAxis) = pVoxelLimit.GetMinExtent(pAxis);
442 
443  if (Inside(pTransform.Inverse().TransformPoint(clipCentre))
444  != kOutside)
445  {
446  pMin = pVoxelLimit.GetMinExtent(pAxis);
447  }
448  else
449  {
450  pMin -= kCarTolerance;
451  }
452  clipCentre(pAxis) = pVoxelLimit.GetMaxExtent(pAxis);
453 
454  if (Inside(pTransform.Inverse().TransformPoint(clipCentre))
455  != kOutside)
456  {
457  pMax = pVoxelLimit.GetMaxExtent(pAxis);
458  }
459  else
460  {
461  pMax += kCarTolerance;
462  }
463  }
464 
465  // Check for case where completely enveloping clipping volume
466  // If point inside then we are confident that the solid completely
467  // envelopes the clipping volume. Hence set min/max extents according
468  // to clipping volume extents along the specified axis.
469 
470  else if (Inside(pTransform.Inverse().TransformPoint(clipCentre))
471  != kOutside)
472  {
473  existsAfterClip = true ;
474  pMin = pVoxelLimit.GetMinExtent(pAxis) ;
475  pMax = pVoxelLimit.GetMaxExtent(pAxis) ;
476  }
477  }
478  delete vertices;
479  return existsAfterClip;
480  }
481 
482 
483 }
484 
486 CreateRotatedVertices(const G4AffineTransform& pTransform) const
487 {
488 
489  G4ThreeVectorList* vertices = new G4ThreeVectorList();
490 
491  if (vertices)
492  {
493  vertices->reserve(8);
494 
495  G4double maxRad = std::sqrt( fDx*fDx + fDy*fDy);
496 
497  G4ThreeVector vertex0(-maxRad,-maxRad,-fDz) ;
498  G4ThreeVector vertex1(maxRad,-maxRad,-fDz) ;
499  G4ThreeVector vertex2(maxRad,maxRad,-fDz) ;
500  G4ThreeVector vertex3(-maxRad,maxRad,-fDz) ;
501  G4ThreeVector vertex4(-maxRad,-maxRad,fDz) ;
502  G4ThreeVector vertex5(maxRad,-maxRad,fDz) ;
503  G4ThreeVector vertex6(maxRad,maxRad,fDz) ;
504  G4ThreeVector vertex7(-maxRad,maxRad,fDz) ;
505 
506  vertices->push_back(pTransform.TransformPoint(vertex0));
507  vertices->push_back(pTransform.TransformPoint(vertex1));
508  vertices->push_back(pTransform.TransformPoint(vertex2));
509  vertices->push_back(pTransform.TransformPoint(vertex3));
510  vertices->push_back(pTransform.TransformPoint(vertex4));
511  vertices->push_back(pTransform.TransformPoint(vertex5));
512  vertices->push_back(pTransform.TransformPoint(vertex6));
513  vertices->push_back(pTransform.TransformPoint(vertex7));
514  }
515  else
516  {
517  DumpInfo();
518  G4Exception("G4VTwistedFaceted::CreateRotatedVertices()",
519  "GeomSolids0003", FatalException,
520  "Error in allocation of vertices. Out of memory !");
521  }
522  return vertices;
523 }
524 
525 //=====================================================================
526 //* Inside ------------------------------------------------------------
527 
529 {
530 
531  G4ThreeVector *tmpp;
532  EInside *tmpin;
533  if (fLastInside.p == p) {
534  return fLastInside.inside;
535  } else {
536  tmpp = const_cast<G4ThreeVector*>(&(fLastInside.p));
537  tmpin = const_cast<EInside*>(&(fLastInside.inside));
538  tmpp->set(p.x(), p.y(), p.z());
539  }
540 
541  *tmpin = kOutside ;
542 
543  G4double phi = p.z()/(2*fDz) * fPhiTwist ; // rotate the point to z=0
544  G4double cphi = std::cos(-phi) ;
545  G4double sphi = std::sin(-phi) ;
546 
547  G4double px = p.x() + fdeltaX * ( -phi/fPhiTwist) ; // shift
548  G4double py = p.y() + fdeltaY * ( -phi/fPhiTwist) ;
549  G4double pz = p.z() ;
550 
551  G4double posx = px * cphi - py * sphi ; // rotation
552  G4double posy = px * sphi + py * cphi ;
553  G4double posz = pz ;
554 
555  G4double xMin = Xcoef(posy,phi,fTAlph) - 2*Xcoef(posy,phi,0.) ;
556  G4double xMax = Xcoef(posy,phi,fTAlph) ;
557 
558  G4double yMax = GetValueB(phi)/2. ; // b(phi)/2 is limit
559  G4double yMin = -yMax ;
560 
561 #ifdef G4TWISTDEBUG
562 
563  G4cout << "inside called: p = " << p << G4endl ;
564  G4cout << "fDx1 = " << fDx1 << G4endl ;
565  G4cout << "fDx2 = " << fDx2 << G4endl ;
566  G4cout << "fDx3 = " << fDx3 << G4endl ;
567  G4cout << "fDx4 = " << fDx4 << G4endl ;
568 
569  G4cout << "fDy1 = " << fDy1 << G4endl ;
570  G4cout << "fDy2 = " << fDy2 << G4endl ;
571 
572  G4cout << "fDz = " << fDz << G4endl ;
573 
574  G4cout << "Tilt angle alpha = " << fAlph << G4endl ;
575  G4cout << "phi,theta = " << fPhi << " , " << fTheta << G4endl ;
576 
577  G4cout << "Twist angle = " << fPhiTwist << G4endl ;
578 
579  G4cout << "posx = " << posx << G4endl ;
580  G4cout << "posy = " << posy << G4endl ;
581  G4cout << "xMin = " << xMin << G4endl ;
582  G4cout << "xMax = " << xMax << G4endl ;
583  G4cout << "yMin = " << yMin << G4endl ;
584  G4cout << "yMax = " << yMax << G4endl ;
585 
586 #endif
587 
588 
589  if ( posx <= xMax - kCarTolerance*0.5
590  && posx >= xMin + kCarTolerance*0.5 )
591  {
592  if ( posy <= yMax - kCarTolerance*0.5
593  && posy >= yMin + kCarTolerance*0.5 )
594  {
595  if (std::fabs(posz) <= fDz - kCarTolerance*0.5 ) *tmpin = kInside ;
596  else if (std::fabs(posz) <= fDz + kCarTolerance*0.5 ) *tmpin = kSurface ;
597  }
598  else if ( posy <= yMax + kCarTolerance*0.5
599  && posy >= yMin - kCarTolerance*0.5 )
600  {
601  if (std::fabs(posz) <= fDz + kCarTolerance*0.5 ) *tmpin = kSurface ;
602  }
603  }
604  else if ( posx <= xMax + kCarTolerance*0.5
605  && posx >= xMin - kCarTolerance*0.5 )
606  {
607  if ( posy <= yMax + kCarTolerance*0.5
608  && posy >= yMin - kCarTolerance*0.5 )
609  {
610  if (std::fabs(posz) <= fDz + kCarTolerance*0.5) *tmpin = kSurface ;
611  }
612  }
613 
614 #ifdef G4TWISTDEBUG
615  G4cout << "inside = " << fLastInside.inside << G4endl ;
616 #endif
617 
618  return fLastInside.inside;
619 
620 }
621 
622 //=====================================================================
623 //* SurfaceNormal -----------------------------------------------------
624 
626 {
627  //
628  // return the normal unit vector to the Hyperbolical Surface at a point
629  // p on (or nearly on) the surface
630  //
631  // Which of the three or four surfaces are we closest to?
632  //
633 
634  if (fLastNormal.p == p)
635  {
636  return fLastNormal.vec;
637  }
638 
639  G4ThreeVector *tmpp = const_cast<G4ThreeVector*>(&(fLastNormal.p));
640  G4ThreeVector *tmpnormal = const_cast<G4ThreeVector*>(&(fLastNormal.vec));
641  G4VTwistSurface **tmpsurface = const_cast<G4VTwistSurface**>(fLastNormal.surface);
642  tmpp->set(p.x(), p.y(), p.z());
643 
644  G4double distance = kInfinity;
645 
646  G4VTwistSurface *surfaces[6];
647 
648  surfaces[0] = fSide0 ;
649  surfaces[1] = fSide90 ;
650  surfaces[2] = fSide180 ;
651  surfaces[3] = fSide270 ;
652  surfaces[4] = fLowerEndcap;
653  surfaces[5] = fUpperEndcap;
654 
655  G4ThreeVector xx;
656  G4ThreeVector bestxx;
657  G4int i;
658  G4int besti = -1;
659  for (i=0; i< 6; i++)
660  {
661  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
662  if (tmpdistance < distance)
663  {
664  distance = tmpdistance;
665  bestxx = xx;
666  besti = i;
667  }
668  }
669 
670  tmpsurface[0] = surfaces[besti];
671  *tmpnormal = tmpsurface[0]->GetNormal(bestxx, true);
672 
673  return fLastNormal.vec;
674 }
675 
676 //=====================================================================
677 //* DistanceToIn (p, v) -----------------------------------------------
678 
680  const G4ThreeVector& v ) const
681 {
682 
683  // DistanceToIn (p, v):
684  // Calculate distance to surface of shape from `outside'
685  // along with the v, allowing for tolerance.
686  // The function returns kInfinity if no intersection or
687  // just grazing within tolerance.
688 
689  //
690  // checking last value
691  //
692 
693  G4ThreeVector *tmpp;
694  G4ThreeVector *tmpv;
695  G4double *tmpdist;
697  {
698  return fLastDistanceToIn.value;
699  }
700  else
701  {
702  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.p));
703  tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.vec));
704  tmpdist = const_cast<G4double*>(&(fLastDistanceToInWithV.value));
705  tmpp->set(p.x(), p.y(), p.z());
706  tmpv->set(v.x(), v.y(), v.z());
707  }
708 
709  //
710  // Calculate DistanceToIn(p,v)
711  //
712 
713  EInside currentside = Inside(p);
714 
715  if (currentside == kInside)
716  {
717  }
718  else if (currentside == kSurface)
719  {
720  // particle is just on a boundary.
721  // if the particle is entering to the volume, return 0
722  //
724  if (normal*v < 0)
725  {
726  *tmpdist = 0;
728  }
729  }
730 
731  // now, we can take smallest positive distance.
732 
733  // Initialize
734  //
735  G4double distance = kInfinity;
736 
737  // Find intersections and choose nearest one
738  //
739  G4VTwistSurface *surfaces[6];
740 
741  surfaces[0] = fSide0;
742  surfaces[1] = fSide90 ;
743  surfaces[2] = fSide180 ;
744  surfaces[3] = fSide270 ;
745  surfaces[4] = fLowerEndcap;
746  surfaces[5] = fUpperEndcap;
747 
748  G4ThreeVector xx;
749  G4ThreeVector bestxx;
750  G4int i;
751  for (i=0; i < 6 ; i++)
752  {
753 #ifdef G4TWISTDEBUG
754  G4cout << G4endl << "surface " << i << ": " << G4endl << G4endl ;
755 #endif
756  G4double tmpdistance = surfaces[i]->DistanceToIn(p, v, xx);
757 #ifdef G4TWISTDEBUG
758  G4cout << "Solid DistanceToIn : distance = " << tmpdistance << G4endl ;
759  G4cout << "intersection point = " << xx << G4endl ;
760 #endif
761  if (tmpdistance < distance)
762  {
763  distance = tmpdistance;
764  bestxx = xx;
765  }
766  }
767 
768 #ifdef G4TWISTDEBUG
769  G4cout << "best distance = " << distance << G4endl ;
770 #endif
771 
772  *tmpdist = distance;
773  // timer.Stop();
775 }
776 
777 
779 {
780  // DistanceToIn(p):
781  // Calculate distance to surface of shape from `outside',
782  // allowing for tolerance
783  //
784 
785  //
786  // checking last value
787  //
788 
789  G4ThreeVector *tmpp;
790  G4double *tmpdist;
791  if (fLastDistanceToIn.p == p)
792  {
793  return fLastDistanceToIn.value;
794  }
795  else
796  {
797  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToIn.p));
798  tmpdist = const_cast<G4double*>(&(fLastDistanceToIn.value));
799  tmpp->set(p.x(), p.y(), p.z());
800  }
801 
802  //
803  // Calculate DistanceToIn(p)
804  //
805 
806  EInside currentside = Inside(p);
807 
808  switch (currentside)
809  {
810  case (kInside) :
811  {
812  }
813 
814  case (kSurface) :
815  {
816  *tmpdist = 0.;
817  return fLastDistanceToIn.value;
818  }
819 
820  case (kOutside) :
821  {
822  // Initialize
823  //
824  G4double distance = kInfinity;
825 
826  // Find intersections and choose nearest one
827  //
828  G4VTwistSurface *surfaces[6];
829 
830  surfaces[0] = fSide0;
831  surfaces[1] = fSide90 ;
832  surfaces[2] = fSide180 ;
833  surfaces[3] = fSide270 ;
834  surfaces[4] = fLowerEndcap;
835  surfaces[5] = fUpperEndcap;
836 
837  G4int i;
838  G4ThreeVector xx;
839  G4ThreeVector bestxx;
840  for (i=0; i< 6; i++)
841  {
842  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
843  if (tmpdistance < distance)
844  {
845  distance = tmpdistance;
846  bestxx = xx;
847  }
848  }
849  *tmpdist = distance;
850  return fLastDistanceToIn.value;
851  }
852 
853  default :
854  {
855  G4Exception("G4VTwistedFaceted::DistanceToIn(p)", "GeomSolids0003",
856  FatalException, "Unknown point location!");
857  }
858  } // switch end
859 
860  return 0;
861 }
862 
863 
864 //=====================================================================
865 //* DistanceToOut (p, v) ----------------------------------------------
866 
867 G4double
869  const G4ThreeVector& v,
870  const G4bool calcNorm,
871  G4bool *validNorm,
872  G4ThreeVector *norm ) const
873 {
874  // DistanceToOut (p, v):
875  // Calculate distance to surface of shape from `inside'
876  // along with the v, allowing for tolerance.
877  // The function returns kInfinity if no intersection or
878  // just grazing within tolerance.
879 
880  //
881  // checking last value
882  //
883 
884  G4ThreeVector *tmpp;
885  G4ThreeVector *tmpv;
886  G4double *tmpdist;
888  {
890  }
891  else
892  {
893  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.p));
894  tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.vec));
895  tmpdist = const_cast<G4double*>(&(fLastDistanceToOutWithV.value));
896  tmpp->set(p.x(), p.y(), p.z());
897  tmpv->set(v.x(), v.y(), v.z());
898  }
899 
900  //
901  // Calculate DistanceToOut(p,v)
902  //
903 
904  EInside currentside = Inside(p);
905 
906  if (currentside == kOutside)
907  {
908  }
909  else if (currentside == kSurface)
910  {
911  // particle is just on a boundary.
912  // if the particle is exiting from the volume, return 0
913  //
915  G4VTwistSurface *blockedsurface = fLastNormal.surface[0];
916  if (normal*v > 0)
917  {
918  if (calcNorm)
919  {
920  *norm = (blockedsurface->GetNormal(p, true));
921  *validNorm = blockedsurface->IsValidNorm();
922  }
923  *tmpdist = 0.;
924  // timer.Stop();
926  }
927  }
928 
929  // now, we can take smallest positive distance.
930 
931  // Initialize
932  G4double distance = kInfinity;
933 
934  // find intersections and choose nearest one.
935  G4VTwistSurface *surfaces[6];
936 
937  surfaces[0] = fSide0;
938  surfaces[1] = fSide90 ;
939  surfaces[2] = fSide180 ;
940  surfaces[3] = fSide270 ;
941  surfaces[4] = fLowerEndcap;
942  surfaces[5] = fUpperEndcap;
943 
944  G4int i;
945  G4int besti = -1;
946  G4ThreeVector xx;
947  G4ThreeVector bestxx;
948  for (i=0; i< 6 ; i++) {
949  G4double tmpdistance = surfaces[i]->DistanceToOut(p, v, xx);
950  if (tmpdistance < distance)
951  {
952  distance = tmpdistance;
953  bestxx = xx;
954  besti = i;
955  }
956  }
957 
958  if (calcNorm)
959  {
960  if (besti != -1)
961  {
962  *norm = (surfaces[besti]->GetNormal(p, true));
963  *validNorm = surfaces[besti]->IsValidNorm();
964  }
965  }
966 
967  *tmpdist = distance;
968  // timer.Stop();
970 }
971 
972 
973 //=====================================================================
974 //* DistanceToOut (p) -------------------------------------------------
975 
977 {
978  // DistanceToOut(p):
979  // Calculate distance to surface of shape from `inside',
980  // allowing for tolerance
981 
982  //
983  // checking last value
984  //
985 
986  G4ThreeVector *tmpp;
987  G4double *tmpdist;
988 
989  if (fLastDistanceToOut.p == p)
990  {
991  return fLastDistanceToOut.value;
992  }
993  else
994  {
995  tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOut.p));
996  tmpdist = const_cast<G4double*>(&(fLastDistanceToOut.value));
997  tmpp->set(p.x(), p.y(), p.z());
998  }
999 
1000  //
1001  // Calculate DistanceToOut(p)
1002  //
1003 
1004  EInside currentside = Inside(p);
1005  G4double retval = kInfinity;
1006 
1007  switch (currentside)
1008  {
1009  case (kOutside) :
1010  {
1011 #ifdef G4SPECSDEBUG
1012  G4int oldprc = G4cout.precision(16) ;
1013  G4cout << G4endl ;
1014  DumpInfo();
1015  G4cout << "Position:" << G4endl << G4endl ;
1016  G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
1017  G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
1018  G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
1019  G4cout.precision(oldprc) ;
1020  G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids1002",
1021  JustWarning, "Point p is outside !?" );
1022 #endif
1023  break;
1024  }
1025  case (kSurface) :
1026  {
1027  *tmpdist = 0.;
1028  retval = fLastDistanceToOut.value;
1029  break;
1030  }
1031 
1032  case (kInside) :
1033  {
1034  // Initialize
1035  //
1036  G4double distance = kInfinity;
1037 
1038  // find intersections and choose nearest one
1039  //
1040  G4VTwistSurface *surfaces[6];
1041 
1042  surfaces[0] = fSide0;
1043  surfaces[1] = fSide90 ;
1044  surfaces[2] = fSide180 ;
1045  surfaces[3] = fSide270 ;
1046  surfaces[4] = fLowerEndcap;
1047  surfaces[5] = fUpperEndcap;
1048 
1049  G4int i;
1050  G4ThreeVector xx;
1051  G4ThreeVector bestxx;
1052  for (i=0; i< 6; i++)
1053  {
1054  G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
1055  if (tmpdistance < distance)
1056  {
1057  distance = tmpdistance;
1058  bestxx = xx;
1059  }
1060  }
1061  *tmpdist = distance;
1062 
1063  retval = fLastDistanceToOut.value;
1064  break;
1065  }
1066 
1067  default :
1068  {
1069  G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids0003",
1070  FatalException, "Unknown point location!");
1071  break;
1072  }
1073  } // switch end
1074 
1075  return retval;
1076 }
1077 
1078 
1079 //=====================================================================
1080 //* StreamInfo --------------------------------------------------------
1081 
1082 std::ostream& G4VTwistedFaceted::StreamInfo(std::ostream& os) const
1083 {
1084  //
1085  // Stream object contents to an output stream
1086  //
1087  G4int oldprc = os.precision(16);
1088  os << "-----------------------------------------------------------\n"
1089  << " *** Dump for solid - " << GetName() << " ***\n"
1090  << " ===================================================\n"
1091  << " Solid type: G4VTwistedFaceted\n"
1092  << " Parameters: \n"
1093  << " polar angle theta = " << fTheta/degree << " deg" << G4endl
1094  << " azimuthal angle phi = " << fPhi/degree << " deg" << G4endl
1095  << " tilt angle alpha = " << fAlph/degree << " deg" << G4endl
1096  << " TWIST angle = " << fPhiTwist/degree << " deg" << G4endl
1097  << " Half length along y (lower endcap) = " << fDy1/cm << " cm"
1098  << G4endl
1099  << " Half length along x (lower endcap, bottom) = " << fDx1/cm << " cm"
1100  << G4endl
1101  << " Half length along x (lower endcap, top) = " << fDx2/cm << " cm"
1102  << G4endl
1103  << " Half length along y (upper endcap) = " << fDy2/cm << " cm"
1104  << G4endl
1105  << " Half length along x (upper endcap, bottom) = " << fDx3/cm << " cm"
1106  << G4endl
1107  << " Half length along x (upper endcap, top) = " << fDx4/cm << " cm"
1108  << G4endl
1109  << "-----------------------------------------------------------\n";
1110  os.precision(oldprc);
1111 
1112  return os;
1113 }
1114 
1115 
1116 //=====================================================================
1117 //* DiscribeYourselfTo ------------------------------------------------
1118 
1120 {
1121  scene.AddSolid (*this);
1122 }
1123 
1124 //=====================================================================
1125 //* GetExtent ---------------------------------------------------------
1126 
1128 {
1129  G4double maxRad = std::sqrt( fDx*fDx + fDy*fDy);
1130 
1131  return G4VisExtent(-maxRad, maxRad ,
1132  -maxRad, maxRad ,
1133  -fDz, fDz );
1134 }
1135 
1136 
1137 //=====================================================================
1138 //* CreateSurfaces ----------------------------------------------------
1139 
1141 {
1142 
1143  // create 6 surfaces of TwistedTub.
1144 
1145  if ( fDx1 == fDx2 && fDx3 == fDx4 ) // special case : Box
1146  {
1147  fSide0 = new G4TwistBoxSide("0deg", fPhiTwist, fDz, fTheta, fPhi,
1148  fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 0.*deg);
1149  fSide180 = new G4TwistBoxSide("180deg", fPhiTwist, fDz, fTheta, fPhi+pi,
1150  fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 180.*deg);
1151  }
1152  else // default general case
1153  {
1155  fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
1156  fSide180 = new G4TwistTrapAlphaSide("180deg", fPhiTwist, fDz, fTheta,
1157  fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
1158  }
1159 
1160  // create parallel sides
1161  //
1163  fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
1165  fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
1166 
1167  // create endcaps
1168  //
1170  fDz, fAlph, fPhi, fTheta, 1 );
1172  fDz, fAlph, fPhi, fTheta, -1 );
1173 
1174  // Set neighbour surfaces
1175 
1182 
1183 }
1184 
1185 
1186 //=====================================================================
1187 //* GetEntityType -----------------------------------------------------
1188 
1190 {
1191  return G4String("G4VTwistedFaceted");
1192 }
1193 
1194 
1195 //=====================================================================
1196 //* GetPolyhedron -----------------------------------------------------
1197 
1199 {
1200  if (!fpPolyhedron ||
1203  fpPolyhedron->GetNumberOfRotationSteps())
1204  {
1205  G4AutoLock l(&polyhedronMutex);
1206  delete fpPolyhedron;
1208  fRebuildPolyhedron = false;
1209  l.unlock();
1210  }
1211 
1212  return fpPolyhedron;
1213 }
1214 
1215 
1216 //=====================================================================
1217 //* GetPointInSolid ---------------------------------------------------
1218 
1220 {
1221 
1222 
1223  // this routine is only used for a test
1224  // can be deleted ...
1225 
1226  if ( z == fDz ) z -= 0.1*fDz ;
1227  if ( z == -fDz ) z += 0.1*fDz ;
1228 
1229  G4double phi = z/(2*fDz)*fPhiTwist ;
1230 
1231  return G4ThreeVector(fdeltaX * phi/fPhiTwist, fdeltaY * phi/fPhiTwist, z ) ;
1232 }
1233 
1234 
1235 //=====================================================================
1236 //* GetPointOnSurface -------------------------------------------------
1237 
1239 {
1240 
1242  G4double u , umin, umax ; // variable for twisted surfaces
1243  G4double y ; // variable for flat surface (top and bottom)
1244 
1245  // Compute the areas. Attention: Only correct for trapezoids
1246  // where the twisting is done along the z-axis. In the general case
1247  // the computed surface area is more difficult. However this simplification
1248  // does not affect the tracking through the solid.
1249 
1256 
1257 #ifdef G4TWISTDEBUG
1258  G4cout << "Surface 0 deg = " << a1 << G4endl ;
1259  G4cout << "Surface 90 deg = " << a2 << G4endl ;
1260  G4cout << "Surface 180 deg = " << a3 << G4endl ;
1261  G4cout << "Surface 270 deg = " << a4 << G4endl ;
1262  G4cout << "Surface Lower = " << a5 << G4endl ;
1263  G4cout << "Surface Upper = " << a6 << G4endl ;
1264 #endif
1265 
1266  G4double chose = G4RandFlat::shoot(0.,a1 + a2 + a3 + a4 + a5 + a6) ;
1267 
1268  if(chose < a1)
1269  {
1270 
1271  umin = fSide0->GetBoundaryMin(phi) ;
1272  umax = fSide0->GetBoundaryMax(phi) ;
1273  u = G4RandFlat::shoot(umin,umax) ;
1274 
1275  return fSide0->SurfacePoint(phi, u, true) ; // point on 0deg surface
1276  }
1277 
1278  else if( (chose >= a1) && (chose < a1 + a2 ) )
1279  {
1280 
1281  umin = fSide90->GetBoundaryMin(phi) ;
1282  umax = fSide90->GetBoundaryMax(phi) ;
1283 
1284  u = G4RandFlat::shoot(umin,umax) ;
1285 
1286  return fSide90->SurfacePoint(phi, u, true); // point on 90deg surface
1287  }
1288 
1289  else if( (chose >= a1 + a2 ) && (chose < a1 + a2 + a3 ) )
1290  {
1291 
1292  umin = fSide180->GetBoundaryMin(phi) ;
1293  umax = fSide180->GetBoundaryMax(phi) ;
1294  u = G4RandFlat::shoot(umin,umax) ;
1295 
1296  return fSide180->SurfacePoint(phi, u, true); // point on 180 deg surface
1297  }
1298 
1299  else if( (chose >= a1 + a2 + a3 ) && (chose < a1 + a2 + a3 + a4 ) )
1300  {
1301 
1302  umin = fSide270->GetBoundaryMin(phi) ;
1303  umax = fSide270->GetBoundaryMax(phi) ;
1304  u = G4RandFlat::shoot(umin,umax) ;
1305 
1306  return fSide270->SurfacePoint(phi, u, true); // point on 270 deg surface
1307  }
1308 
1309  else if( (chose >= a1 + a2 + a3 + a4 ) && (chose < a1 + a2 + a3 + a4 + a5 ) )
1310  {
1311 
1312  y = G4RandFlat::shoot(-fDy1,fDy1) ;
1313  umin = fLowerEndcap->GetBoundaryMin(y) ;
1314  umax = fLowerEndcap->GetBoundaryMax(y) ;
1315  u = G4RandFlat::shoot(umin,umax) ;
1316 
1317  return fLowerEndcap->SurfacePoint(u,y,true); // point on lower endcap
1318  }
1319  else {
1320 
1321  y = G4RandFlat::shoot(-fDy2,fDy2) ;
1322  umin = fUpperEndcap->GetBoundaryMin(y) ;
1323  umax = fUpperEndcap->GetBoundaryMax(y) ;
1324  u = G4RandFlat::shoot(umin,umax) ;
1325 
1326  return fUpperEndcap->SurfacePoint(u,y,true) ; // point on upper endcap
1327 
1328  }
1329 }
1330 
1331 
1332 //=====================================================================
1333 //* CreatePolyhedron --------------------------------------------------
1334 
1336 {
1337  // number of meshes
1338  const G4int k =
1339  G4int(G4Polyhedron::GetNumberOfRotationSteps() * fPhiTwist / twopi) + 2;
1340  const G4int n = k;
1341 
1342  const G4int nnodes = 4*(k-1)*(n-2) + 2*k*k ;
1343  const G4int nfaces = 4*(k-1)*(n-1) + 2*(k-1)*(k-1) ;
1344 
1345  G4Polyhedron *ph=new G4Polyhedron;
1346  typedef G4double G4double3[3];
1347  typedef G4int G4int4[4];
1348  G4double3* xyz = new G4double3[nnodes]; // number of nodes
1349  G4int4* faces = new G4int4[nfaces] ; // number of faces
1350 
1351  fLowerEndcap->GetFacets(k,k,xyz,faces,0) ;
1352  fUpperEndcap->GetFacets(k,k,xyz,faces,1) ;
1353  fSide270->GetFacets(k,n,xyz,faces,2) ;
1354  fSide0->GetFacets(k,n,xyz,faces,3) ;
1355  fSide90->GetFacets(k,n,xyz,faces,4) ;
1356  fSide180->GetFacets(k,n,xyz,faces,5) ;
1357 
1358  ph->createPolyhedron(nnodes,nfaces,xyz,faces);
1359 
1360  return ph;
1361 }
static const double cm
Definition: G4SIunits.hh:118
G4String GetName() const
void ClipCrossSection(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:347
ThreeVector shoot(const G4int Ap, const G4int Af)
virtual void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4ThreeVectorList * CreateRotatedVertices(const G4AffineTransform &pTransform) const
G4VTwistSurface * fSide90
static const G4double kInfinity
Definition: geomdefs.hh:42
G4double GetMinYExtent() const
G4VTwistedFaceted(const G4String &pname, G4double PhiTwist, G4double pDz, G4double pTheta, G4double pPhi, G4double pDy1, G4double pDx1, G4double pDx2, G4double pDy2, G4double pDx3, G4double pDx4, G4double pAlph)
G4VTwistSurface * fUpperEndcap
CLHEP::Hep3Vector G4ThreeVector
LastValueWithDoubleVector fLastDistanceToOutWithV
virtual G4ThreeVector GetNormal(const G4ThreeVector &xx, G4bool isGlobal)=0
G4AffineTransform Inverse() const
G4double z
Definition: TRTMaterials.hh:39
G4bool IsYLimited() const
virtual G4GeometryType GetEntityType() const
static const G4double a1
G4ThreeVector GetPointOnSurface() const
G4ThreeVector GetPointInSolid(G4double z) const
G4bool IsRotated() const
virtual G4Polyhedron * CreatePolyhedron() const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4ThreeVector NetTranslation() const
static const G4double a4
G4bool IsXLimited() const
G4double a
Definition: TRTMaterials.hh:39
G4VTwistedFaceted & operator=(const G4VTwistedFaceted &rhs)
virtual void AddSolid(const G4Box &)=0
int G4int
Definition: G4Types.hh:78
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:175
G4bool IsValidNorm() const
G4VTwistSurface * fSide180
void DumpInfo() const
virtual G4double DistanceToIn(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
G4double Xcoef(G4double u, G4double phi, G4double ftg) const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
G4double GetMaxXExtent() const
virtual G4double GetSurfaceArea()=0
G4VTwistSurface * fSide270
G4double GetMinZExtent() const
G4GLOB_DLL std::ostream G4cout
G4bool IsLimited() const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
void SetNeighbours(G4VTwistSurface *axis0min, G4VTwistSurface *axis1min, G4VTwistSurface *axis0max, G4VTwistSurface *axis1max)
virtual G4double DistanceTo(const G4ThreeVector &gp, G4ThreeVector &gxx)
virtual G4double DistanceToOut(const G4ThreeVector &gp, const G4ThreeVector &gv, G4ThreeVector &gxxbest)
static const double deg
Definition: G4SIunits.hh:151
virtual G4ThreeVector SurfacePoint(G4double, G4double, G4bool isGlobal=false)=0
bool G4bool
Definition: G4Types.hh:79
static const double twopi
Definition: G4SIunits.hh:75
G4double GetValueB(G4double phi) const
G4VTwistSurface * fLowerEndcap
std::vector< G4ThreeVector > G4ThreeVectorList
Definition: G4VSolid.hh:79
const G4int n
LastValueWithDoubleVector fLastDistanceToInWithV
virtual G4double GetBoundaryMax(G4double)=0
static const G4double a3
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual void GetFacets(G4int m, G4int n, G4double xyz[][3], G4int faces[][4], G4int iside)=0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=0, G4ThreeVector *n=0) const
virtual G4double GetBoundaryMin(G4double)=0
G4double GetMinXExtent() const
G4int G4Mutex
Definition: G4Threading.hh:173
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double GetMaxZExtent() const
virtual G4Polyhedron * GetPolyhedron() const
static const double pi
Definition: G4SIunits.hh:74
EInside
Definition: geomdefs.hh:58
virtual G4VisExtent GetExtent() const
virtual EInside Inside(const G4ThreeVector &p) const
EAxis
Definition: geomdefs.hh:54
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
static const double degree
Definition: G4SIunits.hh:143
#define G4endl
Definition: G4ios.hh:61
G4double GetMaxYExtent() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:110
G4double kCarTolerance
Definition: G4VSolid.hh:304
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4VTwistSurface * fSide0
static const G4double a5
double G4double
Definition: G4Types.hh:76
virtual std::ostream & StreamInfo(std::ostream &os) const
G4double GetMaxExtent(const EAxis pAxis) const
G4bool IsZLimited() const
static const double mm
Definition: G4SIunits.hh:114
G4Polyhedron * fpPolyhedron
G4double GetAngularTolerance() const
void ClipBetweenSections(G4ThreeVectorList *pVertices, const G4int pSectionIndex, const G4VoxelLimits &pVoxelLimit, const EAxis pAxis, G4double &pMin, G4double &pMax) const
Definition: G4VSolid.cc:378
G4double GetMinExtent(const EAxis pAxis) const
static const G4double a2
static G4GeometryTolerance * GetInstance()