source: trunk/source/geometry/solids/specific/src/G4PolyconeSide.cc @ 1202

Last change on this file since 1202 was 1058, checked in by garnier, 15 years ago

file release beta

File size: 27.3 KB
Line 
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: G4PolyconeSide.cc,v 1.19 2008/05/15 11:41:59 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//
31// --------------------------------------------------------------------
32// GEANT 4 class source file
33//
34//
35// G4PolyconeSide.cc
36//
37// Implementation of the face representing one conical side of a polycone
38//
39// --------------------------------------------------------------------
40
41#include "G4PolyconeSide.hh"
42#include "G4IntersectingCone.hh"
43#include "G4ClippablePolygon.hh"
44#include "G4AffineTransform.hh"
45#include "meshdefs.hh"
46#include "G4SolidExtentList.hh"
47#include "G4GeometryTolerance.hh"
48
49#include "Randomize.hh"
50
51//
52// Constructor
53//
54// Values for r1,z1 and r2,z2 should be specified in clockwise
55// order in (r,z).
56//
57G4PolyconeSide::G4PolyconeSide( const G4PolyconeSideRZ *prevRZ,
58                                const G4PolyconeSideRZ *tail,
59                                const G4PolyconeSideRZ *head,
60                                const G4PolyconeSideRZ *nextRZ,
61                                      G4double thePhiStart, 
62                                      G4double theDeltaPhi, 
63                                      G4bool thePhiIsOpen, 
64                                      G4bool isAllBehind )
65  : ncorners(0), corners(0)
66{
67  kCarTolerance = G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
68  fSurfaceArea = 0.0;
69
70  //
71  // Record values
72  //
73  r[0] = tail->r; z[0] = tail->z;
74  r[1] = head->r; z[1] = head->z;
75 
76  phiIsOpen = thePhiIsOpen;
77  if (phiIsOpen)
78  {
79    deltaPhi = theDeltaPhi;
80    startPhi = thePhiStart;
81
82    //
83    // Set phi values to our conventions
84    //
85    while (deltaPhi < 0.0) deltaPhi += twopi;
86    while (startPhi < 0.0) startPhi += twopi;
87   
88    //
89    // Calculate corner coordinates
90    //
91    ncorners = 4;
92    corners = new G4ThreeVector[ncorners];
93   
94    corners[0] = G4ThreeVector( tail->r*std::cos(startPhi),
95                                tail->r*std::sin(startPhi), tail->z );
96    corners[1] = G4ThreeVector( head->r*std::cos(startPhi),
97                                head->r*std::sin(startPhi), head->z );
98    corners[2] = G4ThreeVector( tail->r*std::cos(startPhi+deltaPhi),
99                                tail->r*std::sin(startPhi+deltaPhi), tail->z );
100    corners[3] = G4ThreeVector( head->r*std::cos(startPhi+deltaPhi),
101                                head->r*std::sin(startPhi+deltaPhi), head->z );
102  }
103  else
104  {
105    deltaPhi = twopi;
106    startPhi = 0.0;
107  }
108 
109  allBehind = isAllBehind;
110   
111  //
112  // Make our intersecting cone
113  //
114  cone = new G4IntersectingCone( r, z );
115 
116  //
117  // Calculate vectors in r,z space
118  //
119  rS = r[1]-r[0]; zS = z[1]-z[0];
120  length = std::sqrt( rS*rS + zS*zS);
121  rS /= length; zS /= length;
122 
123  rNorm = +zS;
124  zNorm = -rS;
125 
126  G4double lAdj;
127 
128  prevRS = r[0]-prevRZ->r;
129  prevZS = z[0]-prevRZ->z;
130  lAdj = std::sqrt( prevRS*prevRS + prevZS*prevZS );
131  prevRS /= lAdj;
132  prevZS /= lAdj;
133
134  rNormEdge[0] = rNorm + prevZS;
135  zNormEdge[0] = zNorm - prevRS;
136  lAdj = std::sqrt( rNormEdge[0]*rNormEdge[0] + zNormEdge[0]*zNormEdge[0] );
137  rNormEdge[0] /= lAdj;
138  zNormEdge[0] /= lAdj;
139
140  nextRS = nextRZ->r-r[1];
141  nextZS = nextRZ->z-z[1];
142  lAdj = std::sqrt( nextRS*nextRS + nextZS*nextZS );
143  nextRS /= lAdj;
144  nextZS /= lAdj;
145
146  rNormEdge[1] = rNorm + nextZS;
147  zNormEdge[1] = zNorm - nextRS;
148  lAdj = std::sqrt( rNormEdge[1]*rNormEdge[1] + zNormEdge[1]*zNormEdge[1] );
149  rNormEdge[1] /= lAdj;
150  zNormEdge[1] /= lAdj;
151}
152
153
154//
155// Fake default constructor - sets only member data and allocates memory
156//                            for usage restricted to object persistency.
157//
158G4PolyconeSide::G4PolyconeSide( __void__& )
159  : phiIsOpen(false), cone(0), ncorners(0), corners(0)
160{
161}
162
163
164//
165// Destructor
166// 
167G4PolyconeSide::~G4PolyconeSide()
168{
169  delete cone;
170  if (phiIsOpen) delete [] corners;
171}
172
173
174//
175// Copy constructor
176//
177G4PolyconeSide::G4PolyconeSide( const G4PolyconeSide &source )
178  : G4VCSGface()
179{
180  CopyStuff( source );
181}
182
183
184//
185// Assignment operator
186//
187G4PolyconeSide& G4PolyconeSide::operator=( const G4PolyconeSide &source )
188{
189  if (this == &source) return *this;
190
191  delete cone;
192  if (phiIsOpen) delete [] corners;
193 
194  CopyStuff( source );
195 
196  return *this;
197}
198
199
200//
201// CopyStuff
202//
203void G4PolyconeSide::CopyStuff( const G4PolyconeSide &source )
204{
205  r[0]    = source.r[0];
206  r[1]    = source.r[1];
207  z[0]    = source.z[0];
208  z[1]    = source.z[1];
209 
210  startPhi  = source.startPhi;
211  deltaPhi  = source.deltaPhi;
212  phiIsOpen  = source.phiIsOpen;
213  allBehind  = source.allBehind;
214
215  kCarTolerance = source.kCarTolerance;
216  fSurfaceArea = source.fSurfaceArea;
217 
218  cone    = new G4IntersectingCone( *source.cone );
219 
220  rNorm    = source.rNorm;
221  zNorm    = source.zNorm;
222  rS    = source.rS;
223  zS    = source.zS;
224  length    = source.length;
225  prevRS    = source.prevRS;
226  prevZS    = source.prevZS;
227  nextRS    = source.nextRS;
228  nextZS    = source.nextZS;
229 
230  rNormEdge[0]   = source.rNormEdge[0];
231  rNormEdge[1]  = source.rNormEdge[1];
232  zNormEdge[0]  = source.zNormEdge[0];
233  zNormEdge[1]  = source.zNormEdge[1];
234 
235  if (phiIsOpen)
236  {
237    ncorners = 4;
238    corners = new G4ThreeVector[ncorners];
239   
240    corners[0] = source.corners[0];
241    corners[1] = source.corners[1];
242    corners[2] = source.corners[2];
243    corners[3] = source.corners[3];
244  }
245}
246
247
248//
249// Intersect
250//
251G4bool G4PolyconeSide::Intersect( const G4ThreeVector &p,
252                                  const G4ThreeVector &v, 
253                                        G4bool outgoing,
254                                        G4double surfTolerance,
255                                        G4double &distance,
256                                        G4double &distFromSurface,
257                                        G4ThreeVector &normal,
258                                        G4bool &isAllBehind )
259{
260  G4double s1, s2;
261  G4double normSign = outgoing ? +1 : -1;
262 
263  isAllBehind = allBehind;
264
265  //
266  // Check for two possible intersections
267  //
268  G4int nside = cone->LineHitsCone( p, v, &s1, &s2 );
269  if (nside == 0) return false;
270   
271  //
272  // Check the first side first, since it is (supposed to be) closest
273  //
274  G4ThreeVector hit = p + s1*v;
275 
276  if (PointOnCone( hit, normSign, p, v, normal ))
277  {
278    //
279    // Good intersection! What about the normal?
280    //
281    if (normSign*v.dot(normal) > 0)
282    {
283      //
284      // We have a valid intersection, but it could very easily
285      // be behind the point. To decide if we tolerate this,
286      // we have to see if the point p is on the surface near
287      // the intersecting point.
288      //
289      // What does it mean exactly for the point p to be "near"
290      // the intersection? It means that if we draw a line from
291      // p to the hit, the line remains entirely within the
292      // tolerance bounds of the cone. To test this, we can
293      // ask if the normal is correct near p.
294      //
295      G4double pr = p.perp();
296      if (pr < DBL_MIN) pr = DBL_MIN;
297      G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
298      if (normSign*v.dot(pNormal) > 0)
299      {
300        //
301        // p and intersection in same hemisphere
302        //
303        G4double distOutside2;
304        distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
305        if (distOutside2 < surfTolerance*surfTolerance)
306        {
307          if (distFromSurface > -surfTolerance)
308          {
309            //
310            // We are just inside or away from the
311            // surface. Accept *any* value of distance.
312            //
313            distance = s1;
314            return true;
315          }
316        }
317      }
318      else 
319        distFromSurface = s1;
320     
321      //
322      // Accept positive distances
323      //
324      if (s1 > 0)
325      {
326        distance = s1;
327        return true;
328      }
329    }
330  } 
331 
332  if (nside==1) return false;
333 
334  //
335  // Well, try the second hit
336  // 
337  hit = p + s2*v;
338 
339  if (PointOnCone( hit, normSign, p, v, normal ))
340  {
341    //
342    // Good intersection! What about the normal?
343    //
344    if (normSign*v.dot(normal) > 0)
345    {
346      G4double pr = p.perp();
347      if (pr < DBL_MIN) pr = DBL_MIN;
348      G4ThreeVector pNormal( rNorm*p.x()/pr, rNorm*p.y()/pr, zNorm );
349      if (normSign*v.dot(pNormal) > 0)
350      {
351        G4double distOutside2;
352        distFromSurface = -normSign*DistanceAway( p, false, distOutside2 );
353        if (distOutside2 < surfTolerance*surfTolerance)
354        {
355          if (distFromSurface > -surfTolerance)
356          {
357            distance = s2;
358            return true;
359          }
360        }
361      }
362      else 
363        distFromSurface = s2;
364     
365      if (s2 > 0)
366      {
367        distance = s2;
368        return true;
369      }
370    }
371  } 
372
373  //
374  // Better luck next time
375  //
376  return false;
377}
378
379
380G4double G4PolyconeSide::Distance( const G4ThreeVector &p, G4bool outgoing )
381{
382  G4double normSign = outgoing ? -1 : +1;
383  G4double distFrom, distOut2;
384 
385  //
386  // We have two tries for each hemisphere. Try the closest first.
387  //
388  distFrom = normSign*DistanceAway( p, false, distOut2 );
389  if (distFrom > -0.5*kCarTolerance )
390  {
391    //
392    // Good answer
393    //
394    if (distOut2 > 0) 
395      return std::sqrt( distFrom*distFrom + distOut2 );
396    else 
397      return std::fabs(distFrom);
398  }
399 
400  //
401  // Try second side.
402  //
403  distFrom = normSign*DistanceAway( p,  true, distOut2 );
404  if (distFrom > -0.5*kCarTolerance)
405  {
406
407    if (distOut2 > 0) 
408      return std::sqrt( distFrom*distFrom + distOut2 );
409    else
410      return std::fabs(distFrom);
411  }
412 
413  return kInfinity;
414}
415
416
417//
418// Inside
419//
420EInside G4PolyconeSide::Inside( const G4ThreeVector &p,
421                                      G4double tolerance, 
422                                      G4double *bestDistance )
423{
424  //
425  // Check both sides
426  //
427  G4double distFrom[2], distOut2[2], dist2[2];
428  G4double edgeRZnorm[2];
429     
430  distFrom[0] =  DistanceAway( p, false, distOut2[0], edgeRZnorm   );
431  distFrom[1] =  DistanceAway( p,  true, distOut2[1], edgeRZnorm+1 );
432 
433  dist2[0] = distFrom[0]*distFrom[0] + distOut2[0];
434  dist2[1] = distFrom[1]*distFrom[1] + distOut2[1];
435 
436  //
437  // Who's closest?
438  //
439  G4int i = std::fabs(dist2[0]) < std::fabs(dist2[1]) ? 0 : 1;
440 
441  *bestDistance = std::sqrt( dist2[i] );
442 
443  //
444  // Okay then, inside or out?
445  //
446  if ( (std::fabs(edgeRZnorm[i]) < tolerance)
447    && (distOut2[i] < tolerance*tolerance) )
448    return kSurface;
449  else if (edgeRZnorm[i] < 0)
450    return kInside;
451  else
452    return kOutside;
453}
454
455
456//
457// Normal
458//
459G4ThreeVector G4PolyconeSide::Normal( const G4ThreeVector &p,
460                                            G4double *bestDistance )
461{
462  if (p == G4ThreeVector(0.,0.,0.))  { return p; }
463
464  G4ThreeVector dFrom;
465  G4double dOut2;
466 
467  dFrom = DistanceAway( p, false, dOut2 );
468 
469  *bestDistance = std::sqrt( dFrom*dFrom + dOut2 );
470 
471  G4double rad = p.perp();
472  return G4ThreeVector( rNorm*p.x()/rad, rNorm*p.y()/rad, zNorm );
473}
474
475
476//
477// Extent
478//
479G4double G4PolyconeSide::Extent( const G4ThreeVector axis )
480{
481  if (axis.perp2() < DBL_MIN)
482  {
483    //
484    // Special case
485    //
486    return axis.z() < 0 ? -cone->ZLo() : cone->ZHi();
487  }
488
489  //
490  // Is the axis pointing inside our phi gap?
491  //
492  if (phiIsOpen)
493  {
494    G4double phi = axis.phi();
495    while( phi < startPhi ) phi += twopi;
496   
497    if (phi > deltaPhi+startPhi)
498    {
499      //
500      // Yeah, looks so. Make four three vectors defining the phi
501      // opening
502      //
503      G4double cosP = std::cos(startPhi), sinP = std::sin(startPhi);
504      G4ThreeVector a( r[0]*cosP, r[0]*sinP, z[0] );
505      G4ThreeVector b( r[1]*cosP, r[1]*sinP, z[1] );
506      cosP = std::cos(startPhi+deltaPhi); sinP = std::sin(startPhi+deltaPhi);
507      G4ThreeVector c( r[0]*cosP, r[0]*sinP, z[0] );
508      G4ThreeVector d( r[1]*cosP, r[1]*sinP, z[1] );
509     
510      G4double ad = axis.dot(a),
511         bd = axis.dot(b),
512         cd = axis.dot(c),
513         dd = axis.dot(d);
514     
515      if (bd > ad) ad = bd;
516      if (cd > ad) ad = cd;
517      if (dd > ad) ad = dd;
518     
519      return ad;
520    }
521  }
522
523  //
524  // Check either end
525  //
526  G4double aPerp = axis.perp();
527 
528  G4double a = aPerp*r[0] + axis.z()*z[0];
529  G4double b = aPerp*r[1] + axis.z()*z[1];
530 
531  if (b > a) a = b;
532 
533  return a;
534}
535
536
537
538//
539// CalculateExtent
540//
541// See notes in G4VCSGface
542//
543void G4PolyconeSide::CalculateExtent( const EAxis axis, 
544                                      const G4VoxelLimits &voxelLimit,
545                                      const G4AffineTransform &transform,
546                                            G4SolidExtentList &extentList )
547{
548  G4ClippablePolygon polygon;
549 
550  //
551  // Here we will approximate (ala G4Cons) and divide our conical section
552  // into segments, like G4Polyhedra. When doing so, the radius
553  // is extented far enough such that the segments always lie
554  // just outside the surface of the conical section we are
555  // approximating.
556  //
557 
558  //
559  // Choose phi size of our segment(s) based on constants as
560  // defined in meshdefs.hh
561  //
562  G4int numPhi = (G4int)(deltaPhi/kMeshAngleDefault) + 1;
563  if (numPhi < kMinMeshSections) 
564    numPhi = kMinMeshSections;
565  else if (numPhi > kMaxMeshSections)
566    numPhi = kMaxMeshSections;
567   
568  G4double sigPhi = deltaPhi/numPhi;
569 
570  //
571  // Determine radius factor to keep segments outside
572  //
573  G4double rFudge = 1.0/std::cos(0.5*sigPhi);
574 
575  //
576  // Decide which radius to use on each end of the side,
577  // and whether a transition mesh is required
578  //
579  // {r0,z0}  - Beginning of this side
580  // {r1,z1}  - Ending of this side
581  // {r2,z0}  - Beginning of transition piece connecting previous
582  //            side (and ends at beginning of this side)
583  //
584  // So, order is 2 --> 0 --> 1.
585  //                    -------
586  //
587  // r2 < 0 indicates that no transition piece is required
588  //
589  G4double r0, r1, r2, z0, z1;
590 
591  r2 = -1;  // By default: no transition piece
592 
593  if (rNorm < -DBL_MIN)
594  {
595    //
596    // This side faces *inward*, and so our mesh has
597    // the same radius
598    //
599    r1 = r[1];
600    z1 = z[1];
601    z0 = z[0];
602    r0 = r[0];
603   
604    r2 = -1;
605   
606    if (prevZS > DBL_MIN)
607    {
608      //
609      // The previous side is facing outwards
610      //
611      if ( prevRS*zS - prevZS*rS > 0 )
612      {
613        //
614        // Transition was convex: build transition piece
615        //
616        if (r[0] > DBL_MIN) r2 = r[0]*rFudge;
617      }
618      else
619      {
620        //
621        // Transition was concave: short this side
622        //
623        FindLineIntersect( z0, r0, zS, rS,
624                           z0, r0*rFudge, prevZS, prevRS*rFudge, z0, r0 );
625      }
626    }
627   
628    if ( nextZS > DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
629    {
630      //
631      // The next side is facing outwards, forming a
632      // concave transition: short this side
633      //
634      FindLineIntersect( z1, r1, zS, rS,
635                         z1, r1*rFudge, nextZS, nextRS*rFudge, z1, r1 );
636    }
637  }
638  else if (rNorm > DBL_MIN)
639  {
640    //
641    // This side faces *outward* and is given a boost to
642    // it radius
643    //
644    r0 = r[0]*rFudge;
645    z0 = z[0];
646    r1 = r[1]*rFudge;
647    z1 = z[1];
648   
649    if (prevZS < -DBL_MIN)
650    {
651      //
652      // The previous side is facing inwards
653      //
654      if ( prevRS*zS - prevZS*rS > 0 )
655      {
656        //
657        // Transition was convex: build transition piece
658        //
659        if (r[0] > DBL_MIN) r2 = r[0];
660      }
661      else
662      {
663        //
664        // Transition was concave: short this side
665        //
666        FindLineIntersect( z0, r0, zS, rS*rFudge,
667                           z0, r[0], prevZS, prevRS, z0, r0 );
668      }
669    }
670   
671    if ( nextZS < -DBL_MIN && (rS*nextZS - zS*nextRS < 0) )
672    {
673      //
674      // The next side is facing inwards, forming a
675      // concave transition: short this side
676      //
677      FindLineIntersect( z1, r1, zS, rS*rFudge,
678                         z1, r[1], nextZS, nextRS, z1, r1 );
679    }
680  }
681  else
682  {
683    //
684    // This side is perpendicular to the z axis (is a disk)
685    //
686    // Whether or not r0 needs a rFudge factor depends
687    // on the normal of the previous edge. Similar with r1
688    // and the next edge. No transition piece is required.
689    //
690    r0 = r[0];
691    r1 = r[1];
692    z0 = z[0];
693    z1 = z[1];
694   
695    if (prevZS > DBL_MIN) r0 *= rFudge;
696    if (nextZS > DBL_MIN) r1 *= rFudge;
697  }
698 
699  //
700  // Loop
701  //
702  G4double phi = startPhi, 
703           cosPhi = std::cos(phi), 
704           sinPhi = std::sin(phi);
705 
706  G4ThreeVector v0( r0*cosPhi, r0*sinPhi, z0 ),
707                    v1( r1*cosPhi, r1*sinPhi, z1 ),
708  v2, w0, w1, w2;
709  transform.ApplyPointTransform( v0 );
710  transform.ApplyPointTransform( v1 );
711 
712  if (r2 >= 0)
713  {
714    v2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
715    transform.ApplyPointTransform( v2 );
716  }
717
718  do
719  {
720    phi += sigPhi;
721    if (numPhi == 1) phi = startPhi+deltaPhi;  // Try to avoid roundoff
722    cosPhi = std::cos(phi), 
723    sinPhi = std::sin(phi);
724   
725    w0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z0 );
726    w1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z1 );
727    transform.ApplyPointTransform( w0 );
728    transform.ApplyPointTransform( w1 );
729   
730    G4ThreeVector deltaV = r0 > r1 ? w0-v0 : w1-v1;
731   
732    //
733    // Build polygon, taking special care to keep the vertices
734    // in order
735    //
736    polygon.ClearAllVertices();
737
738    polygon.AddVertexInOrder( v0 );
739    polygon.AddVertexInOrder( v1 );
740    polygon.AddVertexInOrder( w1 );
741    polygon.AddVertexInOrder( w0 );
742
743    //
744    // Get extent
745    //
746    if (polygon.PartialClip( voxelLimit, axis ))
747    {
748      //
749      // Get dot product of normal with target axis
750      //
751      polygon.SetNormal( deltaV.cross(v1-v0).unit() );
752     
753      extentList.AddSurface( polygon );
754    }
755   
756    if (r2 >= 0)
757    {
758      //
759      // Repeat, for transition piece
760      //
761      w2 = G4ThreeVector( r2*cosPhi, r2*sinPhi, z0 );
762      transform.ApplyPointTransform( w2 );
763
764      polygon.ClearAllVertices();
765
766      polygon.AddVertexInOrder( v2 );
767      polygon.AddVertexInOrder( v0 );
768      polygon.AddVertexInOrder( w0 );
769      polygon.AddVertexInOrder( w2 );
770
771      if (polygon.PartialClip( voxelLimit, axis ))
772      {
773        polygon.SetNormal( deltaV.cross(v0-v2).unit() );
774       
775        extentList.AddSurface( polygon );
776      }
777     
778      v2 = w2;
779    }
780   
781    //
782    // Next vertex
783    //   
784    v0 = w0;
785    v1 = w1;
786  } while( --numPhi > 0 );
787 
788  //
789  // We are almost done. But, it is important that we leave no
790  // gaps in the surface of our solid. By using rFudge, however,
791  // we've done exactly that, if we have a phi segment.
792  // Add two additional faces if necessary
793  //
794  if (phiIsOpen && rNorm > DBL_MIN)
795  {   
796    G4double cosPhi = std::cos(startPhi),
797       sinPhi = std::sin(startPhi);
798
799    G4ThreeVector a0( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
800                  a1( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
801                  b0( r0*cosPhi, r0*sinPhi, z[0] ),
802                  b1( r1*cosPhi, r1*sinPhi, z[1] );
803 
804    transform.ApplyPointTransform( a0 );
805    transform.ApplyPointTransform( a1 );
806    transform.ApplyPointTransform( b0 );
807    transform.ApplyPointTransform( b1 );
808
809    polygon.ClearAllVertices();
810
811    polygon.AddVertexInOrder( a0 );
812    polygon.AddVertexInOrder( a1 );
813    polygon.AddVertexInOrder( b0 );
814    polygon.AddVertexInOrder( b1 );
815   
816    if (polygon.PartialClip( voxelLimit , axis))
817    {
818      G4ThreeVector normal( sinPhi, -cosPhi, 0 );
819      polygon.SetNormal( transform.TransformAxis( normal ) );
820       
821      extentList.AddSurface( polygon );
822    }
823   
824    cosPhi = std::cos(startPhi+deltaPhi);
825    sinPhi = std::sin(startPhi+deltaPhi);
826   
827    a0 = G4ThreeVector( r[0]*cosPhi, r[0]*sinPhi, z[0] ),
828    a1 = G4ThreeVector( r[1]*cosPhi, r[1]*sinPhi, z[1] ),
829    b0 = G4ThreeVector( r0*cosPhi, r0*sinPhi, z[0] ),
830    b1 = G4ThreeVector( r1*cosPhi, r1*sinPhi, z[1] );
831    transform.ApplyPointTransform( a0 );
832    transform.ApplyPointTransform( a1 );
833    transform.ApplyPointTransform( b0 );
834    transform.ApplyPointTransform( b1 );
835
836    polygon.ClearAllVertices();
837
838    polygon.AddVertexInOrder( a0 );
839    polygon.AddVertexInOrder( a1 );
840    polygon.AddVertexInOrder( b0 );
841    polygon.AddVertexInOrder( b1 );
842   
843    if (polygon.PartialClip( voxelLimit, axis ))
844    {
845      G4ThreeVector normal( -sinPhi, cosPhi, 0 );
846      polygon.SetNormal( transform.TransformAxis( normal ) );
847       
848      extentList.AddSurface( polygon );
849    }
850  }
851   
852  return;
853}
854
855
856//
857// DistanceAway
858//
859// Calculate distance of a point from our conical surface, including the effect
860// of any phi segmentation
861//
862// Arguments:
863//  p             - (in) Point to check
864//  opposite      - (in) If true, check opposite hemisphere (see below)
865//  distOutside   - (out) Additional distance outside the edges of the surface
866//  edgeRZnorm    - (out) if negative, point is inside
867//
868//  return value = distance from the conical plane, if extrapolated beyond edges,
869//                 signed by whether the point is in inside or outside the shape
870//
871// Notes:
872//  * There are two answers, depending on which hemisphere is considered.
873//
874G4double G4PolyconeSide::DistanceAway( const G4ThreeVector &p,
875                                             G4bool opposite,
876                                             G4double &distOutside2,
877                                             G4double *edgeRZnorm  )
878{
879  //
880  // Convert our point to r and z
881  //
882  G4double rx = p.perp(), zx = p.z();
883 
884  //
885  // Change sign of r if opposite says we should
886  //
887  if (opposite) rx = -rx;
888 
889  //
890  // Calculate return value
891  //
892  G4double deltaR  = rx - r[0], deltaZ = zx - z[0];
893  G4double answer = deltaR*rNorm + deltaZ*zNorm;
894 
895  //
896  // Are we off the surface in r,z space?
897  //
898  G4double s = deltaR*rS + deltaZ*zS;
899  if (s < 0)
900  {
901    distOutside2 = s*s;
902    if (edgeRZnorm) *edgeRZnorm = deltaR*rNormEdge[0] + deltaZ*zNormEdge[0];
903  }
904  else if (s > length)
905  {
906    distOutside2 = sqr( s-length );
907    if (edgeRZnorm)
908    {
909      G4double deltaR  = rx - r[1], deltaZ = zx - z[1];
910      *edgeRZnorm = deltaR*rNormEdge[1] + deltaZ*zNormEdge[1];
911    }
912  }
913  else
914  {
915    distOutside2 = 0;
916    if (edgeRZnorm) *edgeRZnorm = answer;
917  }
918
919  if (phiIsOpen)
920  {
921    //
922    // Finally, check phi
923    //
924    G4double phi = p.phi();
925    while( phi < startPhi ) phi += twopi;
926   
927    if (phi > startPhi+deltaPhi)
928    {
929      //
930      // Oops. Are we closer to the start phi or end phi?
931      //
932      G4double d1 = phi-startPhi-deltaPhi;
933      while( phi > startPhi ) phi -= twopi;
934      G4double d2 = startPhi-phi;
935     
936      if (d2 < d1) d1 = d2;
937     
938      //
939      // Add result to our distance
940      //
941      G4double dist = d1*rx;
942     
943      distOutside2 += dist*dist;
944      if (edgeRZnorm)
945      {
946        *edgeRZnorm = std::max(std::fabs(*edgeRZnorm),std::fabs(dist));
947      }
948    }
949  }
950
951  return answer;
952}
953
954
955//
956// PointOnCone
957//
958// Decide if a point is on a cone and return normal if it is
959//
960G4bool G4PolyconeSide::PointOnCone( const G4ThreeVector &hit,
961                                          G4double normSign,
962                                    const G4ThreeVector &p,
963                                    const G4ThreeVector &v,
964                                          G4ThreeVector &normal )
965{
966  G4double rx = hit.perp();
967  //
968  // Check radial/z extent, as appropriate
969  //
970  if (!cone->HitOn( rx, hit.z() )) return false;
971 
972  if (phiIsOpen)
973  {
974    G4double phiTolerant = 2.0*kCarTolerance/(rx+kCarTolerance);
975    //
976    // Check phi segment. Here we have to be careful
977    // to use the standard method consistent with
978    // PolyPhiFace. See PolyPhiFace::InsideEdgesExact
979    //
980    G4double phi = hit.phi();
981    while( phi < startPhi-phiTolerant ) phi += twopi;
982   
983    if (phi > startPhi+deltaPhi+phiTolerant) return false;
984   
985    if (phi > startPhi+deltaPhi-phiTolerant)
986    {
987      //
988      // Exact treatment
989      //
990      G4ThreeVector qx = p + v;
991      G4ThreeVector qa = qx - corners[2],
992              qb = qx - corners[3];
993      G4ThreeVector qacb = qa.cross(qb);
994     
995      if (normSign*qacb.dot(v) < 0) return false;
996    }
997    else if (phi < phiTolerant)
998    {
999      G4ThreeVector qx = p + v;
1000      G4ThreeVector qa = qx - corners[1],
1001              qb = qx - corners[0];
1002      G4ThreeVector qacb = qa.cross(qb);
1003     
1004      if (normSign*qacb.dot(v) < 0) return false;
1005    }
1006  }
1007 
1008  //
1009  // We have a good hit! Calculate normal
1010  //
1011  if (rx < DBL_MIN) 
1012    normal = G4ThreeVector( 0, 0, zNorm < 0 ? -1 : 1 );
1013  else
1014    normal = G4ThreeVector( rNorm*hit.x()/rx, rNorm*hit.y()/rx, zNorm );
1015  return true;
1016}
1017
1018
1019//
1020// FindLineIntersect
1021//
1022// Decide the point at which two 2-dimensional lines intersect
1023//
1024// Equation of line: x = x1 + s*tx1
1025//                   y = y1 + s*ty1
1026//
1027// It is assumed that the lines are *not* parallel
1028//
1029void G4PolyconeSide::FindLineIntersect( G4double x1,  G4double y1,
1030                                        G4double tx1, G4double ty1,
1031                                        G4double x2,  G4double y2,
1032                                        G4double tx2, G4double ty2,
1033                                        G4double &x,  G4double &y )
1034{
1035  //
1036  // The solution is a simple linear equation
1037  //
1038  G4double deter = tx1*ty2 - tx2*ty1;
1039 
1040  G4double s1 = ((x2-x1)*ty2 - tx2*(y2-y1))/deter;
1041  G4double s2 = ((x2-x1)*ty1 - tx1*(y2-y1))/deter;
1042
1043  //
1044  // We want the answer to not depend on which order the
1045  // lines were specified. Take average.
1046  //
1047  x = 0.5*( x1+s1*tx1 + x2+s2*tx2 );
1048  y = 0.5*( y1+s1*ty1 + y2+s2*ty2 );
1049}
1050
1051//
1052// Calculate surface area for GetPointOnSurface()
1053//
1054G4double G4PolyconeSide::SurfaceArea() 
1055{ 
1056  if(fSurfaceArea==0)
1057  {
1058    fSurfaceArea = (r[0]+r[1])* std::sqrt(sqr(r[0]-r[1])+sqr(z[0]-z[1]));
1059    fSurfaceArea *= 0.5*(deltaPhi);
1060  } 
1061  return fSurfaceArea;
1062}
1063
1064//
1065// GetPointOnFace
1066//
1067G4ThreeVector G4PolyconeSide::GetPointOnFace()
1068{
1069  G4double x,y,zz;
1070  G4double rr,phi,dz,dr;
1071  dr=r[1]-r[0];dz=z[1]-z[0];
1072  phi=startPhi+deltaPhi*G4UniformRand();
1073  rr=r[0]+dr*G4UniformRand();
1074 
1075  x=rr*std::cos(phi);
1076  y=rr*std::sin(phi);
1077
1078  // PolyconeSide has a Ring Form
1079  //
1080  if (dz==0.)
1081  {
1082    zz=z[0];
1083  }
1084  else
1085  {
1086    if(dr==0.)  // PolyconeSide has a Tube Form
1087    {
1088      zz = z[0]+dz*G4UniformRand();
1089    }
1090    else
1091    {
1092      zz = z[0]+(rr-r[0])*dz/dr;
1093    }
1094  }
1095
1096  return G4ThreeVector(x,y,zz);
1097}
Note: See TracBrowser for help on using the repository browser.