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