Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 66356 2012-12-18 09:02:32Z gcosmo $
28 //
29 //
30 // class G4NormalNavigation Implementation
31 //
32 // Author: P.Kent, 1996
33 //
34 // --------------------------------------------------------------------
35 
36 #include "G4NormalNavigation.hh"
37 #include "G4AffineTransform.hh"
38 
39 // ********************************************************************
40 // Constructor
41 // ********************************************************************
42 //
44  : fCheck(false)
45 {
46  fLogger = new G4NavigationLogger("G4NormalNavigation");
47 }
48 
49 // ********************************************************************
50 // Destructor
51 // ********************************************************************
52 //
54 {
55  delete fLogger;
56 }
57 
58 // ********************************************************************
59 // ComputeStep
60 // ********************************************************************
61 //
64  const G4ThreeVector &localDirection,
65  const G4double currentProposedStepLength,
66  G4double &newSafety,
67  G4NavigationHistory &history,
68  G4bool &validExitNormal,
69  G4ThreeVector &exitNormal,
70  G4bool &exiting,
71  G4bool &entering,
72  G4VPhysicalVolume *(*pBlockedPhysical),
73  G4int &blockedReplicaNo)
74 {
75  G4VPhysicalVolume *motherPhysical, *samplePhysical, *blockedExitedVol=0;
76  G4LogicalVolume *motherLogical;
77  G4VSolid *motherSolid;
78  G4ThreeVector sampleDirection;
79  G4double ourStep=currentProposedStepLength, motherSafety, ourSafety;
80  G4int localNoDaughters, sampleNo;
81 
82  motherPhysical = history.GetTopVolume();
83  motherLogical = motherPhysical->GetLogicalVolume();
84  motherSolid = motherLogical->GetSolid();
85 
86  // Compute mother safety
87  //
88  motherSafety = motherSolid->DistanceToOut(localPoint);
89  ourSafety = motherSafety; // Working isotropic safety
90 
91 #ifdef G4VERBOSE
92  if ( fCheck )
93  {
94  fLogger->PreComputeStepLog(motherPhysical, motherSafety, localPoint);
95  }
96 #endif
97 
98  //
99  // Compute daughter safeties & intersections
100  //
101 
102  // Exiting normal optimisation
103  //
104  if ( exiting&&validExitNormal )
105  {
106  if ( localDirection.dot(exitNormal)>=kMinExitingNormalCosine )
107  {
108  // Block exited daughter volume
109  //
110  blockedExitedVol =* pBlockedPhysical;
111  ourSafety = 0;
112  }
113  }
114  exiting = false;
115  entering = false;
116 
117  localNoDaughters = motherLogical->GetNoDaughters();
118  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo--)
119  {
120  samplePhysical = motherLogical->GetDaughter(sampleNo);
121  if ( samplePhysical!=blockedExitedVol )
122  {
123  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
124  samplePhysical->GetTranslation());
125  sampleTf.Invert();
126  const G4ThreeVector samplePoint =
127  sampleTf.TransformPoint(localPoint);
128  const G4VSolid *sampleSolid =
129  samplePhysical->GetLogicalVolume()->GetSolid();
130  const G4double sampleSafety =
131  sampleSolid->DistanceToIn(samplePoint);
132 #ifdef G4VERBOSE
133  if( fCheck )
134  {
135  fLogger->PrintDaughterLog(sampleSolid, samplePoint, sampleSafety, 0);
136  }
137 #endif
138  if ( sampleSafety<ourSafety )
139  {
140  ourSafety=sampleSafety;
141  }
142  if ( sampleSafety<=ourStep )
143  {
144  sampleDirection = sampleTf.TransformAxis(localDirection);
145  const G4double sampleStep =
146  sampleSolid->DistanceToIn(samplePoint,sampleDirection);
147 
148 #ifdef G4VERBOSE
149  if( fCheck )
150  {
151  fLogger->PrintDaughterLog(sampleSolid, samplePoint,
152  sampleSafety, sampleStep);
153  }
154 #endif
155  if ( sampleStep<=ourStep )
156  {
157  ourStep = sampleStep;
158  entering = true;
159  exiting = false;
160  *pBlockedPhysical = samplePhysical;
161  blockedReplicaNo = -1;
162 #ifdef G4VERBOSE
163  if( fCheck )
164  {
165  fLogger->AlongComputeStepLog(sampleSolid, samplePoint,
166  sampleDirection, localDirection, sampleSafety, sampleStep);
167  }
168 #endif
169  }
170  }
171  }
172  }
173  if ( currentProposedStepLength<ourSafety )
174  {
175  // Guaranteed physics limited
176  //
177  entering = false;
178  exiting = false;
179  *pBlockedPhysical = 0;
180  ourStep = kInfinity;
181  }
182  else
183  {
184  // Compute mother intersection if required
185  //
186  if ( motherSafety<=ourStep )
187  {
188  G4double motherStep = motherSolid->DistanceToOut(localPoint,
189  localDirection,
190  true,
191  &validExitNormal,
192  &exitNormal);
193 #ifdef G4VERBOSE
194  if ( fCheck )
195  {
196  fLogger->PostComputeStepLog(motherSolid, localPoint, localDirection,
197  motherStep, motherSafety);
198  }
199 #endif
200 
201  if ( motherStep<=ourStep )
202  {
203  ourStep = motherStep;
204  exiting = true;
205  entering = false;
206  if ( validExitNormal )
207  {
208  const G4RotationMatrix *rot = motherPhysical->GetRotation();
209  if (rot)
210  {
211  exitNormal *= rot->inverse();
212  }
213  }
214  }
215  else
216  {
217  validExitNormal = false;
218  }
219  }
220  }
221  newSafety = ourSafety;
222  return ourStep;
223 }
224 
225 // ********************************************************************
226 // ComputeSafety
227 // ********************************************************************
228 //
230  const G4NavigationHistory &history,
231  const G4double)
232 {
233  G4VPhysicalVolume *motherPhysical, *samplePhysical;
234  G4LogicalVolume *motherLogical;
235  G4VSolid *motherSolid;
236  G4double motherSafety, ourSafety;
237  G4int localNoDaughters, sampleNo;
238 
239  motherPhysical = history.GetTopVolume();
240  motherLogical = motherPhysical->GetLogicalVolume();
241  motherSolid = motherLogical->GetSolid();
242 
243  // Compute mother safety
244  //
245  motherSafety = motherSolid->DistanceToOut(localPoint);
246  ourSafety = motherSafety; // Working isotropic safety
247 
248 #ifdef G4VERBOSE
249  if( fCheck )
250  {
251  fLogger->ComputeSafetyLog(motherSolid, localPoint, motherSafety, true);
252  }
253 #endif
254 
255  // Compute daughter safeties
256  //
257  localNoDaughters = motherLogical->GetNoDaughters();
258  for ( sampleNo=localNoDaughters-1; sampleNo>=0; sampleNo-- )
259  {
260  samplePhysical = motherLogical->GetDaughter(sampleNo);
261  G4AffineTransform sampleTf(samplePhysical->GetRotation(),
262  samplePhysical->GetTranslation());
263  sampleTf.Invert();
264  const G4ThreeVector samplePoint =
265  sampleTf.TransformPoint(localPoint);
266  const G4VSolid *sampleSolid =
267  samplePhysical->GetLogicalVolume()->GetSolid();
268  const G4double sampleSafety =
269  sampleSolid->DistanceToIn(samplePoint);
270  if ( sampleSafety<ourSafety )
271  {
272  ourSafety = sampleSafety;
273  }
274 #ifdef G4VERBOSE
275  if(fCheck)
276  {
277  fLogger->ComputeSafetyLog(sampleSolid, samplePoint, sampleSafety, false);
278  }
279 #endif
280  }
281  return ourSafety;
282 }
void AlongComputeStepLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, const G4ThreeVector &sampleDirection, const G4ThreeVector &localDirection, G4double sampleSafety, G4double sampleStep) const
G4VPhysicalVolume * GetTopVolume() const
const G4ThreeVector & GetTranslation() const
double dot(const Hep3Vector &) const
void PrintDaughterLog(const G4VSolid *sampleSolid, const G4ThreeVector &samplePoint, G4double sampleSafety, G4double sampleStep) const
void ComputeSafetyLog(const G4VSolid *solid, const G4ThreeVector &point, G4double safety, G4bool banner) const
G4VPhysicalVolume * GetDaughter(const G4int i) const
int G4int
Definition: G4Types.hh:78
HepRotation inverse() const
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4NavigationHistory &history, const G4double pMaxLength=DBL_MAX)
G4AffineTransform & Invert()
bool G4bool
Definition: G4Types.hh:79
virtual G4double DistanceToIn(const G4ThreeVector &p, const G4ThreeVector &v) const =0
G4int GetNoDaughters() 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 PreComputeStepLog(const G4VPhysicalVolume *motherPhysical, G4double motherSafety, const G4ThreeVector &localPoint) const
const G4RotationMatrix * GetRotation() const
void PostComputeStepLog(const G4VSolid *motherSolid, const G4ThreeVector &localPoint, const G4ThreeVector &localDirection, G4double motherStep, G4double motherSafety) const
double G4double
Definition: G4Types.hh:76
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=0, G4ThreeVector *n=0) const =0
G4VSolid * GetSolid() const