Geant4-11
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TwistedTubs.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4TwistedTubs implementation
27//
28// 01-Aug-2002 - Kotoyo Hoshina (hoshina@hepburn.s.chiba-u.ac.jp), created.
29// 13-Nov-2003 - O.Link (Oliver.Link@cern.ch), Integration in Geant4
30// from original version in Jupiter-2.5.02 application.
31// --------------------------------------------------------------------
32
33#include "G4TwistedTubs.hh"
34
35#include "G4GeomTools.hh"
37#include "G4SystemOfUnits.hh"
39#include "G4VoxelLimits.hh"
40#include "G4AffineTransform.hh"
41#include "G4BoundingEnvelope.hh"
42#include "G4ClippablePolygon.hh"
44#include "meshdefs.hh"
45
46#include "G4VGraphicsScene.hh"
47#include "G4Polyhedron.hh"
48#include "G4VisExtent.hh"
49
50#include "Randomize.hh"
51
52#include "G4AutoLock.hh"
53
54namespace
55{
57}
58
59//=====================================================================
60//* constructors ------------------------------------------------------
61
63 G4double twistedangle,
64 G4double endinnerrad,
65 G4double endouterrad,
66 G4double halfzlen,
67 G4double dphi)
68 : G4VSolid(pname), fDPhi(dphi),
69 fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
70 fFormerTwisted(0), fInnerHype(0), fOuterHype(0)
71{
72 if (endinnerrad < DBL_MIN)
73 {
74 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
75 FatalErrorInArgument, "Invalid end-inner-radius!");
76 }
77
78 G4double sinhalftwist = std::sin(0.5 * twistedangle);
79
80 G4double endinnerradX = endinnerrad * sinhalftwist;
81 G4double innerrad = std::sqrt( endinnerrad * endinnerrad
82 - endinnerradX * endinnerradX );
83
84 G4double endouterradX = endouterrad * sinhalftwist;
85 G4double outerrad = std::sqrt( endouterrad * endouterrad
86 - endouterradX * endouterradX );
87
88 // temporary treatment!!
89 SetFields(twistedangle, innerrad, outerrad, -halfzlen, halfzlen);
91}
92
94 G4double twistedangle,
95 G4double endinnerrad,
96 G4double endouterrad,
97 G4double halfzlen,
98 G4int nseg,
99 G4double totphi)
100 : G4VSolid(pname),
101 fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
102 fFormerTwisted(0), fInnerHype(0), fOuterHype(0)
103{
104
105 if (!nseg)
106 {
107 std::ostringstream message;
108 message << "Invalid number of segments." << G4endl
109 << " nseg = " << nseg;
110 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
111 FatalErrorInArgument, message);
112 }
113 if (totphi == DBL_MIN || endinnerrad < DBL_MIN)
114 {
115 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
116 FatalErrorInArgument, "Invalid total-phi or end-inner-radius!");
117 }
118
119 G4double sinhalftwist = std::sin(0.5 * twistedangle);
120
121 G4double endinnerradX = endinnerrad * sinhalftwist;
122 G4double innerrad = std::sqrt( endinnerrad * endinnerrad
123 - endinnerradX * endinnerradX );
124
125 G4double endouterradX = endouterrad * sinhalftwist;
126 G4double outerrad = std::sqrt( endouterrad * endouterrad
127 - endouterradX * endouterradX );
128
129 // temporary treatment!!
130 fDPhi = totphi / nseg;
131 SetFields(twistedangle, innerrad, outerrad, -halfzlen, halfzlen);
133}
134
136 G4double twistedangle,
137 G4double innerrad,
138 G4double outerrad,
139 G4double negativeEndz,
140 G4double positiveEndz,
141 G4double dphi)
142 : G4VSolid(pname), fDPhi(dphi),
143 fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
144 fFormerTwisted(0), fInnerHype(0), fOuterHype(0)
145{
146 if (innerrad < DBL_MIN)
147 {
148 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
149 FatalErrorInArgument, "Invalid end-inner-radius!");
150 }
151
152 SetFields(twistedangle, innerrad, outerrad, negativeEndz, positiveEndz);
154}
155
157 G4double twistedangle,
158 G4double innerrad,
159 G4double outerrad,
160 G4double negativeEndz,
161 G4double positiveEndz,
162 G4int nseg,
163 G4double totphi)
164 : G4VSolid(pname),
165 fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0),
166 fFormerTwisted(0), fInnerHype(0), fOuterHype(0)
167{
168 if (!nseg)
169 {
170 std::ostringstream message;
171 message << "Invalid number of segments." << G4endl
172 << " nseg = " << nseg;
173 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
174 FatalErrorInArgument, message);
175 }
176 if (totphi == DBL_MIN || innerrad < DBL_MIN)
177 {
178 G4Exception("G4TwistedTubs::G4TwistedTubs()", "GeomSolids0002",
179 FatalErrorInArgument, "Invalid total-phi or end-inner-radius!");
180 }
181
182 fDPhi = totphi / nseg;
183 SetFields(twistedangle, innerrad, outerrad, negativeEndz, positiveEndz);
185}
186
187//=====================================================================
188//* Fake default constructor ------------------------------------------
189
191 : G4VSolid(a), fPhiTwist(0.), fInnerRadius(0.), fOuterRadius(0.), fDPhi(0.),
192 fZHalfLength(0.), fInnerStereo(0.), fOuterStereo(0.), fTanInnerStereo(0.),
193 fTanOuterStereo(0.), fKappa(0.), fInnerRadius2(0.), fOuterRadius2(0.),
194 fTanInnerStereo2(0.), fTanOuterStereo2(0.), fLowerEndcap(0), fUpperEndcap(0),
195 fLatterTwisted(0), fFormerTwisted(0), fInnerHype(0), fOuterHype(0)
196{
197}
198
199//=====================================================================
200//* destructor --------------------------------------------------------
201
203{
204 if (fLowerEndcap) { delete fLowerEndcap; }
205 if (fUpperEndcap) { delete fUpperEndcap; }
206 if (fLatterTwisted) { delete fLatterTwisted; }
207 if (fFormerTwisted) { delete fFormerTwisted; }
208 if (fInnerHype) { delete fInnerHype; }
209 if (fOuterHype) { delete fOuterHype; }
210 if (fpPolyhedron) { delete fpPolyhedron; fpPolyhedron = nullptr; }
211}
212
213//=====================================================================
214//* Copy constructor --------------------------------------------------
215
217 : G4VSolid(rhs), fPhiTwist(rhs.fPhiTwist),
218 fInnerRadius(rhs.fInnerRadius), fOuterRadius(rhs.fOuterRadius),
219 fDPhi(rhs.fDPhi), fZHalfLength(rhs.fZHalfLength),
220 fInnerStereo(rhs.fInnerStereo), fOuterStereo(rhs.fOuterStereo),
221 fTanInnerStereo(rhs.fTanInnerStereo), fTanOuterStereo(rhs.fTanOuterStereo),
222 fKappa(rhs.fKappa), fInnerRadius2(rhs.fInnerRadius2),
223 fOuterRadius2(rhs.fOuterRadius2), fTanInnerStereo2(rhs.fTanInnerStereo2),
224 fTanOuterStereo2(rhs.fTanOuterStereo2),
225 fLowerEndcap(0), fUpperEndcap(0), fLatterTwisted(0), fFormerTwisted(0),
226 fInnerHype(0), fOuterHype(0),
227 fCubicVolume(rhs.fCubicVolume), fSurfaceArea(rhs.fSurfaceArea),
228 fLastInside(rhs.fLastInside), fLastNormal(rhs.fLastNormal),
229 fLastDistanceToIn(rhs.fLastDistanceToIn),
230 fLastDistanceToOut(rhs.fLastDistanceToOut),
231 fLastDistanceToInWithV(rhs.fLastDistanceToInWithV),
232 fLastDistanceToOutWithV(rhs.fLastDistanceToOutWithV)
233{
234 for (auto i=0; i<2; ++i)
235 {
236 fEndZ[i] = rhs.fEndZ[i];
239 fEndPhi[i] = rhs.fEndPhi[i];
240 fEndZ2[i] = rhs.fEndZ2[i];
241 }
243}
244
245
246//=====================================================================
247//* Assignment operator -----------------------------------------------
248
250{
251 // Check assignment to self
252 //
253 if (this == &rhs) { return *this; }
254
255 // Copy base class data
256 //
258
259 // Copy data
260 //
261 fPhiTwist= rhs.fPhiTwist;
277
278 for (auto i=0; i<2; ++i)
279 {
280 fEndZ[i] = rhs.fEndZ[i];
283 fEndPhi[i] = rhs.fEndPhi[i];
284 fEndZ2[i] = rhs.fEndZ2[i];
285 }
286
288 fRebuildPolyhedron = false;
289 delete fpPolyhedron; fpPolyhedron = nullptr;
290
291 return *this;
292}
293
294//=====================================================================
295//* ComputeDimensions -------------------------------------------------
296
298 const G4int /* n */ ,
299 const G4VPhysicalVolume* /* pRep */ )
300{
301 G4Exception("G4TwistedTubs::ComputeDimensions()",
302 "GeomSolids0001", FatalException,
303 "G4TwistedTubs does not support Parameterisation.");
304}
305
306//=====================================================================
307//* BoundingLimits ----------------------------------------------------
308
310 G4ThreeVector& pMax) const
311{
312 // Find bounding tube
313 G4double rmin = GetInnerRadius();
315
316 G4double zmin = std::min(GetEndZ(0), GetEndZ(1));
317 G4double zmax = std::max(GetEndZ(0), GetEndZ(1));
318
319 G4double dphi = 0.5*GetDPhi();
320 G4double sphi = std::min(GetEndPhi(0), GetEndPhi(1)) - dphi;
321 G4double ephi = std::max(GetEndPhi(0), GetEndPhi(1)) + dphi;
322 G4double totalphi = ephi - sphi;
323
324 // Find bounding box
325 if (dphi <= 0 || totalphi >= CLHEP::twopi)
326 {
327 pMin.set(-rmax,-rmax, zmin);
328 pMax.set( rmax, rmax, zmax);
329 }
330 else
331 {
332 G4TwoVector vmin,vmax;
333 G4GeomTools::DiskExtent(rmin, rmax, sphi, totalphi, vmin, vmax);
334 pMin.set(vmin.x(), vmin.y(), zmin);
335 pMax.set(vmax.x(), vmax.y(), zmax);
336 }
337
338 // Check correctness of the bounding box
339 //
340 if (pMin.x() >= pMax.x() || pMin.y() >= pMax.y() || pMin.z() >= pMax.z())
341 {
342 std::ostringstream message;
343 message << "Bad bounding box (min >= max) for solid: "
344 << GetName() << " !"
345 << "\npMin = " << pMin
346 << "\npMax = " << pMax;
347 G4Exception("G4TwistedTubs::BoundingLimits()", "GeomMgt0001",
348 JustWarning, message);
349 DumpInfo();
350 }
351}
352
353//=====================================================================
354//* CalculateExtent ---------------------------------------------------
355
356G4bool
358 const G4VoxelLimits& pVoxelLimit,
359 const G4AffineTransform& pTransform,
360 G4double& pMin, G4double& pMax) const
361{
362 G4ThreeVector bmin, bmax;
363
364 // Get bounding box
365 BoundingLimits(bmin,bmax);
366
367 // Find extent
368 G4BoundingEnvelope bbox(bmin,bmax);
369 return bbox.CalculateExtent(pAxis,pVoxelLimit,pTransform,pMin,pMax);
370}
371
372
373//=====================================================================
374//* Inside ------------------------------------------------------------
375
377{
378
379 const G4double halftol
381 // static G4int timerid = -1;
382 // G4Timer timer(timerid, "G4TwistedTubs", "Inside");
383 // timer.Start();
384
385 G4ThreeVector *tmpp;
386 EInside *tmpinside;
387 if (fLastInside.p == p)
388 {
389 return fLastInside.inside;
390 }
391 else
392 {
393 tmpp = const_cast<G4ThreeVector*>(&(fLastInside.p));
394 tmpinside = const_cast<EInside*>(&(fLastInside.inside));
395 tmpp->set(p.x(), p.y(), p.z());
396 }
397
398 EInside outerhypearea = ((G4TwistTubsHypeSide *)fOuterHype)->Inside(p);
399 G4double innerhyperho = ((G4TwistTubsHypeSide *)fInnerHype)->GetRhoAtPZ(p);
400 G4double distanceToOut = p.getRho() - innerhyperho; // +ve: inside
401
402 if ((outerhypearea == kOutside) || (distanceToOut < -halftol))
403 {
404 *tmpinside = kOutside;
405 }
406 else if (outerhypearea == kSurface)
407 {
408 *tmpinside = kSurface;
409 }
410 else
411 {
412 if (distanceToOut <= halftol)
413 {
414 *tmpinside = kSurface;
415 }
416 else
417 {
418 *tmpinside = kInside;
419 }
420 }
421
422 return fLastInside.inside;
423}
424
425//=====================================================================
426//* SurfaceNormal -----------------------------------------------------
427
429{
430 //
431 // return the normal unit vector to the Hyperbolical Surface at a point
432 // p on (or nearly on) the surface
433 //
434 // Which of the three or four surfaces are we closest to?
435 //
436
437 if (fLastNormal.p == p)
438 {
439 return fLastNormal.vec;
440 }
441 G4ThreeVector *tmpp =
442 const_cast<G4ThreeVector*>(&(fLastNormal.p));
443 G4ThreeVector *tmpnormal =
444 const_cast<G4ThreeVector*>(&(fLastNormal.vec));
445 G4VTwistSurface **tmpsurface =
446 const_cast<G4VTwistSurface**>(fLastNormal.surface);
447 tmpp->set(p.x(), p.y(), p.z());
448
449 G4double distance = kInfinity;
450
451 G4VTwistSurface *surfaces[6];
452 surfaces[0] = fLatterTwisted;
453 surfaces[1] = fFormerTwisted;
454 surfaces[2] = fInnerHype;
455 surfaces[3] = fOuterHype;
456 surfaces[4] = fLowerEndcap;
457 surfaces[5] = fUpperEndcap;
458
459 G4ThreeVector xx;
460 G4ThreeVector bestxx;
461 G4int besti = -1;
462 for (auto i=0; i<6; ++i)
463 {
464 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
465 if (tmpdistance < distance)
466 {
467 distance = tmpdistance;
468 bestxx = xx;
469 besti = i;
470 }
471 }
472
473 tmpsurface[0] = surfaces[besti];
474 *tmpnormal = tmpsurface[0]->GetNormal(bestxx, true);
475
476 return fLastNormal.vec;
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
522 {
523 if (currentside == kSurface)
524 {
525 // particle is just on a boundary.
526 // If the particle is entering to the volume, return 0.
527 //
529 if (normal*v < 0)
530 {
531 *tmpdist = 0.;
533 }
534 }
535 }
536
537 // now, we can take smallest positive distance.
538
539 // Initialize
540 //
541 G4double distance = kInfinity;
542
543 // find intersections and choose nearest one.
544 //
545 G4VTwistSurface* surfaces[6];
546 surfaces[0] = fLowerEndcap;
547 surfaces[1] = fUpperEndcap;
548 surfaces[2] = fLatterTwisted;
549 surfaces[3] = fFormerTwisted;
550 surfaces[4] = fInnerHype;
551 surfaces[5] = fOuterHype;
552
553 G4ThreeVector xx;
554 G4ThreeVector bestxx;
555 for (auto i=0; i<6; ++i)
556 {
557 G4double tmpdistance = surfaces[i]->DistanceToIn(p, v, xx);
558 if (tmpdistance < distance)
559 {
560 distance = tmpdistance;
561 bestxx = xx;
562 }
563 }
564 *tmpdist = distance;
565
567}
568
569//=====================================================================
570//* DistanceToIn (p) --------------------------------------------------
571
573{
574 // DistanceToIn(p):
575 // Calculate distance to surface of shape from `outside',
576 // allowing for tolerance
577
578 //
579 // checking last value
580 //
581
582 G4ThreeVector* tmpp;
583 G4double* tmpdist;
584 if (fLastDistanceToIn.p == p)
585 {
587 }
588 else
589 {
590 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToIn.p));
591 tmpdist = const_cast<G4double*>(&(fLastDistanceToIn.value));
592 tmpp->set(p.x(), p.y(), p.z());
593 }
594
595 //
596 // Calculate DistanceToIn(p)
597 //
598
599 EInside currentside = Inside(p);
600
601 switch (currentside)
602 {
603 case (kInside) :
604 {}
605 case (kSurface) :
606 {
607 *tmpdist = 0.;
609 }
610 case (kOutside) :
611 {
612 // Initialize
613 G4double distance = kInfinity;
614
615 // find intersections and choose nearest one.
616 G4VTwistSurface *surfaces[6];
617 surfaces[0] = fLowerEndcap;
618 surfaces[1] = fUpperEndcap;
619 surfaces[2] = fLatterTwisted;
620 surfaces[3] = fFormerTwisted;
621 surfaces[4] = fInnerHype;
622 surfaces[5] = fOuterHype;
623
624 G4ThreeVector xx;
625 G4ThreeVector bestxx;
626 for (auto i=0; i<6; ++i)
627 {
628 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
629 if (tmpdistance < distance)
630 {
631 distance = tmpdistance;
632 bestxx = xx;
633 }
634 }
635 *tmpdist = distance;
637 }
638 default :
639 {
640 G4Exception("G4TwistedTubs::DistanceToIn(p)", "GeomSolids0003",
641 FatalException, "Unknown point location!");
642 }
643 } // switch end
644
645 return kInfinity;
646}
647
648//=====================================================================
649//* DistanceToOut (p, v) ----------------------------------------------
650
652 const G4ThreeVector& v,
653 const G4bool calcNorm,
654 G4bool *validNorm,
655 G4ThreeVector *norm ) const
656{
657 // DistanceToOut (p, v):
658 // Calculate distance to surface of shape from `inside'
659 // along with the v, allowing for tolerance.
660 // The function returns kInfinity if no intersection or
661 // just grazing within tolerance.
662
663 //
664 // checking last value
665 //
666
667 G4ThreeVector* tmpp;
668 G4ThreeVector* tmpv;
669 G4double* tmpdist;
671 {
673 }
674 else
675 {
676 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.p));
677 tmpv = const_cast<G4ThreeVector*>(&(fLastDistanceToOutWithV.vec));
678 tmpdist = const_cast<G4double*>(&(fLastDistanceToOutWithV.value));
679 tmpp->set(p.x(), p.y(), p.z());
680 tmpv->set(v.x(), v.y(), v.z());
681 }
682
683 //
684 // Calculate DistanceToOut(p,v)
685 //
686
687 EInside currentside = Inside(p);
688
689 if (currentside == kOutside)
690 {
691 }
692 else
693 {
694 if (currentside == kSurface)
695 {
696 // particle is just on a boundary.
697 // If the particle is exiting from the volume, return 0.
698 //
700 G4VTwistSurface *blockedsurface = fLastNormal.surface[0];
701 if (normal*v > 0)
702 {
703 if (calcNorm)
704 {
705 *norm = (blockedsurface->GetNormal(p, true));
706 *validNorm = blockedsurface->IsValidNorm();
707 }
708 *tmpdist = 0.;
710 }
711 }
712 }
713
714 // now, we can take smallest positive distance.
715
716 // Initialize
717 //
718 G4double distance = kInfinity;
719
720 // find intersections and choose nearest one.
721 //
722 G4VTwistSurface* surfaces[6];
723 surfaces[0] = fLatterTwisted;
724 surfaces[1] = fFormerTwisted;
725 surfaces[2] = fInnerHype;
726 surfaces[3] = fOuterHype;
727 surfaces[4] = fLowerEndcap;
728 surfaces[5] = fUpperEndcap;
729
730 G4int besti = -1;
731 G4ThreeVector xx;
732 G4ThreeVector bestxx;
733 for (auto i=0; i<6; ++i)
734 {
735 G4double tmpdistance = surfaces[i]->DistanceToOut(p, v, xx);
736 if (tmpdistance < distance)
737 {
738 distance = tmpdistance;
739 bestxx = xx;
740 besti = i;
741 }
742 }
743
744 if (calcNorm)
745 {
746 if (besti != -1)
747 {
748 *norm = (surfaces[besti]->GetNormal(p, true));
749 *validNorm = surfaces[besti]->IsValidNorm();
750 }
751 }
752
753 *tmpdist = distance;
754
756}
757
758
759//=====================================================================
760//* DistanceToOut (p) ----------------------------------------------
761
763{
764 // DistanceToOut(p):
765 // Calculate distance to surface of shape from `inside',
766 // allowing for tolerance
767
768 //
769 // checking last value
770 //
771
772 G4ThreeVector* tmpp;
773 G4double* tmpdist;
774 if (fLastDistanceToOut.p == p)
775 {
777 }
778 else
779 {
780 tmpp = const_cast<G4ThreeVector*>(&(fLastDistanceToOut.p));
781 tmpdist = const_cast<G4double*>(&(fLastDistanceToOut.value));
782 tmpp->set(p.x(), p.y(), p.z());
783 }
784
785 //
786 // Calculate DistanceToOut(p)
787 //
788
789 EInside currentside = Inside(p);
790
791 switch (currentside)
792 {
793 case (kOutside) :
794 {
795 }
796 case (kSurface) :
797 {
798 *tmpdist = 0.;
800 }
801 case (kInside) :
802 {
803 // Initialize
804 G4double distance = kInfinity;
805
806 // find intersections and choose nearest one.
807 G4VTwistSurface* surfaces[6];
808 surfaces[0] = fLatterTwisted;
809 surfaces[1] = fFormerTwisted;
810 surfaces[2] = fInnerHype;
811 surfaces[3] = fOuterHype;
812 surfaces[4] = fLowerEndcap;
813 surfaces[5] = fUpperEndcap;
814
815 G4ThreeVector xx;
816 G4ThreeVector bestxx;
817 for (auto i=0; i<6; ++i)
818 {
819 G4double tmpdistance = surfaces[i]->DistanceTo(p, xx);
820 if (tmpdistance < distance)
821 {
822 distance = tmpdistance;
823 bestxx = xx;
824 }
825 }
826 *tmpdist = distance;
827
829 }
830 default :
831 {
832 G4Exception("G4TwistedTubs::DistanceToOut(p)", "GeomSolids0003",
833 FatalException, "Unknown point location!");
834 }
835 } // switch end
836
837 return 0.;
838}
839
840//=====================================================================
841//* StreamInfo --------------------------------------------------------
842
843std::ostream& G4TwistedTubs::StreamInfo(std::ostream& os) const
844{
845 //
846 // Stream object contents to an output stream
847 //
848 G4int oldprc = os.precision(16);
849 os << "-----------------------------------------------------------\n"
850 << " *** Dump for solid - " << GetName() << " ***\n"
851 << " ===================================================\n"
852 << " Solid type: G4TwistedTubs\n"
853 << " Parameters: \n"
854 << " -ve end Z : " << fEndZ[0]/mm << " mm \n"
855 << " +ve end Z : " << fEndZ[1]/mm << " mm \n"
856 << " inner end radius(-ve z): " << fEndInnerRadius[0]/mm << " mm \n"
857 << " inner end radius(+ve z): " << fEndInnerRadius[1]/mm << " mm \n"
858 << " outer end radius(-ve z): " << fEndOuterRadius[0]/mm << " mm \n"
859 << " outer end radius(+ve z): " << fEndOuterRadius[1]/mm << " mm \n"
860 << " inner radius (z=0) : " << fInnerRadius/mm << " mm \n"
861 << " outer radius (z=0) : " << fOuterRadius/mm << " mm \n"
862 << " twisted angle : " << fPhiTwist/degree << " degrees \n"
863 << " inner stereo angle : " << fInnerStereo/degree << " degrees \n"
864 << " outer stereo angle : " << fOuterStereo/degree << " degrees \n"
865 << " phi-width of a piece : " << fDPhi/degree << " degrees \n"
866 << "-----------------------------------------------------------\n";
867 os.precision(oldprc);
868
869 return os;
870}
871
872
873//=====================================================================
874//* DiscribeYourselfTo ------------------------------------------------
875
877{
878 scene.AddSolid (*this);
879}
880
881//=====================================================================
882//* GetExtent ---------------------------------------------------------
883
885{
886 // Define the sides of the box into which the G4Tubs instance would fit.
887 //
888 G4ThreeVector pmin,pmax;
889 BoundingLimits(pmin,pmax);
890 return G4VisExtent(pmin.x(),pmax.x(),
891 pmin.y(),pmax.y(),
892 pmin.z(),pmax.z());
893}
894
895//=====================================================================
896//* CreatePolyhedron --------------------------------------------------
897
899{
900 // number of meshes
901 //
902 G4double absPhiTwist = std::abs(fPhiTwist);
903 G4double dA = std::max(fDPhi,absPhiTwist);
904 const G4int k =
906 const G4int n =
908
909 const G4int nnodes = 4*(k-1)*(n-2) + 2*k*k ;
910 const G4int nfaces = 4*(k-1)*(n-1) + 2*(k-1)*(k-1) ;
911
912 G4Polyhedron* ph = new G4Polyhedron;
913 typedef G4double G4double3[3];
914 typedef G4int G4int4[4];
915 G4double3* xyz = new G4double3[nnodes]; // number of nodes
916 G4int4* faces = new G4int4[nfaces] ; // number of faces
917 fLowerEndcap->GetFacets(k,k,xyz,faces,0) ;
918 fUpperEndcap->GetFacets(k,k,xyz,faces,1) ;
919 fInnerHype->GetFacets(k,n,xyz,faces,2) ;
920 fFormerTwisted->GetFacets(k,n,xyz,faces,3) ;
921 fOuterHype->GetFacets(k,n,xyz,faces,4) ;
922 fLatterTwisted->GetFacets(k,n,xyz,faces,5) ;
923
924 ph->createPolyhedron(nnodes,nfaces,xyz,faces);
925
926 delete[] xyz;
927 delete[] faces;
928
929 return ph;
930}
931
932//=====================================================================
933//* GetPolyhedron -----------------------------------------------------
934
936{
937 if (fpPolyhedron == nullptr ||
941 {
943 delete fpPolyhedron;
945 fRebuildPolyhedron = false;
946 l.unlock();
947 }
948 return fpPolyhedron;
949}
950
951//=====================================================================
952//* CreateSurfaces ----------------------------------------------------
953
955{
956 // create 6 surfaces of TwistedTub
957
958 G4ThreeVector x0(0, 0, fEndZ[0]);
959 G4ThreeVector n (0, 0, -1);
960
961 fLowerEndcap = new G4TwistTubsFlatSide("LowerEndcap",
963 fDPhi, fEndPhi, fEndZ, -1) ;
964
965 fUpperEndcap = new G4TwistTubsFlatSide("UpperEndcap",
967 fDPhi, fEndPhi, fEndZ, 1) ;
968
969 G4RotationMatrix rotHalfDPhi;
970 rotHalfDPhi.rotateZ(0.5*fDPhi);
971
972 fLatterTwisted = new G4TwistTubsSide("LatterTwisted",
976 1 ) ;
977 fFormerTwisted = new G4TwistTubsSide("FormerTwisted",
981 -1 ) ;
982
983 fInnerHype = new G4TwistTubsHypeSide("InnerHype",
988 fOuterHype = new G4TwistTubsHypeSide("OuterHype",
993
994
995 // set neighbour surfaces
996 //
1009}
1010
1011
1012//=====================================================================
1013//* GetEntityType -----------------------------------------------------
1014
1016{
1017 return G4String("G4TwistedTubs");
1018}
1019
1020//=====================================================================
1021//* Clone -------------------------------------------------------------
1022
1024{
1025 return new G4TwistedTubs(*this);
1026}
1027
1028//=====================================================================
1029//* GetCubicVolume ----------------------------------------------------
1030
1032{
1033 if (fCubicVolume == 0.)
1034 {
1035 G4double DPhi = GetDPhi();
1036 G4double Z0 = GetEndZ(0);
1037 G4double Z1 = GetEndZ(1);
1038 G4double Ain = GetInnerRadius();
1039 G4double Aout = GetOuterRadius();
1040 G4double R0in = GetEndInnerRadius(0);
1041 G4double R1in = GetEndInnerRadius(1);
1042 G4double R0out = GetEndOuterRadius(0);
1043 G4double R1out = GetEndOuterRadius(1);
1044
1045 // V_hyperboloid = pi*h*(2*a*a + R*R)/3
1046 fCubicVolume = (2.*(Z1 - Z0)*(Aout + Ain)*(Aout - Ain)
1047 + Z1*(R1out + R1in)*(R1out - R1in)
1048 - Z0*(R0out + R0in)*(R0out - R0in))*DPhi/6.;
1049 }
1050 return fCubicVolume;
1051}
1052
1053//=====================================================================
1054//* GetSurfaceArea ----------------------------------------------------
1055
1057{
1059 return fSurfaceArea;
1060}
1061
1062//=====================================================================
1063//* GetPointOnSurface -------------------------------------------------
1064
1066{
1067
1069 G4double phi , phimin, phimax ;
1070 G4double x , xmin, xmax ;
1071 G4double r , rmin, rmax ;
1072
1079
1080 G4double chose = G4RandFlat::shoot(0.,a1 + a2 + a3 + a4 + a5 + a6) ;
1081
1082 if(chose < a1)
1083 {
1084
1085 phimin = fOuterHype->GetBoundaryMin(z) ;
1086 phimax = fOuterHype->GetBoundaryMax(z) ;
1087 phi = G4RandFlat::shoot(phimin,phimax) ;
1088
1089 return fOuterHype->SurfacePoint(phi,z,true) ;
1090
1091 }
1092 else if ( (chose >= a1) && (chose < a1 + a2 ) )
1093 {
1094
1095 phimin = fInnerHype->GetBoundaryMin(z) ;
1096 phimax = fInnerHype->GetBoundaryMax(z) ;
1097 phi = G4RandFlat::shoot(phimin,phimax) ;
1098
1099 return fInnerHype->SurfacePoint(phi,z,true) ;
1100
1101 }
1102 else if ( (chose >= a1 + a2 ) && (chose < a1 + a2 + a3 ) )
1103 {
1104
1105 xmin = fLatterTwisted->GetBoundaryMin(z) ;
1106 xmax = fLatterTwisted->GetBoundaryMax(z) ;
1107 x = G4RandFlat::shoot(xmin,xmax) ;
1108
1109 return fLatterTwisted->SurfacePoint(x,z,true) ;
1110
1111 }
1112 else if ( (chose >= a1 + a2 + a3 ) && (chose < a1 + a2 + a3 + a4 ) )
1113 {
1114
1115 xmin = fFormerTwisted->GetBoundaryMin(z) ;
1116 xmax = fFormerTwisted->GetBoundaryMax(z) ;
1117 x = G4RandFlat::shoot(xmin,xmax) ;
1118
1119 return fFormerTwisted->SurfacePoint(x,z,true) ;
1120 }
1121 else if( (chose >= a1 + a2 + a3 + a4 )&&(chose < a1 + a2 + a3 + a4 + a5 ) )
1122 {
1123 rmin = GetEndInnerRadius(0) ;
1124 rmax = GetEndOuterRadius(0) ;
1125 r = std::sqrt(G4RandFlat::shoot()*(sqr(rmax)-sqr(rmin))+sqr(rmin));
1126
1127 phimin = fLowerEndcap->GetBoundaryMin(r) ;
1128 phimax = fLowerEndcap->GetBoundaryMax(r) ;
1129 phi = G4RandFlat::shoot(phimin,phimax) ;
1130
1131 return fLowerEndcap->SurfacePoint(phi,r,true) ;
1132 }
1133 else
1134 {
1135 rmin = GetEndInnerRadius(1) ;
1136 rmax = GetEndOuterRadius(1) ;
1137 r = rmin + (rmax-rmin)*std::sqrt(G4RandFlat::shoot());
1138
1139 phimin = fUpperEndcap->GetBoundaryMin(r) ;
1140 phimax = fUpperEndcap->GetBoundaryMax(r) ;
1141 phi = G4RandFlat::shoot(phimin,phimax) ;
1142
1143 return fUpperEndcap->SurfacePoint(phi,r,true) ;
1144 }
1145}
@ 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
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
double x() const
double y() const
double z() const
double x() const
double y() const
double getRho() const
void set(double x, double y, double z)
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimits, const G4Transform3D &pTransform3D, G4double &pMin, G4double &pMax) const
static G4bool DiskExtent(G4double rmin, G4double rmax, G4double startPhi, G4double delPhi, G4TwoVector &pmin, G4TwoVector &pmax)
Definition: G4GeomTools.cc:390
G4double GetRadialTolerance() const
static G4GeometryTolerance * GetInstance()
G4int GetNumberOfRotationStepsAtTimeOfCreation() const
G4VTwistSurface ** surface
G4double fPhiTwist
G4double fOuterRadius
G4double fEndZ2[2]
G4double GetOuterRadius() const
G4double fTanInnerStereo2
G4double fOuterStereo
LastValue fLastDistanceToOut
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
G4double fSurfaceArea
G4bool CalculateExtent(const EAxis pAxis, const G4VoxelLimits &pVoxelLimit, const G4AffineTransform &pTransform, G4double &pMin, G4double &pMax) const
G4ThreeVector GetPointOnSurface() const
G4double fTanInnerStereo
LastValue fLastDistanceToIn
G4double fEndZ[2]
G4double GetSurfaceArea()
G4VTwistSurface * fLatterTwisted
G4double fTanOuterStereo
void DescribeYourselfTo(G4VGraphicsScene &scene) const
G4double fCubicVolume
G4Polyhedron * CreatePolyhedron() const
void SetFields(G4double phitwist, G4double innerrad, G4double outerrad, G4double negativeEndz, G4double positiveEndz)
void ComputeDimensions(G4VPVParameterisation *, const G4int, const G4VPhysicalVolume *)
G4VTwistSurface * fFormerTwisted
G4VSolid * Clone() const
G4TwistedTubs & operator=(const G4TwistedTubs &rhs)
G4double GetEndInnerRadius() const
G4double GetEndOuterRadius() const
G4VTwistSurface * fUpperEndcap
void CreateSurfaces()
G4double fEndInnerRadius[2]
G4double fInnerStereo
G4double GetCubicVolume()
G4Polyhedron * GetPolyhedron() const
G4double fTanOuterStereo2
G4VTwistSurface * fLowerEndcap
G4double fInnerRadius2
G4GeometryType GetEntityType() const
virtual ~G4TwistedTubs()
LastValueWithDoubleVector fLastDistanceToOutWithV
G4double GetEndPhi(G4int i) const
G4TwistedTubs(const G4String &pname, G4double twistedangle, G4double endinnerrad, G4double endouterrad, G4double halfzlen, G4double dphi)
LastVector fLastNormal
LastValueWithDoubleVector fLastDistanceToInWithV
G4VTwistSurface * fInnerHype
G4double fInnerRadius
EInside Inside(const G4ThreeVector &p) const
G4double fZHalfLength
G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcnorm=false, G4bool *validnorm=nullptr, G4ThreeVector *n=nullptr) const
G4VisExtent GetExtent() const
G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const
G4bool fRebuildPolyhedron
G4double GetEndZ(G4int i) const
G4Polyhedron * fpPolyhedron
G4double fOuterRadius2
LastState fLastInside
G4double GetInnerRadius() const
void BoundingLimits(G4ThreeVector &pMin, G4ThreeVector &pMax) const
G4double fEndPhi[2]
G4VTwistSurface * fOuterHype
G4double fEndOuterRadius[2]
std::ostream & StreamInfo(std::ostream &os) const
G4double GetDPhi() const
virtual void AddSolid(const G4Box &)=0
G4String GetName() const
void DumpInfo() const
G4VSolid & operator=(const G4VSolid &rhs)
Definition: G4VSolid.cc:107
virtual G4double GetSurfaceArea()
Definition: G4VSolid.cc:250
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
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 constexpr double twopi
Definition: SystemOfUnits.h:56
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:79
ThreeVector shoot(const G4int Ap, const G4int Af)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static const G4double Z1[5]
Definition: paraMaker.cc:41
string pname
Definition: eplot.py:33
T sqr(const T &x)
Definition: templates.hh:128
#define DBL_MIN
Definition: templates.hh:54