Geant4  10.01.p02
G4ReplicaNavigation.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: G4ReplicaNavigation.cc 90836 2015-06-10 09:31:06Z gcosmo $
28 //
29 //
30 // class G4ReplicaNavigation Implementation
31 //
32 // Author: P.Kent, 1996
33 //
34 // --------------------------------------------------------------------
35 
36 #include "G4ReplicaNavigation.hh"
37 
38 #include "G4AffineTransform.hh"
39 #include "G4SmartVoxelProxy.hh"
40 #include "G4SmartVoxelNode.hh"
41 #include "G4VSolid.hh"
42 #include "G4GeometryTolerance.hh"
43 
44 namespace
45 {
46  const G4ThreeVector VecCartAxes[3]=
47  { G4ThreeVector(1.,0.,0.), G4ThreeVector(0.,1.,0.), G4ThreeVector(0.,0.,1.) };
48  const G4ExitNormal::ESide SideCartAxesPlus[3]=
50  const G4ExitNormal::ESide SideCartAxesMinus[3]=
51  { G4ExitNormal::kMX, G4ExitNormal::kMX, G4ExitNormal::kMX };
52 }
53 
54 // ********************************************************************
55 // Constructor
56 // ********************************************************************
57 //
59 : fCheck(false), fVerbose(0)
60 {
64 }
65 
66 // ********************************************************************
67 // Destructor
68 // ********************************************************************
69 //
71 {
72 }
73 
74 // ********************************************************************
75 // Inside
76 // ********************************************************************
77 //
78 EInside
80  const G4int replicaNo,
81  const G4ThreeVector &localPoint) const
82 {
83  EInside in = kOutside;
84 
85  // Replication data
86  //
87  EAxis axis;
88  G4int nReplicas;
89  G4double width, offset;
90  G4bool consuming;
91 
92  G4double coord, rad2, rmin, tolRMax2, rmax, tolRMin2;
93 
94  pVol->GetReplicationData(axis, nReplicas, width, offset, consuming);
95 
96  switch (axis)
97  {
98  case kXAxis:
99  case kYAxis:
100  case kZAxis:
101  coord = std::fabs(localPoint(axis))-width*0.5;
102  if ( coord<=-kCarTolerance*0.5 )
103  {
104  in = kInside;
105  }
106  else if ( coord<=kCarTolerance*0.5 )
107  {
108  in = kSurface;
109  }
110  break;
111  case kPhi:
112  if ( localPoint.y()||localPoint.x() )
113  {
114  coord = std::fabs(std::atan2(localPoint.y(),localPoint.x()))-width*0.5;
115  if ( coord<=-kAngTolerance*0.5 )
116  {
117  in = kInside;
118  }
119  else if ( coord<=kAngTolerance*0.5 )
120  {
121  in = kSurface;
122  }
123  }
124  else
125  {
126  in = kSurface;
127  }
128  break;
129  case kRho:
130  rad2 = localPoint.perp2();
131  rmax = (replicaNo+1)*width+offset;
132  tolRMax2 = rmax-kRadTolerance*0.5;
133  tolRMax2 *= tolRMax2;
134  if ( rad2>tolRMax2 )
135  {
136  tolRMax2 = rmax+kRadTolerance*0.5;
137  tolRMax2 *= tolRMax2;
138  if ( rad2<=tolRMax2 )
139  {
140  in = kSurface;
141  }
142  }
143  else
144  {
145  // Known to be inside outer radius
146  //
147  if ( replicaNo||offset )
148  {
149  rmin = rmax-width;
150  tolRMin2 = rmin-kRadTolerance*0.5;
151  tolRMin2 *= tolRMin2;
152  if ( rad2>tolRMin2 )
153  {
154  tolRMin2 = rmin+kRadTolerance*0.5;
155  tolRMin2 *= tolRMin2;
156  if ( rad2>=tolRMin2 )
157  {
158  in = kInside;
159  }
160  else
161  {
162  in = kSurface;
163  }
164  }
165  }
166  else
167  {
168  in = kInside;
169  }
170  }
171  break;
172  default:
173  G4Exception("G4ReplicaNavigation::Inside()", "GeomNav0002",
174  FatalException, "Unknown axis!");
175  break;
176  }
177  return in;
178 }
179 
180 // ********************************************************************
181 // DistanceToOut
182 // ********************************************************************
183 //
184 G4double
186  const G4int replicaNo,
187  const G4ThreeVector &localPoint) const
188 {
189  // Replication data
190  //
191  EAxis axis;
192  G4int nReplicas;
193  G4double width,offset;
194  G4bool consuming;
195 
196  G4double safety=0.;
197  G4double safe1,safe2;
198  G4double coord, rho, rmin, rmax;
199 
200  pVol->GetReplicationData(axis, nReplicas, width, offset, consuming);
201  switch(axis)
202  {
203  case kXAxis:
204  case kYAxis:
205  case kZAxis:
206  coord = localPoint(axis);
207  safe1 = width*0.5-coord;
208  safe2 = width*0.5+coord;
209  safety = (safe1<=safe2) ? safe1 : safe2;
210  break;
211  case kPhi:
212  if ( localPoint.y()<=0 )
213  {
214  safety = localPoint.x()*std::sin(width*0.5)
215  + localPoint.y()*std::cos(width*0.5);
216  }
217  else
218  {
219  safety = localPoint.x()*std::sin(width*0.5)
220  - localPoint.y()*std::cos(width*0.5);
221  }
222  break;
223  case kRho:
224  rho = localPoint.perp();
225  rmax = width*(replicaNo+1)+offset;
226  if ( replicaNo||offset )
227  {
228  rmin = rmax-width;
229  safe1 = rho-rmin;
230  safe2 = rmax-rho;
231  safety = (safe1<=safe2) ? safe1 : safe2;
232  }
233  else
234  {
235  safety = rmax-rho;
236  }
237  break;
238  default:
239  G4Exception("G4ReplicaNavigation::DistanceToOut()", "GeomNav0002",
240  FatalException, "Unknown axis!");
241  break;
242  }
243  return (safety >= kCarTolerance) ? safety : 0;
244 }
245 
246 // ********************************************************************
247 // DistanceToOut
248 // ********************************************************************
249 //
250 G4double
252  const G4int replicaNo,
253  const G4ThreeVector &localPoint,
254  const G4ThreeVector &localDirection,
255  G4ExitNormal& arExitNormal ) const
256 {
257  // Replication data
258  //
259  EAxis axis;
260  G4int nReplicas;
261  G4double width, offset;
262  G4bool consuming;
263 
264  G4double Dist=kInfinity;
265  G4double coord, Comp, lindist;
266  G4double signC = 0.0;
267  G4ExitNormal candidateNormal;
268 
269  pVol->GetReplicationData(axis, nReplicas, width, offset, consuming);
270  switch(axis)
271  {
272  case kXAxis:
273  case kYAxis:
274  case kZAxis:
275  coord = localPoint(axis);
276  Comp = localDirection(axis);
277  if ( Comp>0 )
278  {
279  lindist = width*0.5-coord;
280  Dist = (lindist>0) ? lindist/Comp : 0;
281  signC= 1.0;
282  }
283  else if ( Comp<0 )
284  {
285  lindist = width*0.5+coord;
286  Dist = (lindist>0) ? -lindist/Comp : 0;
287  signC= -1.0;
288  }
289  else
290  {
291  Dist = kInfinity;
292  }
293  // signC = sign<G4double>(Comp)
294  candidateNormal.exitNormal = ( signC * VecCartAxes[axis]);
295  candidateNormal.calculated = true;
296  candidateNormal.validConvex = true;
297  candidateNormal.exitSide =
298  (Comp>0) ? SideCartAxesPlus[axis] : SideCartAxesMinus[axis];
299  break;
300  case kPhi:
301  Dist = DistanceToOutPhi(localPoint,localDirection,width,candidateNormal);
302  // candidateNormal set in call
303  break;
304  case kRho:
305  Dist = DistanceToOutRad(localPoint,localDirection,width,offset,
306  replicaNo,candidateNormal);
307  // candidateNormal set in call
308  break;
309  default:
310  G4Exception("G4ReplicaNavigation::DistanceToOut()", "GeomNav0002",
311  FatalException, "Unknown axis!");
312  break;
313  }
314 
315  arExitNormal= candidateNormal; // .exitNormal;
316 
317  return Dist;
318 }
319 
320 // ********************************************************************
321 // DistanceToOutPhi
322 // ********************************************************************
323 //
324 G4double
326  const G4ThreeVector &localDirection,
327  const G4double width,
328  G4ExitNormal& foundNormal ) const
329 {
330  // Phi Intersection
331  // NOTE: width<=pi by definition
332  //
333  G4double sinSPhi= -2.0, cosSPhi= -2.0;
334  G4double pDistS, pDistE, compS, compE, Dist, dist2, yi;
336  G4ThreeVector candidateNormal;
337 
338  if ( (localPoint.x()!=0.0) || (localPoint.y()!=0.0) )
339  {
340  sinSPhi = std::sin(-width*0.5); // SIN of starting phi plane
341  cosSPhi = std::cos(width*0.5); // COS of starting phi plane
342 
343  // pDist -ve when inside
344  //
345  pDistS = localPoint.x()*sinSPhi-localPoint.y()*cosSPhi;
346  // Start plane at phi= -S
347  pDistE = localPoint.x()*sinSPhi+localPoint.y()*cosSPhi;
348  // End plane at phi= +S
349 
350  // Comp -ve when in direction of outwards normal
351  //
352  compS = -sinSPhi*localDirection.x()+cosSPhi*localDirection.y();
353  compE = -sinSPhi*localDirection.x()-cosSPhi*localDirection.y();
354 
355  if ( (pDistS<=0)&&(pDistE<=0) )
356  {
357  // Inside both phi *full* planes
358  //
359  if ( compS<0 )
360  {
361  dist2 = pDistS/compS;
362  yi = localPoint.y()+dist2*localDirection.y();
363 
364  // Check intersecting with correct half-plane (no -> no intersect)
365  //
366  if ( yi<=0 )
367  {
368  Dist = (pDistS<=-kCarTolerance*0.5) ? dist2 : 0;
369  sidePhi= G4ExitNormal::kSPhi; // tbc
370  }
371  else
372  {
373  Dist = kInfinity;
374  }
375  }
376  else
377  {
378  Dist = kInfinity;
379  }
380  if ( compE<0 )
381  {
382  dist2 = pDistE/compE;
383 
384  // Only check further if < starting phi intersection
385  //
386  if ( dist2<Dist )
387  {
388  yi = localPoint.y()+dist2*localDirection.y();
389 
390  // Check intersecting with correct half-plane
391  //
392  if ( yi>=0 )
393  {
394  // Leaving via ending phi
395  //
396  Dist = (pDistE<=-kCarTolerance*0.5) ? dist2 : 0;
397  sidePhi = G4ExitNormal::kEPhi;
398  }
399  }
400  }
401  }
402  else if ( (pDistS>=0)&&(pDistE>=0) )
403  {
404  // Outside both *full* phi planes
405  // if towards both >=0 then once inside will remain inside
406  //
407  Dist = ((compS>=0)&&(compE>=0)) ? kInfinity : 0;
408  }
409  else if ( (pDistS>0)&&(pDistE<0) )
410  {
411  // Outside full starting plane, inside full ending plane
412  //
413  if ( compE<0 )
414  {
415  dist2 = pDistE/compE;
416  yi = localPoint.y()+dist2*localDirection.y();
417 
418  // Check intersection in correct half-plane
419  // (if not -> remain in extent)
420  //
421  Dist = (yi>0) ? dist2 : kInfinity;
422  if( yi> 0 ) { sidePhi = G4ExitNormal::kEPhi; }
423  }
424  else // Leaving immediately by starting phi
425  {
426  Dist = kInfinity;
427  }
428  }
429  else
430  {
431  // Must be (pDistS<0)&&(pDistE>0)
432  // Inside full starting plane, outside full ending plane
433  //
434  if ( compE>=0 )
435  {
436  if ( compS<0 )
437  {
438  dist2 = pDistS/compS;
439  yi = localPoint.y()+dist2*localDirection.y();
440 
441  // Check intersection in correct half-plane
442  // (if not -> remain in extent)
443  //
444  Dist = (yi<0) ? dist2 : kInfinity;
445  if(yi<0) { sidePhi = G4ExitNormal::kSPhi; }
446  }
447  else
448  {
449  Dist = kInfinity;
450  }
451  }
452  else
453  {
454  // Leaving immediately by ending phi
455  //
456  Dist = 0;
457  sidePhi= G4ExitNormal::kEPhi;
458  }
459  }
460  }
461  else
462  {
463  // On z axis + travel not || to z axis -> use direction vector
464  //
465  if( (std::fabs(localDirection.phi())<=width*0.5) )
466  {
467  Dist= kInfinity;
468  }
469  else
470  {
471  Dist= 0;
472  sidePhi= G4ExitNormal::kMY;
473  }
474  }
475 
476  if(sidePhi == G4ExitNormal::kSPhi )
477  {
478  candidateNormal = G4ThreeVector(sinSPhi,-cosSPhi,0.) ;
479  }
480  else if (sidePhi == G4ExitNormal::kEPhi)
481  {
482  candidateNormal = G4ThreeVector(sinSPhi,cosSPhi,0.) ;
483  }
484  else if (sidePhi == G4ExitNormal::kMY )
485  {
486  candidateNormal = G4ThreeVector(0., -1.0, 0.); // Split -S and +S 'phi'
487  }
488  foundNormal.calculated= (sidePhi != G4ExitNormal::kNull );
489  foundNormal.exitNormal= candidateNormal;
490 
491  return Dist;
492 }
493 
494 // ********************************************************************
495 // DistanceToOutRad
496 // ********************************************************************
497 //
498 G4double
500  const G4ThreeVector &localDirection,
501  const G4double width,
502  const G4double offset,
503  const G4int replicaNo,
504  G4ExitNormal& foundNormal ) const
505 {
506  G4double rmin, rmax, t1, t2, t3, deltaR;
507  G4double b, c, d2, srd;
509 
510  //
511  // Radial Intersections
512  //
513 
514  // Find intersction with cylinders at rmax/rmin
515  // Intersection point (xi,yi,zi) on line
516  // x=localPoint.x+t*localDirection.x etc.
517  //
518  // Intersects with x^2+y^2=R^2
519  //
520  // Hence (localDirection.x^2+localDirection.y^2)t^2+
521  // 2t(localPoint.x*localDirection.x+localPoint.y*localDirection.y)+
522  // localPoint.x^2+localPoint.y^2-R^2=0
523  //
524  // t1 t2 t3
525 
526  rmin = replicaNo*width+offset;
527  rmax = (replicaNo+1)*width+offset;
528 
529  t1 = 1.0-localDirection.z()*localDirection.z(); // since v normalised
530  t2 = localPoint.x()*localDirection.x()+localPoint.y()*localDirection.y();
531  t3 = localPoint.x()*localPoint.x()+localPoint.y()*localPoint.y();
532 
533  if ( t1>0 ) // Check not parallel
534  {
535  // Calculate srd, r exit distance
536  //
537  if ( t2>=0 )
538  {
539  // Delta r not negative => leaving via rmax
540  //
541  deltaR = t3-rmax*rmax;
542 
543  // NOTE: Should use
544  // rho-rmax<-kRadTolerance*0.5 - [no sqrts for efficiency]
545  //
546  if ( deltaR<-kRadTolerance*0.5 )
547  {
548  b = t2/t1;
549  c = deltaR/t1;
550  srd = -b+std::sqrt(b*b-c);
551  sideR= G4ExitNormal::kRMax;
552  }
553  else
554  {
555  // On tolerant boundary & heading outwards (or locally
556  // perpendicular to) outer radial surface -> leaving immediately
557  //
558  srd = 0;
559  sideR= G4ExitNormal::kRMax;
560  }
561  }
562  else
563  {
564  // Possible rmin intersection
565  //
566  if (rmin)
567  {
568  deltaR = t3-rmin*rmin;
569  b = t2/t1;
570  c = deltaR/t1;
571  d2 = b*b-c;
572  if ( d2>=0 )
573  {
574  // Leaving via rmin
575  // NOTE: Should use
576  // rho-rmin>kRadTolerance*0.5 - [no sqrts for efficiency]
577  //
578  srd = (deltaR>kRadTolerance*0.5) ? -b-std::sqrt(d2) : 0.0;
579  // Is the following more accurate ?
580  // srd = (deltaR>kRadTolerance*0.5) ? c/( -b - std::sqrt(d2)) : 0.0;
581  sideR= G4ExitNormal::kRMin;
582  }
583  else
584  {
585  // No rmin intersect -> must be rmax intersect
586  //
587  deltaR = t3-rmax*rmax;
588  c = deltaR/t1;
589  d2 = b*b-c;
590  srd = (d2 < 0.) ? 0.0 : -b+std::sqrt(d2);
591  sideR= G4ExitNormal::kRMax;
592  }
593  }
594  else
595  {
596  // No rmin intersect -> must be rmax intersect
597  //
598  deltaR = t3-rmax*rmax;
599  b = t2/t1;
600  c = deltaR/t1;
601  d2 = b*b-c;
602  srd = (d2 < 0.) ? 0.0 : -b+std::sqrt(d2);
603  sideR= G4ExitNormal::kRMax;
604  }
605  }
606  }
607  else
608  {
609  srd=kInfinity;
610  sideR= G4ExitNormal::kNull;
611  }
612 
613  if( sideR != G4ExitNormal::kNull ) // if ((side == kRMax) || (side==kRMin))
614  {
615  // Note: returned vector not explicitly normalised
616  // (divided by fRMax for unit vector)
617  G4double xi, yi;
618  xi = localPoint.x() + srd*localDirection.x() ;
619  yi = localPoint.y() + srd*localDirection.y() ;
620  G4ThreeVector normalR = G4ThreeVector(xi,yi,0.0);
621 
622  if( sideR == G4ExitNormal::kRMax ){
623  normalR *= 1.0/rmax;
624  }else{
625  normalR *= (-1.0)/rmin;
626  }
627  foundNormal.exitNormal= normalR;
628  foundNormal.calculated= true;
629  foundNormal.validConvex = (sideR == G4ExitNormal::kRMax);
630  foundNormal.exitSide = sideR;
631  }else{
632  foundNormal.calculated= false;
633  }
634 
635  return srd;
636 }
637 
638 // ********************************************************************
639 // ComputeTransformation
640 //
641 // Setup transformation and transform point into local system
642 // ********************************************************************
643 //
644 void
646  G4VPhysicalVolume* pVol,
647  G4ThreeVector& point) const
648 {
649  G4double val,cosv,sinv,tmpx,tmpy;
650 
651  // Replication data
652  //
653  EAxis axis;
654  G4int nReplicas;
655  G4double width,offset;
656  G4bool consuming;
657 
658  pVol->GetReplicationData(axis, nReplicas, width, offset, consuming);
659 
660  switch (axis)
661  {
662  case kXAxis:
663  val = -width*0.5*(nReplicas-1)+width*replicaNo;
664  pVol->SetTranslation(G4ThreeVector(val,0,0));
665  point.setX(point.x()-val);
666  break;
667  case kYAxis:
668  val = -width*0.5*(nReplicas-1)+width*replicaNo;
669  pVol->SetTranslation(G4ThreeVector(0,val,0));
670  point.setY(point.y()-val);
671  break;
672  case kZAxis:
673  val = -width*0.5*(nReplicas-1)+width*replicaNo;
674  pVol->SetTranslation(G4ThreeVector(0,0,val));
675  point.setZ(point.z()-val);
676  break;
677  case kPhi:
678  val = -(offset+width*(replicaNo+0.5));
679  SetPhiTransformation(val,pVol);
680  cosv = std::cos(val);
681  sinv = std::sin(val);
682  tmpx = point.x()*cosv-point.y()*sinv;
683  tmpy = point.x()*sinv+point.y()*cosv;
684  point.setY(tmpy);
685  point.setX(tmpx);
686  break;
687  case kRho:
688  // No setup required for radial case
689  default:
690  break;
691  }
692 }
693 
694 // ********************************************************************
695 // ComputeTransformation
696 //
697 // Setup transformation into local system
698 // ********************************************************************
699 //
700 void
702  G4VPhysicalVolume* pVol) const
703 {
704  G4double val;
705 
706  // Replication data
707  //
708  EAxis axis;
709  G4int nReplicas;
710  G4double width, offset;
711  G4bool consuming;
712 
713  pVol->GetReplicationData(axis, nReplicas, width, offset, consuming);
714 
715  switch (axis)
716  {
717  case kXAxis:
718  val = -width*0.5*(nReplicas-1)+width*replicaNo;
719  pVol->SetTranslation(G4ThreeVector(val,0,0));
720  break;
721  case kYAxis:
722  val = -width*0.5*(nReplicas-1)+width*replicaNo;
723  pVol->SetTranslation(G4ThreeVector(0,val,0));
724  break;
725  case kZAxis:
726  val = -width*0.5*(nReplicas-1)+width*replicaNo;
727  pVol->SetTranslation(G4ThreeVector(0,0,val));
728  break;
729  case kPhi:
730  val = -(offset+width*(replicaNo+0.5));
731  SetPhiTransformation(val,pVol);
732  break;
733  case kRho:
734  // No setup required for radial case
735  default:
736  break;
737  }
738 }
739 
740 // ********************************************************************
741 // ComputeStep
742 // ********************************************************************
743 //
744 G4double
746  const G4ThreeVector &globalDirection,
747  const G4ThreeVector &localPoint,
748  const G4ThreeVector &localDirection,
749  const G4double currentProposedStepLength,
750  G4double &newSafety,
751  G4NavigationHistory &history,
752  // std::pair<G4bool,G4bool> &validAndCalculated
753  G4bool &validExitNormal,
754  G4bool &calculatedExitNormal,
755  G4ThreeVector &exitNormalVector,
756  G4bool &exiting,
757  G4bool &entering,
758  G4VPhysicalVolume *(*pBlockedPhysical),
759  G4int &blockedReplicaNo )
760 {
761  G4VPhysicalVolume *repPhysical, *motherPhysical;
762  G4VPhysicalVolume *samplePhysical, *blockedExitedVol=0;
763  G4LogicalVolume *repLogical;
764  G4VSolid *motherSolid;
765  G4ThreeVector repPoint, repDirection, sampleDirection;
766  G4double ourStep=currentProposedStepLength;
767  G4double ourSafety=kInfinity;
768  G4double sampleStep, sampleSafety, motherStep, motherSafety;
769  G4int localNoDaughters, sampleNo;
770  G4int depth;
771  G4ExitNormal exitNormalStc;
772  // G4int depthDeterminingStep= -1; // Useful only for debugging - for now
773 
774  calculatedExitNormal= false;
775 
776  // Exiting normal optimisation
777  //
778  if ( exiting&&validExitNormal )
779  {
780  if ( localDirection.dot(exitNormalVector)>=kMinExitingNormalCosine )
781  {
782  // Block exited daughter volume
783  //
784  blockedExitedVol = *pBlockedPhysical;
785  ourSafety = 0;
786  }
787  }
788  exiting = false;
789  entering = false;
790 
791  repPhysical = history.GetTopVolume();
792  repLogical = repPhysical->GetLogicalVolume();
793 
794  //
795  // Compute intersection with replica boundaries & replica safety
796  //
797 
798  sampleSafety = DistanceToOut(repPhysical,
799  history.GetTopReplicaNo(),
800  localPoint);
801  G4ExitNormal normalOutStc;
802  const G4int topDepth= history.GetDepth();
803 
804  ourSafety = std::min( ourSafety, sampleSafety);
805 
806  if ( sampleSafety<ourStep )
807  {
808 
809  sampleStep = DistanceToOut(repPhysical,
810  history.GetTopReplicaNo(),
811  localPoint,
812  localDirection,
813  normalOutStc);
814  if ( sampleStep<ourStep )
815  {
816  ourStep = sampleStep;
817  exiting = true;
818  validExitNormal = normalOutStc.validConvex; // false; -> Old,Conservative
819 
820  exitNormalStc= normalOutStc;
821  exitNormalStc.exitNormal= history.GetTopTransform().Inverse().
822  TransformAxis(normalOutStc.exitNormal);
823  calculatedExitNormal= true;
824  }
825  }
826  const G4int secondDepth= topDepth;
827  depth = secondDepth;
828 
829  while ( history.GetVolumeType(depth)==kReplica )
830  {
831  const G4AffineTransform& GlobalToLocal= history.GetTransform(depth);
832  repPoint = GlobalToLocal.TransformPoint(globalPoint);
833  // repPoint = history.GetTransform(depth).TransformPoint(globalPoint);
834 
835  sampleSafety = DistanceToOut(history.GetVolume(depth),
836  history.GetReplicaNo(depth),
837  repPoint);
838  if ( sampleSafety < ourSafety )
839  {
840  ourSafety = sampleSafety;
841  }
842  if ( sampleSafety < ourStep )
843  {
844  G4ThreeVector newLocalDirection = GlobalToLocal.TransformAxis(globalDirection);
845  sampleStep = DistanceToOut(history.GetVolume(depth),
846  history.GetReplicaNo(depth),
847  repPoint,
848  newLocalDirection,
849  normalOutStc);
850  if ( sampleStep < ourStep )
851  {
852  ourStep = sampleStep;
853  exiting = true;
854 
855  // As step is limited by this level, must set Exit Normal
856  //
857  G4ThreeVector localExitNorm= normalOutStc.exitNormal;
858  G4ThreeVector globalExitNorm=
859  GlobalToLocal.Inverse().TransformAxis(localExitNorm);
860 
861  exitNormalStc= normalOutStc; // Normal, convex, calculated, side
862  exitNormalStc.exitNormal= globalExitNorm;
863  calculatedExitNormal= true;
864  }
865  }
866  depth--;
867  }
868 
869  // Compute mother safety & intersection
870  //
871  G4ThreeVector exitVectorMother;
872  G4bool exitConvex= false; // Value obtained in DistanceToOut(p,v) call
873  G4ExitNormal motherNormalStc;
874 
875  repPoint = history.GetTransform(depth).TransformPoint(globalPoint);
876  motherPhysical = history.GetVolume(depth);
877  motherSolid = motherPhysical->GetLogicalVolume()->GetSolid();
878  motherSafety = motherSolid->DistanceToOut(repPoint);
879  repDirection = history.GetTransform(depth).TransformAxis(globalDirection);
880 
881  motherStep = motherSolid->DistanceToOut(repPoint,repDirection,true,
882  &exitConvex,&exitVectorMother);
883  if( exitConvex )
884  {
885  motherNormalStc = G4ExitNormal( exitVectorMother, true, false,
887  calculatedExitNormal= true;
888  }
889  const G4AffineTransform& globalToLocalTop = history.GetTopTransform();
890 
891  G4bool motherDeterminedStep= (motherStep<ourStep);
892 
893  if( (!exitConvex) && motherDeterminedStep )
894  {
895  exitVectorMother= motherSolid->SurfaceNormal( repPoint );
896  motherNormalStc= G4ExitNormal( exitVectorMother, true, false,
898  // CalculatedExitNormal -> true;
899  // Convex -> false: do not know value
900  // ExitSide -> kMother (or kNull)
901 
902  calculatedExitNormal= true;
903  }
904  if( motherDeterminedStep)
905  {
906  G4ThreeVector globalExitNormalTop=
907  globalToLocalTop.Inverse().TransformAxis(exitVectorMother);
908 
909  exitNormalStc= motherNormalStc;
910  exitNormalStc.exitNormal= globalExitNormalTop;
911  }
912 
913  // Push in principle no longer necessary. G4Navigator now takes care of ...
914  // Removing this will however generate warnings for pushed particles from
915  // G4Navigator, particularly for the case of 3D replicas (Cartesian or
916  // combined Radial/Phi cases).
917  // Requires further investigation and eventually reimplementation of
918  // LevelLocate() to take into account point and direction ...
919  //
920  if ( ( !ourStep && (sampleSafety<0.5*kCarTolerance) )
921  && ( repLogical->GetSolid()->Inside(localPoint)==kSurface ) )
922  {
923  ourStep += kCarTolerance;
924  }
925 
926  if ( motherSafety<ourSafety )
927  {
928  ourSafety = motherSafety;
929  }
930 
931 #ifdef G4VERBOSE
932  if ( fCheck )
933  {
934  if( motherSolid->Inside(localPoint)==kOutside )
935  {
936  std::ostringstream message;
937  message << "Point outside volume !" << G4endl
938  << " Point " << localPoint
939  << " is outside current volume " << motherPhysical->GetName()
940  << G4endl;
941  G4double estDistToSolid= motherSolid->DistanceToIn(localPoint);
942  message << " Estimated isotropic distance to solid (distToIn)= "
943  << estDistToSolid << G4endl;
944  if( estDistToSolid > 100.0 * kCarTolerance )
945  {
946  motherSolid->DumpInfo();
947  G4Exception("G4ReplicaNavigation::ComputeStep()",
948  "GeomNav0003", FatalException, message,
949  "Point is far outside Current Volume !" );
950  }
951  else
952  G4Exception("G4ReplicaNavigation::ComputeStep()",
953  "GeomNav1002", JustWarning, message,
954  "Point is a little outside Current Volume.");
955  }
956  }
957 #endif
958 
959  // Comparison of steps may need precision protection
960  //
961 #if 1
962  if( motherDeterminedStep)
963  {
964  ourStep = motherStep;
965  exiting = true;
966  }
967 
968  // Transform it to the Grand-Mother Reference Frame (current convention)
969  //
970  if ( calculatedExitNormal )
971  {
972  if ( motherDeterminedStep )
973  {
974  exitNormalVector= motherNormalStc.exitNormal;
975  }else{
976  G4ThreeVector exitNormalGlobal= exitNormalStc.exitNormal;
977  exitNormalVector= globalToLocalTop.TransformAxis(exitNormalGlobal);
978  // exitNormalVector= globalToLocal2nd.TransformAxis(exitNormalGlobal);
979  // Alt Make it in one go to Grand-Mother, avoiding transform below
980  }
981  // Transform to Grand-mother reference frame
982  const G4RotationMatrix* rot = motherPhysical->GetRotation();
983  if ( rot )
984  {
985  exitNormalVector *= rot->inverse();
986  }
987 
988  }
989  else
990  {
991  validExitNormal = false;
992  }
993 
994 #else
995  if ( motherSafety<=ourStep )
996  {
997  if ( motherStep<=ourStep )
998  {
999  ourStep = motherStep;
1000  exiting = true;
1001  if ( validExitNormal )
1002  {
1003  const G4RotationMatrix* rot = motherPhysical->GetRotation();
1004  if ( rot )
1005  {
1006  exitNormal *= rot->inverse();
1007  }
1008  }
1009  }
1010  else
1011  {
1012  validExitNormal = false;
1013  // calculatedExitNormal= false;
1014  }
1015  }
1016 #endif
1017 
1018 
1019  G4bool daughterDeterminedStep=false;
1020  G4ThreeVector daughtNormRepCrd;
1021  // Exit normal of daughter transformed to
1022  // the coordinate system of Replica (i.e. last depth)
1023 
1024  //
1025  // Compute daughter safeties & intersections
1026  //
1027  localNoDaughters = repLogical->GetNoDaughters();
1028  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
1029  {
1030  samplePhysical = repLogical->GetDaughter(sampleNo);
1031  if ( samplePhysical!=blockedExitedVol )
1032  {
1033  G4ThreeVector localExitNorm;
1034  G4ThreeVector normReplicaCoord;
1035 
1036  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
1037  samplePhysical->GetTranslation());
1038  sampleTf.Invert();
1039  const G4ThreeVector samplePoint =
1040  sampleTf.TransformPoint(localPoint);
1041  const G4VSolid* sampleSolid =
1042  samplePhysical->GetLogicalVolume()->GetSolid();
1043  const G4double sampleSafetyDistance =
1044  sampleSolid->DistanceToIn(samplePoint);
1045  if ( sampleSafetyDistance<ourSafety )
1046  {
1047  ourSafety = sampleSafetyDistance;
1048  }
1049  if ( sampleSafetyDistance<=ourStep )
1050  {
1051  sampleDirection = sampleTf.TransformAxis(localDirection);
1052  const G4double sampleStepDistance =
1053  sampleSolid->DistanceToIn(samplePoint,sampleDirection);
1054  if ( sampleStepDistance<=ourStep )
1055  {
1056  daughterDeterminedStep= true;
1057 
1058  ourStep = sampleStepDistance;
1059  entering = true;
1060  exiting = false;
1061  *pBlockedPhysical = samplePhysical;
1062  blockedReplicaNo = -1;
1063 
1064 #ifdef DAUGHTER_NORMAL_ALSO
1065  // This norm can be calculated later, if needed daughter is available
1066  localExitNorm = sampleSolid->SurfaceNormal(samplePoint);
1067  daughtNormRepCrd = sampleTf.Inverse().TransformAxis(localExitNorm);
1068 #endif
1069 
1070 #ifdef G4VERBOSE
1071  // Check to see that the resulting point is indeed in/on volume.
1072  // This check could eventually be made only for successful candidate.
1073 
1074  if ( ( fCheck ) && ( sampleStepDistance < kInfinity ) )
1075  {
1076  G4ThreeVector intersectionPoint;
1077  intersectionPoint= samplePoint
1078  + sampleStepDistance * sampleDirection;
1079  EInside insideIntPt= sampleSolid->Inside(intersectionPoint);
1080  if ( insideIntPt != kSurface )
1081  {
1082  G4int oldcoutPrec = G4cout.precision(16);
1083  std::ostringstream message;
1084  message << "Navigator gets conflicting response from Solid."
1085  << G4endl
1086  << " Inaccurate DistanceToIn for solid "
1087  << sampleSolid->GetName() << G4endl
1088  << " Solid gave DistanceToIn = "
1089  << sampleStepDistance << " yet returns " ;
1090  if ( insideIntPt == kInside )
1091  message << "-kInside-";
1092  else if ( insideIntPt == kOutside )
1093  message << "-kOutside-";
1094  else
1095  message << "-kSurface-";
1096  message << " for this point !" << G4endl
1097  << " Point = " << intersectionPoint << G4endl;
1098  if ( insideIntPt != kInside )
1099  message << " DistanceToIn(p) = "
1100  << sampleSolid->DistanceToIn(intersectionPoint)
1101  << G4endl;
1102  if ( insideIntPt != kOutside )
1103  message << " DistanceToOut(p) = "
1104  << sampleSolid->DistanceToOut(intersectionPoint);
1105  G4Exception("G4ReplicaNavigation::ComputeStep()",
1106  "GeomNav1002", JustWarning, message);
1107  G4cout.precision(oldcoutPrec);
1108  }
1109  }
1110 #endif
1111  }
1112  }
1113  }
1114  }
1115 
1116  calculatedExitNormal &= (!daughterDeterminedStep);
1117 
1118 #ifdef DAUGHTER_NORMAL_ALSO
1119  if( daughterDeterminedStep )
1120  {
1121  // G4ThreeVector daughtNormGlobal =
1122  // GlobalToLastDepth.Inverse().TransformAxis(daughtNormRepCrd);
1123  // ==> Can calculate it, but have no way to transmit it to caller (for now)
1124 
1125  exitNormalVector=globalToLocalTop.Inverse().TransformAxis(daughtNormGlobal);
1126  validExitNormal= false; // Entering daughter - never convex for parent
1127 
1128  calculatedExitNormal= true;
1129  }
1130  // calculatedExitNormal= true; // Force it to true -- dubious
1131 #endif
1132 
1133  newSafety = ourSafety;
1134  return ourStep;
1135 }
1136 
1137 // ********************************************************************
1138 // ComputeSafety
1139 //
1140 // Compute the isotropic distance to current volume's boundaries
1141 // and to daughter volumes.
1142 // ********************************************************************
1143 //
1144 G4double
1146  const G4ThreeVector &localPoint,
1147  G4NavigationHistory &history,
1148  const G4double )
1149 {
1150  G4VPhysicalVolume *repPhysical, *motherPhysical;
1151  G4VPhysicalVolume *samplePhysical, *blockedExitedVol=0;
1152  G4LogicalVolume *repLogical;
1153  G4VSolid *motherSolid;
1154  G4ThreeVector repPoint;
1155  G4double ourSafety=kInfinity;
1156  G4double sampleSafety;
1157  G4int localNoDaughters, sampleNo;
1158  G4int depth;
1159 
1160  repPhysical = history.GetTopVolume();
1161  repLogical = repPhysical->GetLogicalVolume();
1162 
1163  //
1164  // Compute intersection with replica boundaries & replica safety
1165  //
1166 
1167  sampleSafety = DistanceToOut(history.GetTopVolume(),
1168  history.GetTopReplicaNo(),
1169  localPoint);
1170  if ( sampleSafety<ourSafety )
1171  {
1172  ourSafety = sampleSafety;
1173  }
1174 
1175  depth = history.GetDepth()-1;
1176  while ( history.GetVolumeType(depth)==kReplica )
1177  {
1178  repPoint = history.GetTransform(depth).TransformPoint(globalPoint);
1179  sampleSafety = DistanceToOut(history.GetVolume(depth),
1180  history.GetReplicaNo(depth),
1181  repPoint);
1182  if ( sampleSafety<ourSafety )
1183  {
1184  ourSafety = sampleSafety;
1185  }
1186  depth--;
1187  }
1188 
1189  // Compute mother safety & intersection
1190  //
1191  repPoint = history.GetTransform(depth).TransformPoint(globalPoint);
1192  motherPhysical = history.GetVolume(depth);
1193  motherSolid = motherPhysical->GetLogicalVolume()->GetSolid();
1194  sampleSafety = motherSolid->DistanceToOut(repPoint);
1195 
1196  if ( sampleSafety<ourSafety )
1197  {
1198  ourSafety = sampleSafety;
1199  }
1200 
1201  // Compute daughter safeties & intersections
1202  //
1203  localNoDaughters = repLogical->GetNoDaughters();
1204  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
1205  {
1206  samplePhysical = repLogical->GetDaughter(sampleNo);
1207  if ( samplePhysical!=blockedExitedVol )
1208  {
1209  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
1210  samplePhysical->GetTranslation());
1211  sampleTf.Invert();
1212  const G4ThreeVector samplePoint =
1213  sampleTf.TransformPoint(localPoint);
1214  const G4VSolid *sampleSolid =
1215  samplePhysical->GetLogicalVolume()->GetSolid();
1216  const G4double sampleSafetyDistance =
1217  sampleSolid->DistanceToIn(samplePoint);
1218  if ( sampleSafetyDistance<ourSafety )
1219  {
1220  ourSafety = sampleSafetyDistance;
1221  }
1222  }
1223  }
1224  return ourSafety;
1225 }
1226 
1227 // ********************************************************************
1228 // BackLocate
1229 // ********************************************************************
1230 //
1231 EInside
1233  const G4ThreeVector &globalPoint,
1234  G4ThreeVector &localPoint,
1235  const G4bool &exiting,
1236  G4bool &notKnownInside ) const
1237 {
1238  G4VPhysicalVolume *pNRMother=0;
1239  G4VSolid *motherSolid;
1240  G4ThreeVector repPoint, goodPoint;
1241  G4int mdepth, depth, cdepth;
1242  EInside insideCode;
1243 
1244  cdepth = history.GetDepth();
1245 
1246  // Find non replicated mother
1247  //
1248  for ( mdepth=cdepth-1; mdepth>=0; mdepth-- )
1249  {
1250  if ( history.GetVolumeType(mdepth)!=kReplica )
1251  {
1252  pNRMother = history.GetVolume(mdepth);
1253  break;
1254  }
1255  }
1256 
1257  if( pNRMother==0 )
1258  {
1259  // All the tree of mother volumes were Replicas.
1260  // This is an error, as the World volume must be a Placement
1261  //
1262  G4Exception("G4ReplicaNavigation::BackLocate()", "GeomNav0002",
1263  FatalException, "The World volume must be a Placement!");
1264  return kInside;
1265  }
1266 
1267  motherSolid = pNRMother->GetLogicalVolume()->GetSolid();
1268  goodPoint = history.GetTransform(mdepth).TransformPoint(globalPoint);
1269  insideCode = motherSolid->Inside(goodPoint);
1270  if ( (insideCode==kOutside)||((insideCode==kSurface)&&exiting) )
1271  {
1272  // Outside mother -> back up to mother level
1273  // Locate.. in Navigator will back up one more level
1274  // localPoint not required
1275  //
1276  history.BackLevel(cdepth-mdepth);
1277  // localPoint = goodPoint;
1278  }
1279  else
1280  {
1281  notKnownInside = false;
1282 
1283  // Still within replications
1284  // Check down: if on outside stop at this level
1285  //
1286  for ( depth=mdepth+1; depth<cdepth; depth++)
1287  {
1288  repPoint = history.GetTransform(depth).TransformPoint(globalPoint);
1289  insideCode = Inside(history.GetVolume(depth),
1290  history.GetReplicaNo(depth),
1291  repPoint);
1292  if ( (insideCode==kOutside)||((insideCode==kSurface)&&exiting) )
1293  {
1294  localPoint = goodPoint;
1295  history.BackLevel(cdepth-depth);
1296  return insideCode;
1297  }
1298  else
1299  {
1300  goodPoint = repPoint;
1301  }
1302  }
1303  localPoint = history.GetTransform(depth).TransformPoint(globalPoint);
1304  insideCode = Inside(history.GetVolume(depth),
1305  history.GetReplicaNo(depth),
1306  localPoint);
1307  // If outside level, set localPoint = coordinates in reference system
1308  // of *previous* level - location code in navigator will back up one
1309  // level [And also manage blocking]
1310  //
1311  if ( (insideCode==kOutside)||((insideCode==kSurface)&&exiting) )
1312  {
1313  localPoint = goodPoint;
1314  }
1315  }
1316  return insideCode;
1317 }
Definition: geomdefs.hh:54
G4String GetName() const
G4VPhysicalVolume * GetTopVolume() const
const G4ThreeVector & GetTranslation() const
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4AffineTransform Inverse() const
EInside Inside(const G4VPhysicalVolume *pVol, const G4int replicaNo, const G4ThreeVector &localPoint) const
G4int GetDepth() const
G4double GetSurfaceTolerance() const
G4VPhysicalVolume * GetDaughter(const G4int i) const
#define width
EVolume GetVolumeType(G4int n) const
G4double ComputeStep(const G4ThreeVector &globalPoint, const G4ThreeVector &globalDirection, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4bool &calculatedExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
G4double DistanceToOut(const G4VPhysicalVolume *pVol, const G4int replicaNo, const G4ThreeVector &localPoint) const
void ComputeTransformation(const G4int replicaNo, G4VPhysicalVolume *pVol, G4ThreeVector &point) const
int G4int
Definition: G4Types.hh:78
void DumpInfo() const
G4int GetTopReplicaNo() const
G4ThreeVector exitNormal
G4AffineTransform & Invert()
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
virtual EInside Inside(const G4ThreeVector &p) const =0
G4double GetRadialTolerance() const
bool G4bool
Definition: G4Types.hh:79
virtual G4ThreeVector SurfaceNormal(const G4ThreeVector &p) const =0
G4int GetReplicaNo(G4int n) const
void SetTranslation(const G4ThreeVector &v)
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
static const double kMinExitingNormalCosine
Definition: geomdefs.hh:46
G4int GetNoDaughters() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4LogicalVolume * GetLogicalVolume() const
const G4AffineTransform & GetTransform(G4int n) const
EInside
Definition: geomdefs.hh:58
G4double DistanceToOutRad(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double width, const G4double offset, const G4int replicaNo, G4ExitNormal &foundNormal) const
EAxis
Definition: geomdefs.hh:54
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
const G4RotationMatrix * GetRotation() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
G4double DistanceToOutPhi(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double width, G4ExitNormal &foundNormal) const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0
#define G4endl
Definition: G4ios.hh:61
const G4AffineTransform & GetTopTransform() const
G4double ComputeSafety(const G4ThreeVector &globalPoint, const G4ThreeVector &localPoint, G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX)
double G4double
Definition: G4Types.hh:76
Definition: geomdefs.hh:54
G4VPhysicalVolume * GetVolume(G4int n) const
void SetPhiTransformation(const G4double ang, G4VPhysicalVolume *pVol=0) const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
static const G4double d2
EInside BackLocate(G4NavigationHistory &history, const G4ThreeVector &globalPoint, G4ThreeVector &localPoint, const G4bool &exiting, G4bool &notKnownInside) const
G4double GetAngularTolerance() const
static G4GeometryTolerance * GetInstance()
G4VSolid * GetSolid() const