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

Last change on this file since 834 was 831, checked in by garnier, 17 years ago

import all except CVS

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