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

Last change on this file since 1315 was 1315, checked in by garnier, 14 years ago

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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