Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 (Oliver.Link@cern.ch)
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{
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(0), fUpperEndcap(0), fSide0(0),
75 fSide90(0), fSide180(0), fSide270(0)
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 }
184}
185
186
187//=====================================================================
188//* Fake default constructor ------------------------------------------
189
191 : G4VSolid(a),
192 fTheta(0.), fPhi(0.), fDy1(0.),
193 fDx1(0.), fDx2(0.), fDy2(0.), fDx3(0.), fDx4(0.),
194 fDz(0.), fDx(0.), fDy(0.), fAlph(0.),
195 fTAlph(0.), fdeltaX(0.), fdeltaY(0.), fPhiTwist(0.),
196 fLowerEndcap(0), fUpperEndcap(0), fSide0(0), fSide90(0), fSide180(0),
197 fSide270(0)
198{
199}
200
201
202//=====================================================================
203//* destructor --------------------------------------------------------
204
206{
207 if (fLowerEndcap) { delete fLowerEndcap ; }
208 if (fUpperEndcap) { delete fUpperEndcap ; }
209
210 if (fSide0) { delete fSide0 ; }
211 if (fSide90) { delete fSide90 ; }
212 if (fSide180) { delete fSide180 ; }
213 if (fSide270) { delete fSide270 ; }
214 if (fpPolyhedron) { delete fpPolyhedron; fpPolyhedron = nullptr; }
215}
216
217
218//=====================================================================
219//* Copy constructor --------------------------------------------------
220
222 : G4VSolid(rhs),
223 fTheta(rhs.fTheta), fPhi(rhs.fPhi),
224 fDy1(rhs.fDy1), fDx1(rhs.fDx1), fDx2(rhs.fDx2), fDy2(rhs.fDy2),
225 fDx3(rhs.fDx3), fDx4(rhs.fDx4), fDz(rhs.fDz), fDx(rhs.fDx), fDy(rhs.fDy),
226 fAlph(rhs.fAlph), fTAlph(rhs.fTAlph), fdeltaX(rhs.fdeltaX),
227 fdeltaY(rhs.fdeltaY), fPhiTwist(rhs.fPhiTwist), fLowerEndcap(0),
228 fUpperEndcap(0), fSide0(0), fSide90(0), fSide180(0), fSide270(0),
229 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
230 fLastInside(rhs.fLastInside), fLastNormal(rhs.fLastNormal),
231 fLastDistanceToIn(rhs.fLastDistanceToIn),
232 fLastDistanceToOut(rhs.fLastDistanceToOut),
233 fLastDistanceToInWithV(rhs.fLastDistanceToInWithV),
234 fLastDistanceToOutWithV(rhs.fLastDistanceToOutWithV)
235{
237}
238
239
240//=====================================================================
241//* Assignment operator -----------------------------------------------
242
244{
245 // Check assignment to self
246 //
247 if (this == &rhs) { return *this; }
248
249 // Copy base class data
250 //
252
253 // Copy data
254 //
255 fTheta = rhs.fTheta; fPhi = rhs.fPhi;
256 fDy1= rhs.fDy1; fDx1= rhs.fDx1; fDx2= rhs.fDx2; fDy2= rhs.fDy2;
257 fDx3= rhs.fDx3; fDx4= rhs.fDx4; fDz= rhs.fDz; fDx= rhs.fDx; fDy= rhs.fDy;
258 fAlph= rhs.fAlph; fTAlph= rhs.fTAlph; fdeltaX= rhs.fdeltaX;
260 fUpperEndcap= 0; fSide0= 0; fSide90= 0; fSide180= 0; fSide270= 0;
262 fRebuildPolyhedron = false;
263 delete fpPolyhedron; fpPolyhedron = nullptr;
269
271
272 return *this;
273}
274
275
276//=====================================================================
277//* ComputeDimensions -------------------------------------------------
278
280 const G4int ,
281 const G4VPhysicalVolume* )
282{
283 G4Exception("G4VTwistedFaceted::ComputeDimensions()",
284 "GeomSolids0001", FatalException,
285 "G4VTwistedFaceted does not support Parameterisation.");
286}
287
288
289//=====================================================================
290//* Extent ------------------------------------------------------------
291
293 G4ThreeVector& pMax) const
294{
295 G4double maxRad = std::sqrt(fDx*fDx + fDy*fDy);
296 pMin.set(-maxRad,-maxRad,-fDz);
297 pMax.set( maxRad, maxRad, fDz);
298}
299
300
301//=====================================================================
302//* CalculateExtent ---------------------------------------------------
303
304G4bool
306 const G4VoxelLimits& pVoxelLimit,
307 const G4AffineTransform& pTransform,
308 G4double& pMin,
309 G4double& pMax ) const
310{
311 G4ThreeVector bmin, bmax;
312
313 // Get bounding box
314 BoundingLimits(bmin,bmax);
315
316 // Find extent
317 G4BoundingEnvelope bbox(bmin,bmax);
318 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
319}
320
321
322//=====================================================================
323//* Inside ------------------------------------------------------------
324
326{
327
328 G4ThreeVector *tmpp;
329 EInside *tmpin;
330 if (fLastInside.p == p)
331 {
332 return fLastInside.inside;
333 }
334 else
335 {
336 tmpp = const_cast<G4ThreeVector*>(&(fLastInside.p));
337 tmpin = const_cast<EInside*>(&(fLastInside.inside));
338 tmpp->set(p.x(), p.y(), p.z());
339 }
340
341 *tmpin = kOutside ;
342
343 G4double phi = p.z()/(2*fDz) * fPhiTwist ; // rotate the point to z=0
344 G4double cphi = std::cos(-phi) ;
345 G4double sphi = std::sin(-phi) ;
346
347 G4double px = p.x() + fdeltaX * ( -phi/fPhiTwist) ; // shift
348 G4double py = p.y() + fdeltaY * ( -phi/fPhiTwist) ;
349 G4double pz = p.z() ;
350
351 G4double posx = px * cphi - py * sphi ; // rotation
352 G4double posy = px * sphi + py * cphi ;
353 G4double posz = pz ;
354
355 G4double xMin = Xcoef(posy,phi,fTAlph) - 2*Xcoef(posy,phi,0.) ;
356 G4double xMax = Xcoef(posy,phi,fTAlph) ;
357
358 G4double yMax = GetValueB(phi)/2. ; // b(phi)/2 is limit
359 G4double yMin = -yMax ;
360
361#ifdef G4TWISTDEBUG
362
363 G4cout << "inside called: p = " << p << G4endl ;
364 G4cout << "fDx1 = " << fDx1 << G4endl ;
365 G4cout << "fDx2 = " << fDx2 << G4endl ;
366 G4cout << "fDx3 = " << fDx3 << G4endl ;
367 G4cout << "fDx4 = " << fDx4 << G4endl ;
368
369 G4cout << "fDy1 = " << fDy1 << G4endl ;
370 G4cout << "fDy2 = " << fDy2 << G4endl ;
371
372 G4cout << "fDz = " << fDz << G4endl ;
373
374 G4cout << "Tilt angle alpha = " << fAlph << G4endl ;
375 G4cout << "phi,theta = " << fPhi << " , " << fTheta << G4endl ;
376
377 G4cout << "Twist angle = " << fPhiTwist << G4endl ;
378
379 G4cout << "posx = " << posx << G4endl ;
380 G4cout << "posy = " << posy << G4endl ;
381 G4cout << "xMin = " << xMin << G4endl ;
382 G4cout << "xMax = " << xMax << G4endl ;
383 G4cout << "yMin = " << yMin << G4endl ;
384 G4cout << "yMax = " << yMax << G4endl ;
385
386#endif
387
388
389 if ( posx <= xMax - kCarTolerance*0.5
390 && posx >= xMin + kCarTolerance*0.5 )
391 {
392 if ( posy <= yMax - kCarTolerance*0.5
393 && posy >= yMin + kCarTolerance*0.5 )
394 {
395 if (std::fabs(posz) <= fDz - kCarTolerance*0.5 ) *tmpin = kInside ;
396 else if (std::fabs(posz) <= fDz + kCarTolerance*0.5 ) *tmpin = kSurface ;
397 }
398 else if ( posy <= yMax + kCarTolerance*0.5
399 && posy >= yMin - kCarTolerance*0.5 )
400 {
401 if (std::fabs(posz) <= fDz + kCarTolerance*0.5 ) *tmpin = kSurface ;
402 }
403 }
404 else if ( posx <= xMax + kCarTolerance*0.5
405 && posx >= xMin - kCarTolerance*0.5 )
406 {
407 if ( posy <= yMax + kCarTolerance*0.5
408 && posy >= yMin - kCarTolerance*0.5 )
409 {
410 if (std::fabs(posz) <= fDz + kCarTolerance*0.5) *tmpin = kSurface ;
411 }
412 }
413
414#ifdef G4TWISTDEBUG
415 G4cout << "inside = " << fLastInside.inside << G4endl ;
416#endif
417
418 return fLastInside.inside;
419
420}
421
422
423//=====================================================================
424//* SurfaceNormal -----------------------------------------------------
425
427{
428 //
429 // return the normal unit vector to the Hyperbolical Surface at a point
430 // p on (or nearly on) the surface
431 //
432 // Which of the three or four surfaces are we closest to?
433 //
434
435 if (fLastNormal.p == p)
436 {
437 return fLastNormal.vec;
438 }
439
440 G4ThreeVector* tmpp = const_cast<G4ThreeVector*>(&(fLastNormal.p));
441 G4ThreeVector* tmpnormal = const_cast<G4ThreeVector*>(&(fLastNormal.vec));
442 G4VTwistSurface** tmpsurface
443 = const_cast<G4VTwistSurface**>(fLastNormal.surface);
444 tmpp->set(p.x(), p.y(), p.z());
445
446 G4double distance = kInfinity;
447
448 G4VTwistSurface* surfaces[6];
449
450 surfaces[0] = fSide0 ;
451 surfaces[1] = fSide90 ;
452 surfaces[2] = fSide180 ;
453 surfaces[3] = fSide270 ;
454 surfaces[4] = fLowerEndcap;
455 surfaces[5] = fUpperEndcap;
456
457 G4ThreeVector xx;
458 G4ThreeVector bestxx;
459 G4int i;
460 G4int besti = -1;
461 for (i=0; i< 6; i++)
462 {
463 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
464 if (tmpdistance < distance)
465 {
466 distance = tmpdistance;
467 bestxx = xx;
468 besti = i;
469 }
470 }
471
472 tmpsurface[0] = surfaces[besti];
473 *tmpnormal = tmpsurface[0]->GetNormal(bestxx, true);
474
475 return fLastNormal.vec;
476}
477
478
479//=====================================================================
480//* DistanceToIn (p, v) -----------------------------------------------
481
483 const G4ThreeVector& v ) const
484{
485
486 // DistanceToIn (p, v):
487 // Calculate distance to surface of shape from `outside'
488 // along with the v, allowing for tolerance.
489 // The function returns kInfinity if no intersection or
490 // just grazing within tolerance.
491
492 //
493 // checking last value
494 //
495
496 G4ThreeVector* tmpp;
497 G4ThreeVector* tmpv;
498 G4double* tmpdist;
500 {
502 }
503 else
504 {
505 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.p));
506 tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToInWithV.vec));
507 tmpdist = const_cast<G4double*>(&(fLastDistanceToInWithV.value));
508 tmpp->set(p.x(), p.y(), p.z());
509 tmpv->set(v.x(), v.y(), v.z());
510 }
511
512 //
513 // Calculate DistanceToIn(p,v)
514 //
515
516 EInside currentside = Inside(p);
517
518 if (currentside == kInside)
519 {
520 }
521 else if (currentside == kSurface)
522 {
523 // particle is just on a boundary.
524 // if the particle is entering to the volume, return 0
525 //
527 if (normal*v < 0)
528 {
529 *tmpdist = 0.;
531 }
532 }
533
534 // now, we can take smallest positive distance.
535
536 // Initialize
537 //
538 G4double distance = kInfinity;
539
540 // Find intersections and choose nearest one
541 //
542 G4VTwistSurface *surfaces[6];
543
544 surfaces[0] = fSide0;
545 surfaces[1] = fSide90 ;
546 surfaces[2] = fSide180 ;
547 surfaces[3] = fSide270 ;
548 surfaces[4] = fLowerEndcap;
549 surfaces[5] = fUpperEndcap;
550
551 G4ThreeVector xx;
552 G4ThreeVector bestxx;
553 for (auto i=0; i < 6 ; ++i)
554 {
555#ifdef G4TWISTDEBUG
556 G4cout << G4endl << "surface " << i << ": " << G4endl << G4endl ;
557#endif
558 G4double tmpdistance = surfaces[i]->DistanceToIn(p, v, xx);
559#ifdef G4TWISTDEBUG
560 G4cout << "Solid DistanceToIn : distance = " << tmpdistance << G4endl ;
561 G4cout << "intersection point = " << xx << G4endl ;
562#endif
563 if (tmpdistance < distance)
564 {
565 distance = tmpdistance;
566 bestxx = xx;
567 }
568 }
569
570#ifdef G4TWISTDEBUG
571 G4cout << "best distance = " << distance << G4endl ;
572#endif
573
574 *tmpdist = distance;
575 // timer.Stop();
577}
578
579
580//=====================================================================
581//* DistanceToIn (p) --------------------------------------------------
582
584{
585 // DistanceToIn(p):
586 // Calculate distance to surface of shape from `outside',
587 // allowing for tolerance
588 //
589
590 //
591 // checking last value
592 //
593
594 G4ThreeVector* tmpp;
595 G4double* tmpdist;
596 if (fLastDistanceToIn.p == p)
597 {
599 }
600 else
601 {
602 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToIn.p));
603 tmpdist = const_cast<G4double*>(&(fLastDistanceToIn.value));
604 tmpp->set(p.x(), p.y(), p.z());
605 }
606
607 //
608 // Calculate DistanceToIn(p)
609 //
610
611 EInside currentside = Inside(p);
612
613 switch (currentside)
614 {
615 case (kInside) :
616 {
617 }
618
619 case (kSurface) :
620 {
621 *tmpdist = 0.;
623 }
624
625 case (kOutside) :
626 {
627 // Initialize
628 //
629 G4double distance = kInfinity;
630
631 // Find intersections and choose nearest one
632 //
633 G4VTwistSurface* surfaces[6];
634
635 surfaces[0] = fSide0;
636 surfaces[1] = fSide90 ;
637 surfaces[2] = fSide180 ;
638 surfaces[3] = fSide270 ;
639 surfaces[4] = fLowerEndcap;
640 surfaces[5] = fUpperEndcap;
641
642 G4ThreeVector xx;
643 G4ThreeVector bestxx;
644 for (auto i=0; i< 6; ++i)
645 {
646 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
647 if (tmpdistance < distance)
648 {
649 distance = tmpdistance;
650 bestxx = xx;
651 }
652 }
653 *tmpdist = distance;
655 }
656
657 default:
658 {
659 G4Exception("G4VTwistedFaceted::DistanceToIn(p)", "GeomSolids0003",
660 FatalException, "Unknown point location!");
661 }
662 } // switch end
663
664 return 0.;
665}
666
667
668//=====================================================================
669//* DistanceToOut (p, v) ----------------------------------------------
670
673 const G4ThreeVector& v,
674 const G4bool calcNorm,
675 G4bool* validNorm,
676 G4ThreeVector* norm ) const
677{
678 // DistanceToOut (p, v):
679 // Calculate distance to surface of shape from `inside'
680 // along with the v, allowing for tolerance.
681 // The function returns kInfinity if no intersection or
682 // just grazing within tolerance.
683
684 //
685 // checking last value
686 //
687
688 G4ThreeVector* tmpp;
689 G4ThreeVector* tmpv;
690 G4double* tmpdist;
692 {
694 }
695 else
696 {
697 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.p));
698 tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.vec));
699 tmpdist = const_cast<G4double*>(&(fLastDistanceToOutWithV.value));
700 tmpp->set(p.x(), p.y(), p.z());
701 tmpv->set(v.x(), v.y(), v.z());
702 }
703
704 //
705 // Calculate DistanceToOut(p,v)
706 //
707
708 EInside currentside = Inside(p);
709
710 if (currentside == kOutside)
711 {
712 }
713 else if (currentside == kSurface)
714 {
715 // particle is just on a boundary.
716 // if the particle is exiting from the volume, return 0
717 //
719 G4VTwistSurface *blockedsurface = fLastNormal.surface[0];
720 if (normal*v > 0)
721 {
722 if (calcNorm)
723 {
724 *norm = (blockedsurface->GetNormal(p, true));
725 *validNorm = blockedsurface->IsValidNorm();
726 }
727 *tmpdist = 0.;
728 // timer.Stop();
730 }
731 }
732
733 // now, we can take smallest positive distance.
734
735 // Initialize
736 G4double distance = kInfinity;
737
738 // find intersections and choose nearest one.
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 G4int besti = -1;
749 G4ThreeVector xx;
750 G4ThreeVector bestxx;
751 for (auto i=0; i<6 ; ++i)
752 {
753 G4double tmpdistance = surfaces[i]->DistanceToOut(p, v, xx);
754 if (tmpdistance < distance)
755 {
756 distance = tmpdistance;
757 bestxx = xx;
758 besti = i;
759 }
760 }
761
762 if (calcNorm)
763 {
764 if (besti != -1)
765 {
766 *norm = (surfaces[besti]->GetNormal(p, true));
767 *validNorm = surfaces[besti]->IsValidNorm();
768 }
769 }
770
771 *tmpdist = distance;
773}
774
775
776//=====================================================================
777//* DistanceToOut (p) -------------------------------------------------
778
780{
781 // DistanceToOut(p):
782 // Calculate distance to surface of shape from `inside',
783 // allowing for tolerance
784
785 //
786 // checking last value
787 //
788
789 G4ThreeVector* tmpp;
790 G4double* tmpdist;
791
792 if (fLastDistanceToOut.p == p)
793 {
795 }
796 else
797 {
798 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOut.p));
799 tmpdist = const_cast<G4double*>(&(fLastDistanceToOut.value));
800 tmpp->set(p.x(), p.y(), p.z());
801 }
802
803 //
804 // Calculate DistanceToOut(p)
805 //
806
807 EInside currentside = Inside(p);
808 G4double retval = kInfinity;
809
810 switch (currentside)
811 {
812 case (kOutside) :
813 {
814#ifdef G4SPECSDEBUG
815 G4int oldprc = G4cout.precision(16) ;
816 G4cout << G4endl ;
817 DumpInfo();
818 G4cout << "Position:" << G4endl << G4endl ;
819 G4cout << "p.x() = " << p.x()/mm << " mm" << G4endl ;
820 G4cout << "p.y() = " << p.y()/mm << " mm" << G4endl ;
821 G4cout << "p.z() = " << p.z()/mm << " mm" << G4endl << G4endl ;
822 G4cout.precision(oldprc) ;
823 G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids1002",
824 JustWarning, "Point p is outside !?" );
825#endif
826 break;
827 }
828 case (kSurface) :
829 {
830 *tmpdist = 0.;
831 retval = fLastDistanceToOut.value;
832 break;
833 }
834
835 case (kInside) :
836 {
837 // Initialize
838 //
839 G4double distance = kInfinity;
840
841 // find intersections and choose nearest one
842 //
843 G4VTwistSurface* surfaces[6];
844
845 surfaces[0] = fSide0;
846 surfaces[1] = fSide90 ;
847 surfaces[2] = fSide180 ;
848 surfaces[3] = fSide270 ;
849 surfaces[4] = fLowerEndcap;
850 surfaces[5] = fUpperEndcap;
851
852 G4ThreeVector xx;
853 G4ThreeVector bestxx;
854 for (auto i=0; i<6; ++i)
855 {
856 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
857 if (tmpdistance < distance)
858 {
859 distance = tmpdistance;
860 bestxx = xx;
861 }
862 }
863 *tmpdist = distance;
864
865 retval = fLastDistanceToOut.value;
866 break;
867 }
868
869 default :
870 {
871 G4Exception("G4VTwistedFaceted::DistanceToOut(p)", "GeomSolids0003",
872 FatalException, "Unknown point location!");
873 break;
874 }
875 } // switch end
876
877 return retval;
878}
879
880
881//=====================================================================
882//* StreamInfo --------------------------------------------------------
883
884std::ostream& G4VTwistedFaceted::StreamInfo(std::ostream& os) const
885{
886 //
887 // Stream object contents to an output stream
888 //
889 G4int oldprc = os.precision(16);
890 os << "-----------------------------------------------------------\n"
891 << " *** Dump for solid - " << GetName() << " ***\n"
892 << " ===================================================\n"
893 << " Solid type: G4VTwistedFaceted\n"
894 << " Parameters: \n"
895 << " polar angle theta = " << fTheta/degree << " deg" << G4endl
896 << " azimuthal angle phi = " << fPhi/degree << " deg" << G4endl
897 << " tilt angle alpha = " << fAlph/degree << " deg" << G4endl
898 << " TWIST angle = " << fPhiTwist/degree << " deg" << G4endl
899 << " Half length along y (lower endcap) = " << fDy1/cm << " cm"
900 << G4endl
901 << " Half length along x (lower endcap, bottom) = " << fDx1/cm << " cm"
902 << G4endl
903 << " Half length along x (lower endcap, top) = " << fDx2/cm << " cm"
904 << G4endl
905 << " Half length along y (upper endcap) = " << fDy2/cm << " cm"
906 << G4endl
907 << " Half length along x (upper endcap, bottom) = " << fDx3/cm << " cm"
908 << G4endl
909 << " Half length along x (upper endcap, top) = " << fDx4/cm << " cm"
910 << G4endl
911 << "-----------------------------------------------------------\n";
912 os.precision(oldprc);
913
914 return os;
915}
916
917
918//=====================================================================
919//* DiscribeYourselfTo ------------------------------------------------
920
922{
923 scene.AddSolid (*this);
924}
925
926
927//=====================================================================
928//* GetExtent ---------------------------------------------------------
929
931{
932 G4double maxRad = std::sqrt( fDx*fDx + fDy*fDy);
933
934 return G4VisExtent(-maxRad, maxRad ,
935 -maxRad, maxRad ,
936 -fDz, fDz );
937}
938
939
940//=====================================================================
941//* CreateSurfaces ----------------------------------------------------
942
944{
945
946 // create 6 surfaces of TwistedTub.
947
948 if ( fDx1 == fDx2 && fDx3 == fDx4 ) // special case : Box
949 {
950 fSide0 = new G4TwistBoxSide("0deg", fPhiTwist, fDz, fTheta, fPhi,
951 fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 0.*deg);
952 fSide180 = new G4TwistBoxSide("180deg", fPhiTwist, fDz, fTheta, fPhi+pi,
953 fDy1, fDx1, fDx1, fDy2, fDx3, fDx3, fAlph, 180.*deg);
954 }
955 else // default general case
956 {
958 fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
960 fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
961 }
962
963 // create parallel sides
964 //
966 fPhi, fDy1, fDx1, fDx2, fDy2, fDx3, fDx4, fAlph, 0.*deg);
968 fPhi+pi, fDy1, fDx2, fDx1, fDy2, fDx4, fDx3, fAlph, 180.*deg);
969
970 // create endcaps
971 //
973 fDz, fAlph, fPhi, fTheta, 1 );
975 fDz, fAlph, fPhi, fTheta, -1 );
976
977 // Set neighbour surfaces
978
985
986}
987
988
989//=====================================================================
990//* GetEntityType -----------------------------------------------------
991
993{
994 return G4String("G4VTwistedFaceted");
995}
996
997
998//=====================================================================
999//* GetPolyhedron -----------------------------------------------------
1000
1002{
1003 if (fpPolyhedron == nullptr ||
1007 {
1009 delete fpPolyhedron;
1011 fRebuildPolyhedron = false;
1012 l.unlock();
1013 }
1014
1015 return fpPolyhedron;
1016}
1017
1018
1019//=====================================================================
1020//* GetPointInSolid ---------------------------------------------------
1021
1023{
1024
1025
1026 // this routine is only used for a test
1027 // can be deleted ...
1028
1029 if ( z == fDz ) z -= 0.1*fDz ;
1030 if ( z == -fDz ) z += 0.1*fDz ;
1031
1032 G4double phi = z/(2*fDz)*fPhiTwist ;
1033
1034 return G4ThreeVector(fdeltaX * phi/fPhiTwist, fdeltaY * phi/fPhiTwist, z ) ;
1035}
1036
1037
1038//=====================================================================
1039//* GetPointOnSurface -------------------------------------------------
1040
1042{
1043
1045 G4double u , umin, umax ; // variable for twisted surfaces
1046 G4double y ; // variable for flat surface (top and bottom)
1047
1048 // Compute the areas. Attention: Only correct for trapezoids
1049 // where the twisting is done along the z-axis. In the general case
1050 // the computed surface area is more difficult. However this simplification
1051 // does not affect the tracking through the solid.
1052
1059
1060#ifdef G4TWISTDEBUG
1061 G4cout << "Surface 0 deg = " << a1 << G4endl ;
1062 G4cout << "Surface 90 deg = " << a2 << G4endl ;
1063 G4cout << "Surface 180 deg = " << a3 << G4endl ;
1064 G4cout << "Surface 270 deg = " << a4 << G4endl ;
1065 G4cout << "Surface Lower = " << a5 << G4endl ;
1066 G4cout << "Surface Upper = " << a6 << G4endl ;
1067#endif
1068
1069 G4double chose = G4RandFlat::shoot(0.,a1 + a2 + a3 + a4 + a5 + a6) ;
1070
1071 if(chose < a1)
1072 {
1073 umin = fSide0->GetBoundaryMin(phi) ;
1074 umax = fSide0->GetBoundaryMax(phi) ;
1075 u = G4RandFlat::shoot(umin,umax) ;
1076
1077 return fSide0->SurfacePoint(phi, u, true) ; // point on 0deg surface
1078 }
1079
1080 else if( (chose >= a1) && (chose < a1 + a2 ) )
1081 {
1082 umin = fSide90->GetBoundaryMin(phi) ;
1083 umax = fSide90->GetBoundaryMax(phi) ;
1084
1085 u = G4RandFlat::shoot(umin,umax) ;
1086
1087 return fSide90->SurfacePoint(phi, u, true); // point on 90deg surface
1088 }
1089 else if( (chose >= a1 + a2 ) && (chose < a1 + a2 + a3 ) )
1090 {
1091 umin = fSide180->GetBoundaryMin(phi) ;
1092 umax = fSide180->GetBoundaryMax(phi) ;
1093 u = G4RandFlat::shoot(umin,umax) ;
1094
1095 return fSide180->SurfacePoint(phi, u, true); // point on 180 deg surface
1096 }
1097 else if( (chose >= a1 + a2 + a3 ) && (chose < a1 + a2 + a3 + a4 ) )
1098 {
1099 umin = fSide270->GetBoundaryMin(phi) ;
1100 umax = fSide270->GetBoundaryMax(phi) ;
1101 u = G4RandFlat::shoot(umin,umax) ;
1102 return fSide270->SurfacePoint(phi, u, true); // point on 270 deg surface
1103 }
1104 else if( (chose >= a1 + a2 + a3 + a4 ) && (chose < a1 + a2 + a3 + a4 + a5 ) )
1105 {
1107 umin = fLowerEndcap->GetBoundaryMin(y) ;
1108 umax = fLowerEndcap->GetBoundaryMax(y) ;
1109 u = G4RandFlat::shoot(umin,umax) ;
1110
1111 return fLowerEndcap->SurfacePoint(u,y,true); // point on lower endcap
1112 }
1113 else
1114 {
1116 umin = fUpperEndcap->GetBoundaryMin(y) ;
1117 umax = fUpperEndcap->GetBoundaryMax(y) ;
1118 u = G4RandFlat::shoot(umin,umax) ;
1119
1120 return fUpperEndcap->SurfacePoint(u,y,true) ; // point on upper endcap
1121
1122 }
1123}
1124
1125
1126//=====================================================================
1127//* CreatePolyhedron --------------------------------------------------
1128
1130{
1131 // number of meshes
1132 const G4int k =
1134 std::abs(fPhiTwist) / twopi) + 2;
1135 const G4int n = k;
1136
1137 const G4int nnodes = 4*(k-1)*(n-2) + 2*k*k ;
1138 const G4int nfaces = 4*(k-1)*(n-1) + 2*(k-1)*(k-1) ;
1139
1140 G4Polyhedron* ph = new G4Polyhedron;
1141 typedef G4double G4double3[3];
1142 typedef G4int G4int4[4];
1143 G4double3* xyz = new G4double3[nnodes]; // number of nodes
1144 G4int4* faces = new G4int4[nfaces] ; // number of faces
1145
1146 fLowerEndcap->GetFacets(k,k,xyz,faces,0) ;
1147 fUpperEndcap->GetFacets(k,k,xyz,faces,1) ;
1148 fSide270->GetFacets(k,n,xyz,faces,2) ;
1149 fSide0->GetFacets(k,n,xyz,faces,3) ;
1150 fSide90->GetFacets(k,n,xyz,faces,4) ;
1151 fSide180->GetFacets(k,n,xyz,faces,5) ;
1152
1153 ph->createPolyhedron(nnodes,nfaces,xyz,faces);
1154
1155 return ph;
1156}
@ JustWarning
@ FatalException
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static const G4double pMax
static const G4double pMin
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double mm
Definition: G4SIunits.hh:95
static constexpr double degree
Definition: G4SIunits.hh:124
static constexpr double pi
Definition: G4SIunits.hh:55
static constexpr double cm
Definition: G4SIunits.hh:99
static constexpr double deg
Definition: G4SIunits.hh:132
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
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 G4double DistanceToIn(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
G4VTwistSurface * fSide180
virtual void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4VTwistedFaceted & operator=(const G4VTwistedFaceted &rhs)
G4Polyhedron * fpPolyhedron
virtual G4Polyhedron * GetPolyhedron() const
G4VTwistSurface * fSide0
G4ThreeVector GetPointOnSurface() const
virtual G4GeometryType GetEntityType() const
G4VTwistSurface * fSide90
G4ThreeVector GetPointInSolid(G4double z) const
virtual G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4VTwistSurface * fLowerEndcap
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const
G4VTwistSurface * fSide270
virtual EInside Inside(const G4ThreeVector &p) const
virtual std::ostream & StreamInfo(std::ostream &os) const
virtual void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
virtual G4VisExtent GetExtent() const
virtual G4Polyhedron * CreatePolyhedron() const
G4VTwistSurface * fUpperEndcap
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)
LastValueWithDoubleVector fLastDistanceToOutWithV
virtual void DescribeYourselfTo(G4VGraphicsScene &scene) const
LastValueWithDoubleVector fLastDistanceToInWithV
G4double GetValueB(G4double phi) const
G4double Xcoef(G4double u, G4double phi, G4double ftg) const
static G4int GetNumberOfRotationSteps()
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
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
static const G4double kInfinity
Definition: geomdefs.hh:41
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
ThreeVector shoot(const G4int Ap, const G4int Af)
string pname
Definition: eplot.py:33