Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4Polycone.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27// $Id$
28//
29//
30// --------------------------------------------------------------------
31// GEANT 4 class source file
32//
33//
34// G4Polycone.cc
35//
36// Implementation of a CSG polycone
37//
38// --------------------------------------------------------------------
39
40#include "G4Polycone.hh"
41
42#include "G4PolyconeSide.hh"
43#include "G4PolyPhiFace.hh"
44
45#include "Randomize.hh"
46
47#include "G4Polyhedron.hh"
49#include "G4ReduciblePolygon.hh"
51
52using namespace CLHEP;
53
54//
55// Constructor (GEANT3 style parameters)
56//
58 G4double phiStart,
59 G4double phiTotal,
60 G4int numZPlanes,
61 const G4double zPlane[],
62 const G4double rInner[],
63 const G4double rOuter[] )
64 : G4VCSGfaceted( name ), genericPcon(false)
65{
66 //
67 // Some historical ugliness
68 //
70
74 original_parameters->Z_values = new G4double[numZPlanes];
75 original_parameters->Rmin = new G4double[numZPlanes];
76 original_parameters->Rmax = new G4double[numZPlanes];
77
78 G4int i;
79 for (i=0; i<numZPlanes; i++)
80 {
81 if (( i < numZPlanes-1) && ( zPlane[i] == zPlane[i+1] ))
82 {
83 if( (rInner[i] > rOuter[i+1])
84 ||(rInner[i+1] > rOuter[i]) )
85 {
86 DumpInfo();
87 std::ostringstream message;
88 message << "Cannot create a Polycone with no contiguous segments."
89 << G4endl
90 << " Segments are not contiguous !" << G4endl
91 << " rMin[" << i << "] = " << rInner[i]
92 << " -- rMax[" << i+1 << "] = " << rOuter[i+1] << G4endl
93 << " rMin[" << i+1 << "] = " << rInner[i+1]
94 << " -- rMax[" << i << "] = " << rOuter[i];
95 G4Exception("G4Polycone::G4Polycone()", "GeomSolids0002",
96 FatalErrorInArgument, message);
97 }
98 }
99 original_parameters->Z_values[i] = zPlane[i];
100 original_parameters->Rmin[i] = rInner[i];
101 original_parameters->Rmax[i] = rOuter[i];
102 }
103
104 //
105 // Build RZ polygon using special PCON/PGON GEANT3 constructor
106 //
108 new G4ReduciblePolygon( rInner, rOuter, zPlane, numZPlanes );
109
110 //
111 // Do the real work
112 //
113 Create( phiStart, phiTotal, rz );
114
115 delete rz;
116}
117
118
119//
120// Constructor (generic parameters)
121//
123 G4double phiStart,
124 G4double phiTotal,
125 G4int numRZ,
126 const G4double r[],
127 const G4double z[] )
128 : G4VCSGfaceted( name ), genericPcon(true)
129{
130 G4ReduciblePolygon *rz = new G4ReduciblePolygon( r, z, numRZ );
131
132 Create( phiStart, phiTotal, rz );
133
134 // Set original_parameters struct for consistency
135 //
137
138 delete rz;
139}
140
141
142//
143// Create
144//
145// Generic create routine, called by each constructor after
146// conversion of arguments
147//
149 G4double phiTotal,
151{
152 //
153 // Perform checks of rz values
154 //
155 if (rz->Amin() < 0.0)
156 {
157 std::ostringstream message;
158 message << "Illegal input parameters - " << GetName() << G4endl
159 << " All R values must be >= 0 !";
160 G4Exception("G4Polycone::Create()", "GeomSolids0002",
161 FatalErrorInArgument, message);
162 }
163
164 G4double rzArea = rz->Area();
165 if (rzArea < -kCarTolerance)
166 rz->ReverseOrder();
167
168 else if (rzArea < -kCarTolerance)
169 {
170 std::ostringstream message;
171 message << "Illegal input parameters - " << GetName() << G4endl
172 << " R/Z cross section is zero or near zero: " << rzArea;
173 G4Exception("G4Polycone::Create()", "GeomSolids0002",
174 FatalErrorInArgument, message);
175 }
176
179 {
180 std::ostringstream message;
181 message << "Illegal input parameters - " << GetName() << G4endl
182 << " Too few unique R/Z values !";
183 G4Exception("G4Polycone::Create()", "GeomSolids0002",
184 FatalErrorInArgument, message);
185 }
186
187 if (rz->CrossesItself(1/kInfinity))
188 {
189 std::ostringstream message;
190 message << "Illegal input parameters - " << GetName() << G4endl
191 << " R/Z segments cross !";
192 G4Exception("G4Polycone::Create()", "GeomSolids0002",
193 FatalErrorInArgument, message);
194 }
195
196 numCorner = rz->NumVertices();
197
198 //
199 // Phi opening? Account for some possible roundoff, and interpret
200 // nonsense value as representing no phi opening
201 //
202 if (phiTotal <= 0 || phiTotal > twopi-1E-10)
203 {
204 phiIsOpen = false;
205 startPhi = 0;
206 endPhi = twopi;
207 }
208 else
209 {
210 phiIsOpen = true;
211
212 //
213 // Convert phi into our convention
214 //
215 startPhi = phiStart;
216 while( startPhi < 0 ) startPhi += twopi;
217
218 endPhi = phiStart+phiTotal;
219 while( endPhi < startPhi ) endPhi += twopi;
220 }
221
222 //
223 // Allocate corner array.
224 //
226
227 //
228 // Copy corners
229 //
231
233 iterRZ.Begin();
234 do
235 {
236 next->r = iterRZ.GetA();
237 next->z = iterRZ.GetB();
238 } while( ++next, iterRZ.Next() );
239
240 //
241 // Allocate face pointer array
242 //
244 faces = new G4VCSGface*[numFace];
245
246 //
247 // Construct conical faces
248 //
249 // But! Don't construct a face if both points are at zero radius!
250 //
251 G4PolyconeSideRZ *corner = corners,
252 *prev = corners + numCorner-1,
253 *nextNext;
254 G4VCSGface **face = faces;
255 do
256 {
257 next = corner+1;
258 if (next >= corners+numCorner) next = corners;
259 nextNext = next+1;
260 if (nextNext >= corners+numCorner) nextNext = corners;
261
262 if (corner->r < 1/kInfinity && next->r < 1/kInfinity) continue;
263
264 //
265 // We must decide here if we can dare declare one of our faces
266 // as having a "valid" normal (i.e. allBehind = true). This
267 // is never possible if the face faces "inward" in r.
268 //
269 G4bool allBehind;
270 if (corner->z > next->z)
271 {
272 allBehind = false;
273 }
274 else
275 {
276 //
277 // Otherwise, it is only true if the line passing
278 // through the two points of the segment do not
279 // split the r/z cross section
280 //
281 allBehind = !rz->BisectedBy( corner->r, corner->z,
282 next->r, next->z, kCarTolerance );
283 }
284
285 *face++ = new G4PolyconeSide( prev, corner, next, nextNext,
286 startPhi, endPhi-startPhi, phiIsOpen, allBehind );
287 } while( prev=corner, corner=next, corner > corners );
288
289 if (phiIsOpen)
290 {
291 //
292 // Construct phi open edges
293 //
294 *face++ = new G4PolyPhiFace( rz, startPhi, 0, endPhi );
295 *face++ = new G4PolyPhiFace( rz, endPhi, 0, startPhi );
296 }
297
298 //
299 // We might have dropped a face or two: recalculate numFace
300 //
301 numFace = face-faces;
302
303 //
304 // Make enclosingCylinder
305 //
307 new G4EnclosingCylinder( rz, phiIsOpen, phiStart, phiTotal );
308}
309
310
311//
312// Fake default constructor - sets only member data and allocates memory
313// for usage restricted to object persistency.
314//
316 : G4VCSGfaceted(a), startPhi(0.), endPhi(0.), phiIsOpen(false),
317 genericPcon(false), numCorner(0), corners(0),
318 original_parameters(0), enclosingCylinder(0)
319{
320}
321
322
323//
324// Destructor
325//
327{
328 delete [] corners;
329 delete original_parameters;
330 delete enclosingCylinder;
331}
332
333
334//
335// Copy constructor
336//
338 : G4VCSGfaceted( source )
339{
340 CopyStuff( source );
341}
342
343
344//
345// Assignment operator
346//
348{
349 if (this == &source) return *this;
350
351 G4VCSGfaceted::operator=( source );
352
353 delete [] corners;
355
356 delete enclosingCylinder;
357
358 CopyStuff( source );
359
360 return *this;
361}
362
363
364//
365// CopyStuff
366//
368{
369 //
370 // Simple stuff
371 //
372 startPhi = source.startPhi;
373 endPhi = source.endPhi;
374 phiIsOpen = source.phiIsOpen;
375 numCorner = source.numCorner;
376 genericPcon= source.genericPcon;
377
378 //
379 // The corner array
380 //
382
384 *sourceCorn = source.corners;
385 do
386 {
387 *corn = *sourceCorn;
388 } while( ++sourceCorn, ++corn < corners+numCorner );
389
390 //
391 // Original parameters
392 //
393 if (source.original_parameters)
394 {
397 }
398
399 //
400 // Enclosing cylinder
401 //
403}
404
405
406//
407// Reset
408//
410{
411 if (genericPcon)
412 {
413 std::ostringstream message;
414 message << "Solid " << GetName() << " built using generic construct."
415 << G4endl << "Not applicable to the generic construct !";
416 G4Exception("G4Polycone::Reset()", "GeomSolids1001",
417 JustWarning, message, "Parameters NOT resetted.");
418 return 1;
419 }
420
421 //
422 // Clear old setup
423 //
425 delete [] corners;
426 delete enclosingCylinder;
427
428 //
429 // Rebuild polycone
430 //
438 delete rz;
439
440 return 0;
441}
442
443
444//
445// Inside
446//
447// This is an override of G4VCSGfaceted::Inside, created in order
448// to speed things up by first checking with G4EnclosingCylinder.
449//
451{
452 //
453 // Quick test
454 //
456
457 //
458 // Long answer
459 //
460 return G4VCSGfaceted::Inside(p);
461}
462
463
464//
465// DistanceToIn
466//
467// This is an override of G4VCSGfaceted::Inside, created in order
468// to speed things up by first checking with G4EnclosingCylinder.
469//
471 const G4ThreeVector &v ) const
472{
473 //
474 // Quick test
475 //
477 return kInfinity;
478
479 //
480 // Long answer
481 //
482 return G4VCSGfaceted::DistanceToIn( p, v );
483}
484
485
486//
487// DistanceToIn
488//
490{
492}
493
494
495//
496// ComputeDimensions
497//
499 const G4int n,
500 const G4VPhysicalVolume* pRep )
501{
502 p->ComputeDimensions(*this,n,pRep);
503}
504
505//
506// GetEntityType
507//
509{
510 return G4String("G4Polycone");
511}
512
513//
514// Make a clone of the object
515//
517{
518 return new G4Polycone(*this);
519}
520
521//
522// Stream object contents to an output stream
523//
524std::ostream& G4Polycone::StreamInfo( std::ostream& os ) const
525{
526 G4int oldprc = os.precision(16);
527 os << "-----------------------------------------------------------\n"
528 << " *** Dump for solid - " << GetName() << " ***\n"
529 << " ===================================================\n"
530 << " Solid type: G4Polycone\n"
531 << " Parameters: \n"
532 << " starting phi angle : " << startPhi/degree << " degrees \n"
533 << " ending phi angle : " << endPhi/degree << " degrees \n";
534 G4int i=0;
535 if (!genericPcon)
536 {
538 os << " number of Z planes: " << numPlanes << "\n"
539 << " Z values: \n";
540 for (i=0; i<numPlanes; i++)
541 {
542 os << " Z plane " << i << ": "
543 << original_parameters->Z_values[i] << "\n";
544 }
545 os << " Tangent distances to inner surface (Rmin): \n";
546 for (i=0; i<numPlanes; i++)
547 {
548 os << " Z plane " << i << ": "
549 << original_parameters->Rmin[i] << "\n";
550 }
551 os << " Tangent distances to outer surface (Rmax): \n";
552 for (i=0; i<numPlanes; i++)
553 {
554 os << " Z plane " << i << ": "
555 << original_parameters->Rmax[i] << "\n";
556 }
557 }
558 os << " number of RZ points: " << numCorner << "\n"
559 << " RZ values (corners): \n";
560 for (i=0; i<numCorner; i++)
561 {
562 os << " "
563 << corners[i].r << ", " << corners[i].z << "\n";
564 }
565 os << "-----------------------------------------------------------\n";
566 os.precision(oldprc);
567
568 return os;
569}
570
571
572//
573// GetPointOnCone
574//
575// Auxiliary method for Get Point On Surface
576//
578 G4double fRmin2, G4double fRmax2,
579 G4double zOne, G4double zTwo,
580 G4double& totArea) const
581{
582 // declare working variables
583 //
584 G4double Aone, Atwo, Afive, phi, zRand, fDPhi, cosu, sinu;
585 G4double rRand1, rmin, rmax, chose, rone, rtwo, qone, qtwo,
586 fDz = std::fabs((zTwo-zOne)/2.);
587 G4ThreeVector point, offset;
588 offset = G4ThreeVector(0.,0.,0.5*(zTwo+zOne));
589 fDPhi = endPhi - startPhi;
590 rone = (fRmax1-fRmax2)/(2.*fDz);
591 rtwo = (fRmin1-fRmin2)/(2.*fDz);
592 if(fRmax1==fRmax2){qone=0.;}
593 else{
594 qone = fDz*(fRmax1+fRmax2)/(fRmax1-fRmax2);
595 }
596 if(fRmin1==fRmin2){qtwo=0.;}
597 else{
598 qtwo = fDz*(fRmin1+fRmin2)/(fRmin1-fRmin2);
599 }
600 Aone = 0.5*fDPhi*(fRmax2 + fRmax1)*(sqr(fRmin1-fRmin2)+sqr(zTwo-zOne));
601 Atwo = 0.5*fDPhi*(fRmin2 + fRmin1)*(sqr(fRmax1-fRmax2)+sqr(zTwo-zOne));
602 Afive = fDz*(fRmax1-fRmin1+fRmax2-fRmin2);
603 totArea = Aone+Atwo+2.*Afive;
604
606 cosu = std::cos(phi);
607 sinu = std::sin(phi);
608
609
610 if( (startPhi == 0) && (endPhi == twopi) ) { Afive = 0; }
611 chose = RandFlat::shoot(0.,Aone+Atwo+2.*Afive);
612 if( (chose >= 0) && (chose < Aone) )
613 {
614 if(fRmax1 != fRmax2)
615 {
616 zRand = RandFlat::shoot(-1.*fDz,fDz);
617 point = G4ThreeVector (rone*cosu*(qone-zRand),
618 rone*sinu*(qone-zRand), zRand);
619
620
621 }
622 else
623 {
624 point = G4ThreeVector(fRmax1*cosu, fRmax1*sinu,
625 RandFlat::shoot(-1.*fDz,fDz));
626
627 }
628 }
629 else if(chose >= Aone && chose < Aone + Atwo)
630 {
631 if(fRmin1 != fRmin2)
632 {
633 zRand = RandFlat::shoot(-1.*fDz,fDz);
634 point = G4ThreeVector (rtwo*cosu*(qtwo-zRand),
635 rtwo*sinu*(qtwo-zRand), zRand);
636
637 }
638 else
639 {
640 point = G4ThreeVector(fRmin1*cosu, fRmin1*sinu,
641 RandFlat::shoot(-1.*fDz,fDz));
642
643 }
644 }
645 else if( (chose >= Aone + Atwo + Afive) && (chose < Aone + Atwo + 2.*Afive) )
646 {
647 zRand = RandFlat::shoot(-1.*fDz,fDz);
648 rmin = fRmin2-((zRand-fDz)/(2.*fDz))*(fRmin1-fRmin2);
649 rmax = fRmax2-((zRand-fDz)/(2.*fDz))*(fRmax1-fRmax2);
650 rRand1 = std::sqrt(RandFlat::shoot()*(sqr(rmax)-sqr(rmin))+sqr(rmin));
651 point = G4ThreeVector (rRand1*std::cos(startPhi),
652 rRand1*std::sin(startPhi), zRand);
653 }
654 else
655 {
656 zRand = RandFlat::shoot(-1.*fDz,fDz);
657 rmin = fRmin2-((zRand-fDz)/(2.*fDz))*(fRmin1-fRmin2);
658 rmax = fRmax2-((zRand-fDz)/(2.*fDz))*(fRmax1-fRmax2);
659 rRand1 = std::sqrt(RandFlat::shoot()*(sqr(rmax)-sqr(rmin))+sqr(rmin));
660 point = G4ThreeVector (rRand1*std::cos(endPhi),
661 rRand1*std::sin(endPhi), zRand);
662
663 }
664 return point+offset;
665}
666
667
668//
669// GetPointOnTubs
670//
671// Auxiliary method for GetPoint On Surface
672//
674 G4double zOne, G4double zTwo,
675 G4double& totArea) const
676{
677 G4double xRand,yRand,zRand,phi,cosphi,sinphi,chose,
678 aOne,aTwo,aFou,rRand,fDz,fSPhi,fDPhi;
679 fDz = std::fabs(0.5*(zTwo-zOne));
680 fSPhi = startPhi;
681 fDPhi = endPhi-startPhi;
682
683 aOne = 2.*fDz*fDPhi*fRMax;
684 aTwo = 2.*fDz*fDPhi*fRMin;
685 aFou = 2.*fDz*(fRMax-fRMin);
686 totArea = aOne+aTwo+2.*aFou;
688 cosphi = std::cos(phi);
689 sinphi = std::sin(phi);
690 rRand = fRMin + (fRMax-fRMin)*std::sqrt(RandFlat::shoot());
691
692 if(startPhi == 0 && endPhi == twopi)
693 aFou = 0;
694
695 chose = RandFlat::shoot(0.,aOne+aTwo+2.*aFou);
696 if( (chose >= 0) && (chose < aOne) )
697 {
698 xRand = fRMax*cosphi;
699 yRand = fRMax*sinphi;
700 zRand = RandFlat::shoot(-1.*fDz,fDz);
701 return G4ThreeVector(xRand, yRand, zRand+0.5*(zTwo+zOne));
702 }
703 else if( (chose >= aOne) && (chose < aOne + aTwo) )
704 {
705 xRand = fRMin*cosphi;
706 yRand = fRMin*sinphi;
707 zRand = RandFlat::shoot(-1.*fDz,fDz);
708 return G4ThreeVector(xRand, yRand, zRand+0.5*(zTwo+zOne));
709 }
710 else if( (chose >= aOne+aTwo) && (chose <aOne+aTwo+aFou) )
711 {
712 xRand = rRand*std::cos(fSPhi+fDPhi);
713 yRand = rRand*std::sin(fSPhi+fDPhi);
714 zRand = RandFlat::shoot(-1.*fDz,fDz);
715 return G4ThreeVector(xRand, yRand, zRand+0.5*(zTwo+zOne));
716 }
717
718 // else
719
720 xRand = rRand*std::cos(fSPhi+fDPhi);
721 yRand = rRand*std::sin(fSPhi+fDPhi);
722 zRand = RandFlat::shoot(-1.*fDz,fDz);
723 return G4ThreeVector(xRand, yRand, zRand+0.5*(zTwo+zOne));
724}
725
726
727//
728// GetPointOnRing
729//
730// Auxiliary method for GetPoint On Surface
731//
733 G4double fRMin2,G4double fRMax2,
734 G4double zOne) const
735{
736 G4double xRand,yRand,phi,cosphi,sinphi,rRand1,rRand2,A1,Atot,rCh;
737
739 cosphi = std::cos(phi);
740 sinphi = std::sin(phi);
741
742 if(fRMin1==fRMin2)
743 {
744 rRand1 = fRMin1; A1=0.;
745 }
746 else
747 {
748 rRand1 = RandFlat::shoot(fRMin1,fRMin2);
749 A1=std::fabs(fRMin2*fRMin2-fRMin1*fRMin1);
750 }
751 if(fRMax1==fRMax2)
752 {
753 rRand2=fRMax1; Atot=A1;
754 }
755 else
756 {
757 rRand2 = RandFlat::shoot(fRMax1,fRMax2);
758 Atot = A1+std::fabs(fRMax2*fRMax2-fRMax1*fRMax1);
759 }
760 rCh = RandFlat::shoot(0.,Atot);
761
762 if(rCh>A1) { rRand1=rRand2; }
763
764 xRand = rRand1*cosphi;
765 yRand = rRand1*sinphi;
766
767 return G4ThreeVector(xRand, yRand, zOne);
768}
769
770
771//
772// GetPointOnCut
773//
774// Auxiliary method for Get Point On Surface
775//
777 G4double fRMin2, G4double fRMax2,
778 G4double zOne, G4double zTwo,
779 G4double& totArea) const
780{ if(zOne==zTwo)
781 {
782 return GetPointOnRing(fRMin1, fRMax1,fRMin2,fRMax2,zOne);
783 }
784 if( (fRMin1 == fRMin2) && (fRMax1 == fRMax2) )
785 {
786 return GetPointOnTubs(fRMin1, fRMax1,zOne,zTwo,totArea);
787 }
788 return GetPointOnCone(fRMin1,fRMax1,fRMin2,fRMax2,zOne,zTwo,totArea);
789}
790
791
792//
793// GetPointOnSurface
794//
796{
797 if (!genericPcon) // Polycone by faces
798 {
799 G4double Area=0,totArea=0,Achose1=0,Achose2=0,phi,cosphi,sinphi,rRand;
800 G4int i=0;
802
804 cosphi = std::cos(phi);
805 sinphi = std::sin(phi);
806
807 rRand = original_parameters->Rmin[0] +
809 * std::sqrt(RandFlat::shoot()) );
810
811 std::vector<G4double> areas; // (numPlanes+1);
812 std::vector<G4ThreeVector> points; // (numPlanes-1);
813
814 areas.push_back(pi*(sqr(original_parameters->Rmax[0])
816
817 for(i=0; i<numPlanes-1; i++)
818 {
820 * std::sqrt(sqr(original_parameters->Rmin[i]
824
826 * std::sqrt(sqr(original_parameters->Rmax[i]
830
831 Area *= 0.5*(endPhi-startPhi);
832
833 if(startPhi==0.&& endPhi == twopi)
834 {
835 Area += std::fabs(original_parameters->Z_values[i+1]
841 }
842 areas.push_back(Area);
843 totArea += Area;
844 }
845
846 areas.push_back(pi*(sqr(original_parameters->Rmax[numPlanes-1])-
847 sqr(original_parameters->Rmin[numPlanes-1])));
848
849 totArea += (areas[0]+areas[numPlanes]);
850 G4double chose = RandFlat::shoot(0.,totArea);
851
852 if( (chose>=0.) && (chose<areas[0]) )
853 {
854 return G4ThreeVector(rRand*cosphi, rRand*sinphi,
856 }
857
858 for (i=0; i<numPlanes-1; i++)
859 {
860 Achose1 += areas[i];
861 Achose2 = (Achose1+areas[i+1]);
862 if(chose>=Achose1 && chose<Achose2)
863 {
869 original_parameters->Z_values[i+1], Area);
870 }
871 }
872
873 rRand = original_parameters->Rmin[numPlanes-1] +
874 ( (original_parameters->Rmax[numPlanes-1]-original_parameters->Rmin[numPlanes-1])
875 * std::sqrt(RandFlat::shoot()) );
876
877
878 return G4ThreeVector(rRand*cosphi,rRand*sinphi,
879 original_parameters->Z_values[numPlanes-1]);
880
881 }
882 else // Generic Polycone
883 {
884 return GetPointOnSurfaceGeneric();
885 }
886}
887
888//
889// CreatePolyhedron
890//
892{
893 //
894 // This has to be fixed in visualization. Fake it for the moment.
895 //
896 if (!genericPcon)
897 {
904 }
905 else
906 {
907 // The following code prepares for:
908 // HepPolyhedron::createPolyhedron(int Nnodes, int Nfaces,
909 // const double xyz[][3],
910 // const int faces_vec[][4])
911 // Here is an extract from the header file HepPolyhedron.h:
912 /**
913 * Creates user defined polyhedron.
914 * This function allows to the user to define arbitrary polyhedron.
915 * The faces of the polyhedron should be either triangles or planar
916 * quadrilateral. Nodes of a face are defined by indexes pointing to
917 * the elements in the xyz array. Numeration of the elements in the
918 * array starts from 1 (like in fortran). The indexes can be positive
919 * or negative. Negative sign means that the corresponding edge is
920 * invisible. The normal of the face should be directed to exterior
921 * of the polyhedron.
922 *
923 * @param Nnodes number of nodes
924 * @param Nfaces number of faces
925 * @param xyz nodes
926 * @param faces_vec faces (quadrilaterals or triangles)
927 * @return status of the operation - is non-zero in case of problem
928 */
929 const G4int numSide =
931 * (endPhi - startPhi) / twopi) + 1;
932 G4int nNodes;
933 G4int nFaces;
934 typedef G4double double3[3];
935 double3* xyz;
936 typedef G4int int4[4];
937 int4* faces_vec;
938 if (phiIsOpen)
939 {
940 // Triangulate open ends. Simple ear-chopping algorithm...
941 // I'm not sure how robust this algorithm is (J.Allison).
942 //
943 std::vector<G4bool> chopped(numCorner, false);
944 std::vector<G4int*> triQuads;
945 G4int remaining = numCorner;
946 G4int iStarter = 0;
947 while (remaining >= 3)
948 {
949 // Find unchopped corners...
950 //
951 G4int A = -1, B = -1, C = -1;
952 G4int iStepper = iStarter;
953 do
954 {
955 if (A < 0) { A = iStepper; }
956 else if (B < 0) { B = iStepper; }
957 else if (C < 0) { C = iStepper; }
958 do
959 {
960 if (++iStepper >= numCorner) { iStepper = 0; }
961 }
962 while (chopped[iStepper]);
963 }
964 while (C < 0 && iStepper != iStarter);
965
966 // Check triangle at B is pointing outward (an "ear").
967 // Sign of z cross product determines...
968 //
969 G4double BAr = corners[A].r - corners[B].r;
970 G4double BAz = corners[A].z - corners[B].z;
971 G4double BCr = corners[C].r - corners[B].r;
972 G4double BCz = corners[C].z - corners[B].z;
973 if (BAr * BCz - BAz * BCr < kCarTolerance)
974 {
975 G4int* tq = new G4int[3];
976 tq[0] = A + 1;
977 tq[1] = B + 1;
978 tq[2] = C + 1;
979 triQuads.push_back(tq);
980 chopped[B] = true;
981 --remaining;
982 }
983 else
984 {
985 do
986 {
987 if (++iStarter >= numCorner) { iStarter = 0; }
988 }
989 while (chopped[iStarter]);
990 }
991 }
992 // Transfer to faces...
993 //
994 nNodes = (numSide + 1) * numCorner;
995 nFaces = numSide * numCorner + 2 * triQuads.size();
996 faces_vec = new int4[nFaces];
997 G4int iface = 0;
998 G4int addition = numCorner * numSide;
999 G4int d = numCorner - 1;
1000 for (G4int iEnd = 0; iEnd < 2; ++iEnd)
1001 {
1002 for (size_t i = 0; i < triQuads.size(); ++i)
1003 {
1004 // Negative for soft/auxiliary/normally invisible edges...
1005 //
1006 G4int a, b, c;
1007 if (iEnd == 0)
1008 {
1009 a = triQuads[i][0];
1010 b = triQuads[i][1];
1011 c = triQuads[i][2];
1012 }
1013 else
1014 {
1015 a = triQuads[i][0] + addition;
1016 b = triQuads[i][2] + addition;
1017 c = triQuads[i][1] + addition;
1018 }
1019 G4int ab = std::abs(b - a);
1020 G4int bc = std::abs(c - b);
1021 G4int ca = std::abs(a - c);
1022 faces_vec[iface][0] = (ab == 1 || ab == d)? a: -a;
1023 faces_vec[iface][1] = (bc == 1 || bc == d)? b: -b;
1024 faces_vec[iface][2] = (ca == 1 || ca == d)? c: -c;
1025 faces_vec[iface][3] = 0;
1026 ++iface;
1027 }
1028 }
1029
1030 // Continue with sides...
1031
1032 xyz = new double3[nNodes];
1033 const G4double dPhi = (endPhi - startPhi) / numSide;
1034 G4double phi = startPhi;
1035 G4int ixyz = 0;
1036 for (G4int iSide = 0; iSide < numSide; ++iSide)
1037 {
1038 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1039 {
1040 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1041 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1042 xyz[ixyz][2] = corners[iCorner].z;
1043 if (iSide == 0) // startPhi
1044 {
1045 if (iCorner < numCorner - 1)
1046 {
1047 faces_vec[iface][0] = ixyz + 1;
1048 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1049 faces_vec[iface][2] = ixyz + numCorner + 2;
1050 faces_vec[iface][3] = ixyz + 2;
1051 }
1052 else
1053 {
1054 faces_vec[iface][0] = ixyz + 1;
1055 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1056 faces_vec[iface][2] = ixyz + 2;
1057 faces_vec[iface][3] = ixyz - numCorner + 2;
1058 }
1059 }
1060 else if (iSide == numSide - 1) // endPhi
1061 {
1062 if (iCorner < numCorner - 1)
1063 {
1064 faces_vec[iface][0] = ixyz + 1;
1065 faces_vec[iface][1] = ixyz + numCorner + 1;
1066 faces_vec[iface][2] = ixyz + numCorner + 2;
1067 faces_vec[iface][3] = -(ixyz + 2);
1068 }
1069 else
1070 {
1071 faces_vec[iface][0] = ixyz + 1;
1072 faces_vec[iface][1] = ixyz + numCorner + 1;
1073 faces_vec[iface][2] = ixyz + 2;
1074 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1075 }
1076 }
1077 else
1078 {
1079 if (iCorner < numCorner - 1)
1080 {
1081 faces_vec[iface][0] = ixyz + 1;
1082 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1083 faces_vec[iface][2] = ixyz + numCorner + 2;
1084 faces_vec[iface][3] = -(ixyz + 2);
1085 }
1086 else
1087 {
1088 faces_vec[iface][0] = ixyz + 1;
1089 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1090 faces_vec[iface][2] = ixyz + 2;
1091 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1092 }
1093 }
1094 ++iface;
1095 ++ixyz;
1096 }
1097 phi += dPhi;
1098 }
1099
1100 // Last corners...
1101
1102 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1103 {
1104 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1105 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1106 xyz[ixyz][2] = corners[iCorner].z;
1107 ++ixyz;
1108 }
1109 }
1110 else // !phiIsOpen - i.e., a complete 360 degrees.
1111 {
1112 nNodes = numSide * numCorner;
1113 nFaces = numSide * numCorner;;
1114 xyz = new double3[nNodes];
1115 faces_vec = new int4[nFaces];
1116 const G4double dPhi = (endPhi - startPhi) / numSide;
1117 G4double phi = startPhi;
1118 G4int ixyz = 0, iface = 0;
1119 for (G4int iSide = 0; iSide < numSide; ++iSide)
1120 {
1121 for (G4int iCorner = 0; iCorner < numCorner; ++iCorner)
1122 {
1123 xyz[ixyz][0] = corners[iCorner].r * std::cos(phi);
1124 xyz[ixyz][1] = corners[iCorner].r * std::sin(phi);
1125 xyz[ixyz][2] = corners[iCorner].z;
1126
1127 if (iSide < numSide - 1)
1128 {
1129 if (iCorner < numCorner - 1)
1130 {
1131 faces_vec[iface][0] = ixyz + 1;
1132 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1133 faces_vec[iface][2] = ixyz + numCorner + 2;
1134 faces_vec[iface][3] = -(ixyz + 2);
1135 }
1136 else
1137 {
1138 faces_vec[iface][0] = ixyz + 1;
1139 faces_vec[iface][1] = -(ixyz + numCorner + 1);
1140 faces_vec[iface][2] = ixyz + 2;
1141 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1142 }
1143 }
1144 else // Last side joins ends...
1145 {
1146 if (iCorner < numCorner - 1)
1147 {
1148 faces_vec[iface][0] = ixyz + 1;
1149 faces_vec[iface][1] = -(ixyz + numCorner - nFaces + 1);
1150 faces_vec[iface][2] = ixyz + numCorner - nFaces + 2;
1151 faces_vec[iface][3] = -(ixyz + 2);
1152 }
1153 else
1154 {
1155 faces_vec[iface][0] = ixyz + 1;
1156 faces_vec[iface][1] = -(ixyz - nFaces + numCorner + 1);
1157 faces_vec[iface][2] = ixyz - nFaces + 2;
1158 faces_vec[iface][3] = -(ixyz - numCorner + 2);
1159 }
1160 }
1161 ++ixyz;
1162 ++iface;
1163 }
1164 phi += dPhi;
1165 }
1166 }
1167 G4Polyhedron* polyhedron = new G4Polyhedron;
1168 G4int problem = polyhedron->createPolyhedron(nNodes, nFaces, xyz, faces_vec);
1169 delete [] faces_vec;
1170 delete [] xyz;
1171 if (problem)
1172 {
1173 std::ostringstream message;
1174 message << "Problem creating G4Polyhedron for: " << GetName();
1175 G4Exception("G4Polycone::CreatePolyhedron()", "GeomSolids1002",
1176 JustWarning, message);
1177 delete polyhedron;
1178 return 0;
1179 }
1180 else
1181 {
1182 return polyhedron;
1183 }
1184 }
1185}
1186
1187
1188//
1189// CreateNURBS
1190//
1192{
1193 return 0;
1194}
1195
1196
1197//
1198// G4PolyconeHistorical stuff
1199//
1200
1202 : Start_angle(0.), Opening_angle(0.), Num_z_planes(0),
1203 Z_values(0), Rmin(0), Rmax(0)
1204{
1205}
1206
1208{
1209 delete [] Z_values;
1210 delete [] Rmin;
1211 delete [] Rmax;
1212}
1213
1216{
1217 Start_angle = source.Start_angle;
1219 Num_z_planes = source.Num_z_planes;
1220
1222 Rmin = new G4double[Num_z_planes];
1223 Rmax = new G4double[Num_z_planes];
1224
1225 for( G4int i = 0; i < Num_z_planes; i++)
1226 {
1227 Z_values[i] = source.Z_values[i];
1228 Rmin[i] = source.Rmin[i];
1229 Rmax[i] = source.Rmax[i];
1230 }
1231}
1232
1235{
1236 if ( &right == this ) return *this;
1237
1238 if (&right)
1239 {
1240 Start_angle = right.Start_angle;
1242 Num_z_planes = right.Num_z_planes;
1243
1244 delete [] Z_values;
1245 delete [] Rmin;
1246 delete [] Rmax;
1248 Rmin = new G4double[Num_z_planes];
1249 Rmax = new G4double[Num_z_planes];
1250
1251 for( G4int i = 0; i < Num_z_planes; i++)
1252 {
1253 Z_values[i] = right.Z_values[i];
1254 Rmin[i] = right.Rmin[i];
1255 Rmax[i] = right.Rmax[i];
1256 }
1257 }
1258 return *this;
1259}
@ JustWarning
@ FatalErrorInArgument
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
static double shoot()
Definition: RandFlat.cc:59
G4bool ShouldMiss(const G4ThreeVector &p, const G4ThreeVector &v) const
G4bool MustBeOutside(const G4ThreeVector &p) const
G4double Opening_angle
Definition: G4Polycone.hh:78
G4double * Z_values
Definition: G4Polycone.hh:80
G4PolyconeHistorical & operator=(const G4PolyconeHistorical &right)
Definition: G4Polycone.cc:1234
G4bool phiIsOpen
Definition: G4Polycone.hh:190
G4EnclosingCylinder * enclosingCylinder
Definition: G4Polycone.hh:198
G4NURBS * CreateNURBS() const
Definition: G4Polycone.cc:1191
G4ThreeVector GetPointOnSurface() const
Definition: G4Polycone.cc:795
void SetOriginalParameters()
virtual ~G4Polycone()
Definition: G4Polycone.cc:326
G4GeometryType GetEntityType() const
Definition: G4Polycone.cc:508
G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
Definition: G4Polycone.cc:470
G4Polycone(const G4String &name, G4double phiStart, G4double phiTotal, G4int numZPlanes, const G4double zPlane[], const G4double rInner[], const G4double rOuter[])
Definition: G4Polycone.cc:57
G4ThreeVector GetPointOnTubs(G4double fRMin, G4double fRMax, G4double zOne, G4double zTwo, G4double &totArea) const
Definition: G4Polycone.cc:673
G4ThreeVector GetPointOnCone(G4double fRmin1, G4double fRmax1, G4double fRmin2, G4double fRmax2, G4double zOne, G4double zTwo, G4double &totArea) const
Definition: G4Polycone.cc:577
const G4Polycone & operator=(const G4Polycone &source)
Definition: G4Polycone.cc:347
G4Polyhedron * CreatePolyhedron() const
Definition: G4Polycone.cc:891
void Create(G4double phiStart, G4double phiTotal, G4ReduciblePolygon *rz)
Definition: G4Polycone.cc:148
EInside Inside(const G4ThreeVector &p) const
Definition: G4Polycone.cc:450
G4double endPhi
Definition: G4Polycone.hh:189
G4VSolid * Clone() const
Definition: G4Polycone.cc:516
void CopyStuff(const G4Polycone &source)
Definition: G4Polycone.cc:367
G4bool genericPcon
Definition: G4Polycone.hh:191
G4PolyconeSideRZ * corners
Definition: G4Polycone.hh:193
void ComputeDimensions(G4VPVParameterisation *p, const G4int n, const G4VPhysicalVolume *pRep)
Definition: G4Polycone.cc:498
G4bool Reset()
Definition: G4Polycone.cc:409
G4int numCorner
Definition: G4Polycone.hh:192
G4PolyconeHistorical * original_parameters
Definition: G4Polycone.hh:194
G4ThreeVector GetPointOnRing(G4double fRMin, G4double fRMax, G4double fRMin2, G4double fRMax2, G4double zOne) const
Definition: G4Polycone.cc:732
std::ostream & StreamInfo(std::ostream &os) const
Definition: G4Polycone.cc:524
G4ThreeVector GetPointOnCut(G4double fRMin1, G4double fRMax1, G4double fRMin2, G4double fRMax2, G4double zOne, G4double zTwo, G4double &totArea) const
Definition: G4Polycone.cc:776
G4double startPhi
Definition: G4Polycone.hh:188
G4bool BisectedBy(G4double a1, G4double b1, G4double a2, G4double b2, G4double tolerance)
G4double Amin() const
G4bool RemoveDuplicateVertices(G4double tolerance)
G4int NumVertices() const
G4bool RemoveRedundantVertices(G4double tolerance)
G4bool CrossesItself(G4double tolerance)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const
virtual EInside Inside(const G4ThreeVector &p) const
G4VCSGface ** faces
const G4VCSGfaceted & operator=(const G4VCSGfaceted &source)
G4ThreeVector GetPointOnSurfaceGeneric() const
virtual void ComputeDimensions(G4Box &, const G4int, const G4VPhysicalVolume *) const
G4String GetName() const
void DumpInfo() const
G4double kCarTolerance
Definition: G4VSolid.hh:307
static G4int GetNumberOfRotationSteps()
G4int createPolyhedron(G4int Nnodes, G4int Nfaces, const G4double xyz[][3], const G4int faces[][4])
EInside
Definition: geomdefs.hh:58
@ kOutside
Definition: geomdefs.hh:58
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Definition: DoubConv.h:17
T sqr(const T &x)
Definition: templates.hh:145