Geant4  10.02.p01
G4NormalNavigation.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: G4NormalNavigation.cc 93289 2015-10-15 10:01:15Z gcosmo $
28 //
29 //
30 // class G4NormalNavigation Implementation
31 //
32 // Author: P.Kent, 1996
33 //
34 // --------------------------------------------------------------------
35 
36 #include "G4NormalNavigation.hh"
37 #include "G4NavigationLogger.hh"
38 #include "G4AffineTransform.hh"
39 
40 // ********************************************************************
41 // Constructor
42 // ********************************************************************
43 //
45  : fCheck(false)
46 {
47  fLogger = new G4NavigationLogger("G4NormalNavigation");
48 }
49 
50 // ********************************************************************
51 // Destructor
52 // ********************************************************************
53 //
55 {
56  delete fLogger;
57 }
58 
59 // ********************************************************************
60 // ComputeStep
61 // ********************************************************************
62 //
63 // On entry
64 // exitNormal, validExitNormal: for previous exited volume (daughter)
65 //
66 // On exit
67 // exitNormal, validExitNormal: for mother, if exiting it (else unchanged)
70  const G4ThreeVector &localDirection,
71  const G4double currentProposedStepLength,
72  G4double &newSafety,
73  G4NavigationHistory &history,
74  G4bool &validExitNormal,
75  G4ThreeVector &exitNormal,
76  G4bool &exiting,
77  G4bool &entering,
78  G4VPhysicalVolume *(*pBlockedPhysical),
79  G4int &blockedReplicaNo)
80 {
81  G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=0;
82  G4LogicalVolume *motherLogical;
83  G4VSolid *motherSolid;
84  G4ThreeVector sampleDirection;
85  G4double ourStep=currentProposedStepLength, ourSafety;
86  G4double motherSafety, motherStep=DBL_MAX;
87  G4int localNoDaughters, sampleNo;
88  G4bool motherValidExitNormal=false;
89  G4ThreeVector motherExitNormal;
90 
91  motherPhysical = history.GetTopVolume();
92  motherLogical = motherPhysical->GetLogicalVolume();
93  motherSolid = motherLogical->GetSolid();
94 
95  // Compute mother safety
96  //
97  motherSafety = motherSolid->DistanceToOut(localPoint);
98  ourSafety = motherSafety; // Working isotropic safety
99 
100  localNoDaughters = motherLogical->GetNoDaughters();
101 
102 #ifdef G4VERBOSE
103  if ( fCheck && ( (localNoDaughters>0) || (ourStep < motherSafety) ) )
104  {
105  fLogger->PreComputeStepLog(motherPhysical, motherSafety, localPoint);
106  }
107 #endif
108  // Compute daughter safeties & intersections
109  //
110 
111  // Exiting normal optimisation
112  //
113  if ( exiting&&validExitNormal )
114  {
115  if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
116  {
117  // Block exited daughter volume
118  //
119  blockedExitedVol = (*pBlockedPhysical);
120  ourSafety = 0;
121  }
122  }
123  exiting = false;
124  entering = false;
125 
126 #ifdef G4VERBOSE
127  if ( fCheck )
128  {
129  // Compute early:
130  // a) to check whether point is (wrongly) outside
131  // (signaled if step < 0 or step == kInfinity )
132  // b) to check value against answer of daughters!
133 
134  motherStep = motherSolid->DistanceToOut(localPoint,
135  localDirection,
136  true,
137  &motherValidExitNormal,
138  &motherExitNormal);
139 
140  if( (motherStep >= kInfinity) || (motherStep < 0.0) )
141  {
142  // Error - indication of being outside solid !!
143  fLogger->ReportOutsideMother(localPoint, localDirection, motherPhysical);
144 
145  ourStep = motherStep = 0.0;
146 
147  exiting= true;
148  entering= false;
149 
150  // If we are outside the solid does the normal make sense?
151  validExitNormal= motherValidExitNormal;
152  exitNormal= motherExitNormal;
153 
154  *pBlockedPhysical= 0; // or motherPhysical ?
155  blockedReplicaNo= 0; // or motherReplicaNumber ?
156 
157  newSafety= 0.0;
158  return ourStep;
159  }
160  }
161 #endif
162 
163  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo--)
164  {
165  samplePhysical = motherLogical->GetDaughter(sampleNo);
166  if ( samplePhysical!=blockedExitedVol )
167  {
168  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
169  samplePhysical->GetTranslation());
170  sampleTf.Invert();
171  const G4ThreeVector samplePoint = sampleTf.TransformPoint(localPoint);
172  const G4VSolid *sampleSolid =
173  samplePhysical->GetLogicalVolume()->GetSolid();
174  const G4double sampleSafety =
175  sampleSolid->DistanceToIn(samplePoint);
176 
177  if ( sampleSafety<ourSafety )
178  {
179  ourSafety=sampleSafety;
180  }
181 
182  if ( sampleSafety<=ourStep )
183  {
184  sampleDirection = sampleTf.TransformAxis(localDirection);
185  const G4double sampleStep =
186  sampleSolid->DistanceToIn(samplePoint,sampleDirection);
187 #ifdef G4VERBOSE
188  if( fCheck )
189  {
190  fLogger->PrintDaughterLog(sampleSolid, samplePoint,
191  sampleSafety, true,
192  sampleDirection, sampleStep);
193  }
194 #endif
195  if ( sampleStep<=ourStep )
196  {
197  ourStep = sampleStep;
198  entering = true;
199  exiting = false;
200  *pBlockedPhysical = samplePhysical;
201  blockedReplicaNo = -1;
202 #ifdef G4VERBOSE
203  if( fCheck )
204  {
205  fLogger->AlongComputeStepLog(sampleSolid, samplePoint,
206  sampleDirection, localDirection, sampleSafety, sampleStep);
207  }
208 #endif
209  }
210 
211 #ifdef G4VERBOSE
212  if( fCheck && (sampleStep < kInfinity) && (sampleStep >= motherStep) )
213  {
214  // The intersection point with the daughter is at or after the exit
215  // point from the mother volume. Double check!
216  fLogger->CheckDaughterEntryPoint(sampleSolid,
217  samplePoint, sampleDirection,
218  motherSolid,
219  localPoint, localDirection,
220  motherStep, sampleStep);
221  }
222 #endif
223  } // end of if ( sampleSafety <= ourStep )
224 #ifdef G4VERBOSE
225  else if( fCheck )
226  {
227  fLogger->PrintDaughterLog(sampleSolid, samplePoint,
228  sampleSafety, false,
229  G4ThreeVector(0.,0.,0.), -1.0 );
230  }
231 #endif
232  }
233  }
234  if ( currentProposedStepLength<ourSafety )
235  {
236  // Guaranteed physics limited
237  //
238  entering = false;
239  exiting = false;
240  *pBlockedPhysical = 0;
241  ourStep = kInfinity;
242  }
243  else
244  {
245  // Consider intersection with mother solid
246  //
247  if ( motherSafety<=ourStep )
248  {
249  if ( !fCheck ) // The call is moved above when running in check_mode
250  {
251  motherStep = motherSolid->DistanceToOut(localPoint,
252  localDirection,
253  true,
254  &motherValidExitNormal,
255  &motherExitNormal);
256  }
257 #ifdef G4VERBOSE
258  else // check_mode
259  {
260  fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
261  motherStep, motherSafety);
262  if( motherValidExitNormal )
263  {
264  fLogger->CheckAndReportBadNormal(motherExitNormal,
265  localPoint,
266  localDirection,
267  motherStep,
268  motherSolid,
269  "From motherSolid::DistanceToOut" );
270  }
271  }
272 #endif
273 
274  if( (motherStep >= kInfinity) || (motherStep < 0.0) )
275  {
276  // Clearly outside the mother solid!
277 #ifdef G4VERBOSE
278  fLogger->ReportOutsideMother(localPoint, localDirection, motherPhysical);
279 #endif
280  ourStep = motherStep = 0.0;
281  exiting = true;
282  entering = false;
283  // validExitNormal= motherValidExitNormal;
284  // exitNormal= motherExitNormal;
285  // The normal could be useful - but only if near the mother
286  // But it could be unreliable!
287  validExitNormal = false;
288  *pBlockedPhysical= 0; // or motherPhysical ?
289  blockedReplicaNo= 0; // or motherReplicaNumber ?
290  newSafety= 0.0;
291  return ourStep;
292  }
293 
294  if ( motherStep<=ourStep )
295  {
296  ourStep = motherStep;
297  exiting = true;
298  entering = false;
299  validExitNormal= motherValidExitNormal;
300  exitNormal= motherExitNormal;
301 
302  if ( motherValidExitNormal )
303  {
304  const G4RotationMatrix *rot = motherPhysical->GetRotation();
305  if (rot)
306  {
307  exitNormal *= rot->inverse();
308 #ifdef G4VERBOSE
309  if( fCheck )
310  fLogger->CheckAndReportBadNormal(exitNormal, // rotated
311  motherExitNormal, // original
312  *rot,
313  "From RotationMatrix" );
314 #endif
315  }
316  }
317  }
318  else
319  {
320  validExitNormal = false;
321  }
322  }
323  }
324  newSafety = ourSafety;
325  return ourStep;
326 }
327 
328 // ********************************************************************
329 // ComputeSafety
330 // ********************************************************************
331 //
333  const G4NavigationHistory &history,
334  const G4double)
335 {
336  G4VPhysicalVolume *motherPhysical, *samplePhysical;
337  G4LogicalVolume *motherLogical;
338  G4VSolid *motherSolid;
339  G4double motherSafety, ourSafety;
340  G4int localNoDaughters, sampleNo;
341 
342  motherPhysical = history.GetTopVolume();
343  motherLogical = motherPhysical->GetLogicalVolume();
344  motherSolid = motherLogical->GetSolid();
345 
346  // Compute mother safety
347  //
348  motherSafety = motherSolid->DistanceToOut(localPoint);
349  ourSafety = motherSafety; // Working isotropic safety
350 
351 #ifdef G4VERBOSE
352  if( fCheck )
353  {
354  fLogger->ComputeSafetyLog(motherSolid,localPoint,motherSafety,true,true);
355  }
356 #endif
357 
358  // Compute daughter safeties
359  //
360  localNoDaughters = motherLogical->GetNoDaughters();
361  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
362  {
363  samplePhysical = motherLogical->GetDaughter(sampleNo);
364  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
365  samplePhysical->GetTranslation());
366  sampleTf.Invert();
367  const G4ThreeVector samplePoint =
368  sampleTf.TransformPoint(localPoint);
369  const G4VSolid *sampleSolid =
370  samplePhysical->GetLogicalVolume()->GetSolid();
371  const G4double sampleSafety =
372  sampleSolid->DistanceToIn(samplePoint);
373  if ( sampleSafety<ourSafety )
374  {
375  ourSafety = sampleSafety;
376  }
377 #ifdef G4VERBOSE
378  if(fCheck)
379  {
380  fLogger->ComputeSafetyLog(sampleSolid,samplePoint,
381  sampleSafety,false,false);
382  // Not mother, no banner
383  }
384 #endif
385  }
386  return ourSafety;
387 }
388 
389 // The following methods have been imported to this source file
390 // in order to avoid dependency of the header file on the
391 // header implementation of G4NavigationLogger.
392 
393 // ********************************************************************
394 // GetVerboseLevel
395 // ********************************************************************
396 //
398 {
399  return fLogger->GetVerboseLevel();
400 }
401 
402 // ********************************************************************
403 // SetVerboseLevel
404 // ********************************************************************
405 //
407 {
408  fLogger->SetVerboseLevel(level);
409 }
410 
void SetVerboseLevel(G4int level)
void AlongComputeStepLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4ThreeVector &localDirection, G4double sampleSafety, G4double sampleStep) const
G4VPhysicalVolume * GetTopVolume() const
void PrintDaughterLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, G4double sampleSafety, G4bool onlySafety, const G4ThreeVector &sampleDirection, G4double sampleStep) const
const G4ThreeVector & GetTranslation() const
static const G4double kInfinity
Definition: geomdefs.hh:42
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
G4VPhysicalVolume * GetDaughter(const G4int i) const
int G4int
Definition: G4Types.hh:78
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
void ReportOutsideMother(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4VPhysicalVolume *motherPV, G4double tDist=30.0 *CLHEP::cm) const
G4AffineTransform & Invert()
bool G4bool
Definition: G4Types.hh:79
G4bool CheckAndReportBadNormal(const G4ThreeVector &unitNormal, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double step, const G4VSolid *solid, const char *msg) const
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
static const double kMinExitingNormalCosine
Definition: geomdefs.hh:46
G4int GetNoDaughters() const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4double ComputeStep(const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo)
G4LogicalVolume * GetLogicalVolume() const
void ComputeSafetyLog(const G4VSolid *solid, const G4ThreeVector &point, G4double safety, G4bool isMotherVolume, G4int banner=-1) const
void PreComputeStepLog(const G4VPhysicalVolume *motherPhysical, G4double motherSafety, const G4ThreeVector &localPoint) const
const G4RotationMatrix * GetRotation() const
G4int GetVerboseLevel() const
G4NavigationLogger * fLogger
void PostComputeStepLog(const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double motherSafety) const
void CheckDaughterEntryPoint(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double sampleStep) const
double G4double
Definition: G4Types.hh:76
G4int GetVerboseLevel() const
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
#define DBL_MAX
Definition: templates.hh:83
void SetVerboseLevel(G4int level)
G4VSolid * GetSolid() const