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