Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RegularNavigation.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 G4RegularNavigation implementation
27//
28// Author: Pedro Arce (CIEMAT), May 2007
29// --------------------------------------------------------------------
30
32#include "G4TouchableHistory.hh"
34#include "G4Material.hh"
35#include "G4NormalNavigation.hh"
36#include "G4Navigator.hh"
39
40//------------------------------------------------------------------
42{
44 fMinStep = 101*kCarTolerance;
45}
46
47
48//------------------------------------------------------------------
50
51
52//------------------------------------------------------------------
54 ComputeStep(const G4ThreeVector& localPoint,
55 const G4ThreeVector& localDirection,
56 const G4double currentProposedStepLength,
57 G4double& newSafety,
58 G4NavigationHistory& history,
59 G4bool& validExitNormal,
60 G4ThreeVector& exitNormal,
61 G4bool& exiting,
62 G4bool& entering,
63 G4VPhysicalVolume *(*pBlockedPhysical),
64 G4int& blockedReplicaNo)
65{
66 // This method is never called because to be called the daughter has to be
67 // a regular structure. This would only happen if the track is in the mother
68 // of voxels volume. But the voxels fill completely their mother, so when a
69 // track enters the mother it automatically enters a voxel. Only precision
70 // problems would make this method to be called
71
72 G4ThreeVector globalPoint =
73 history.GetTopTransform().InverseTransformPoint(localPoint);
74 G4ThreeVector globalDirection =
75 history.GetTopTransform().InverseTransformAxis(localDirection);
76
77 G4ThreeVector localPoint2 = localPoint; // take away constantness
78
79 LevelLocate( history, *pBlockedPhysical, blockedReplicaNo,
80 globalPoint, &globalDirection, true, localPoint2 );
81
82
83 // Get in which voxel it is
84 //
85 G4VPhysicalVolume *motherPhysical, *daughterPhysical;
86 G4LogicalVolume *motherLogical;
87 motherPhysical = history.GetTopVolume();
88 motherLogical = motherPhysical->GetLogicalVolume();
89 daughterPhysical = motherLogical->GetDaughter(0);
90
91 auto daughterParam =
92 (G4PhantomParameterisation*)(daughterPhysical->GetParameterisation());
93 G4int copyNo = daughterParam ->GetReplicaNo(localPoint,localDirection);
94
95 G4ThreeVector voxelTranslation = daughterParam->GetTranslation( copyNo );
96 G4ThreeVector daughterPoint = localPoint - voxelTranslation;
97
98
99 // Compute step in voxel
100 //
101 return fnormalNav->ComputeStep(daughterPoint,
102 localDirection,
103 currentProposedStepLength,
104 newSafety,
105 history,
106 validExitNormal,
107 exitNormal,
108 exiting,
109 entering,
110 pBlockedPhysical,
111 blockedReplicaNo);
112}
113
114
115//------------------------------------------------------------------
117 G4ThreeVector& localPoint,
118 const G4ThreeVector& localDirection,
119 const G4double currentProposedStepLength,
120 G4double& newSafety,
121 G4NavigationHistory& history,
122 G4bool& validExitNormal,
123 G4ThreeVector& exitNormal,
124 G4bool& exiting,
125 G4bool& entering,
126 G4VPhysicalVolume *(*pBlockedPhysical),
127 G4int& blockedReplicaNo,
128 G4VPhysicalVolume* pCurrentPhysical)
129{
131
132 auto param =
133 (G4PhantomParameterisation*)(pCurrentPhysical->GetParameterisation());
134
135 if( !param->SkipEqualMaterials() )
136 {
137 return fnormalNav->ComputeStep(localPoint,
138 localDirection,
139 currentProposedStepLength,
140 newSafety,
141 history,
142 validExitNormal,
143 exitNormal,
144 exiting,
145 entering,
146 pBlockedPhysical,
147 blockedReplicaNo);
148 }
149
150
151 G4double ourStep = 0.;
152
153 // To get replica No: transform local point to the reference system of the
154 // param container volume
155 //
156 auto ide = (G4int)history.GetDepth();
157 G4ThreeVector containerPoint = history.GetTransform(ide)
158 .InverseTransformPoint(localPoint);
159
160 // Point in global frame
161 //
162 containerPoint = history.GetTransform(ide).InverseTransformPoint(localPoint);
163
164 // Point in voxel parent volume frame
165 //
166 containerPoint = history.GetTransform(ide-1).TransformPoint(containerPoint);
167
168 // Store previous voxel translation to move localPoint by the difference
169 // with the new one
170 //
171 G4ThreeVector prevVoxelTranslation = containerPoint - localPoint;
172
173 // Do not use the expression below: There are cases where the
174 // fLastLocatedPointLocal does not give the correct answer
175 // (particle reaching a wall and bounced back, particle travelling through
176 // the wall that is deviated in an step, ...; these are pathological cases
177 // that give wrong answers in G4PhantomParameterisation::GetReplicaNo()
178 //
179 // G4ThreeVector prevVoxelTranslation = param->GetTranslation( copyNo );
180
181 G4int copyNo = param->GetReplicaNo(containerPoint,localDirection);
182
183 G4Material* currentMate = param->ComputeMaterial( copyNo, nullptr, nullptr );
184 G4VSolid* voxelBox = pCurrentPhysical->GetLogicalVolume()->GetSolid();
185
186 G4VSolid* containerSolid = param->GetContainerSolid();
187 G4Material* nextMate;
188 G4bool bFirstStep = true;
189 G4double newStep;
190 G4double totalNewStep = 0.;
191
192 // Loop while same material is found
193 //
194 //
195 fNumberZeroSteps = 0;
196 for( G4int ii = 0; ii < fNoStepsAllowed+1; ++ii )
197 {
198 if( ii == fNoStepsAllowed ) {
199 // Must kill this stuck track
200 //
201 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
202 .InverseTransformPoint(localPoint);
203 std::ostringstream message;
204 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials()"
205 << "Stuck Track: potential geometry or navigation problem."
206 << G4endl
207 << " Track stuck, moving for more than "
208 << ii << " steps" << G4endl
209 << "- at point " << pGlobalpoint << G4endl
210 << " local direction: " << localDirection << G4endl;
211 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
212 "GeomRegNav1001",
214 message);
215 }
216 newStep = voxelBox->DistanceToOut( localPoint, localDirection );
217 fLastStepWasZero = (newStep<fMinStep);
218 if( fLastStepWasZero )
219 {
220 ++fNumberZeroSteps;
221#ifdef G4DEBUG_NAVIGATION
222 if( fNumberZeroSteps > 1 )
223 {
224 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
225 .InverseTransformPoint(localPoint);
226 std::ostringstream message;
227 message.precision(16);
228 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials(): another 'zero' step, # "
229 << fNumberZeroSteps
230 << ", at " << pGlobalpoint
231 << ", nav-comp-step calls # " << ii
232 << ", Step= " << newStep;
233 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
234 "GeomRegNav1002", JustWarning, message,
235 "Potential overlap in geometry!");
236 }
237#endif
238 if( fNumberZeroSteps > fActionThreshold_NoZeroSteps-1 )
239 {
240 // Act to recover this stuck track. Pushing it along direction
241 //
242 newStep = std::min(101*kCarTolerance*std::pow(10,fNumberZeroSteps-2),0.1);
243#ifdef G4DEBUG_NAVIGATION
244 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
245 .InverseTransformPoint(localPoint);
246 std::ostringstream message;
247 message.precision(16);
248 message << "Track stuck or not moving." << G4endl
249 << " Track stuck, not moving for "
250 << fNumberZeroSteps << " steps" << G4endl
251 << "- at point " << pGlobalpoint
252 << " (local point " << localPoint << ")" << G4endl
253 << " local direction: " << localDirection
254 << " Potential geometry or navigation problem !"
255 << G4endl
256 << " Trying pushing it of " << newStep << " mm ...";
257 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
258 "GeomRegNav1003", JustWarning, message,
259 "Potential overlap in geometry!");
260#endif
261 }
262 if( fNumberZeroSteps > fAbandonThreshold_NoZeroSteps-1 )
263 {
264 // Must kill this stuck track
265 //
266 G4ThreeVector pGlobalpoint = history.GetTransform(ide)
267 .InverseTransformPoint(localPoint);
268 std::ostringstream message;
269 message << "G4RegularNavigation::ComputeStepSkippingEqualMaterials()"
270 << "Stuck Track: potential geometry or navigation problem."
271 << G4endl
272 << " Track stuck, not moving for "
273 << fNumberZeroSteps << " steps" << G4endl
274 << "- at point " << pGlobalpoint << G4endl
275 << " local direction: " << localDirection << G4endl;
276 G4Exception("G4RegularNavigation::ComputeStepSkippingEqualMaterials()",
277 "GeomRegNav1004",
279 message);
280 }
281 }
282 else
283 {
284 // reset the zero step counter when a non-zero step was performed
285 fNumberZeroSteps = 0;
286 }
287 if( (bFirstStep) && (newStep < currentProposedStepLength) )
288 {
289 exiting = true;
290 }
291 bFirstStep = false;
292
293 newStep += kCarTolerance; // Avoid precision problems
294 ourStep += newStep;
295 totalNewStep += newStep;
296
297 // Physical process is limiting the step, don't continue
298 //
299 if(std::fabs(totalNewStep-currentProposedStepLength) < kCarTolerance)
300 {
301 return currentProposedStepLength;
302 }
303 if(totalNewStep > currentProposedStepLength)
304 {
306 AddStepLength(copyNo, newStep-totalNewStep+currentProposedStepLength);
307 return currentProposedStepLength;
308 }
310
311 // Move container point until wall of voxel
312 //
313 containerPoint += newStep*localDirection;
314 if( containerSolid->Inside( containerPoint ) != kInside )
315 {
316 break;
317 }
318
319 // Get copyNo and translation of new voxel
320 //
321 copyNo = param->GetReplicaNo(containerPoint, localDirection);
322 G4ThreeVector voxelTranslation = param->GetTranslation( copyNo );
323
324 // Move local point until wall of voxel and then put it in the new voxel
325 // local coordinates
326 //
327 localPoint += newStep*localDirection;
328 localPoint += prevVoxelTranslation - voxelTranslation;
329
330 prevVoxelTranslation = voxelTranslation;
331
332 // Check if material of next voxel is the same as that of the current voxel
333 nextMate = param->ComputeMaterial( copyNo, nullptr, nullptr );
334
335 if( currentMate != nextMate ) { break; }
336 }
337
338 return ourStep;
339}
340
341
342//------------------------------------------------------------------
345 const G4NavigationHistory& history,
346 const G4double pMaxLength)
347{
348 // This method is never called because to be called the daughter has to be a
349 // regular structure. This would only happen if the track is in the mother of
350 // voxels volume. But the voxels fill completely their mother, so when a
351 // track enters the mother it automatically enters a voxel. Only precision
352 // problems would make this method to be called
353
354 // Compute step in voxel
355 //
356 return fnormalNav->ComputeSafety(localPoint,
357 history,
358 pMaxLength );
359}
360
361
362//------------------------------------------------------------------
363G4bool
365 const G4VPhysicalVolume* ,
366 const G4int ,
367 const G4ThreeVector& globalPoint,
368 const G4ThreeVector* globalDirection,
369 const G4bool, // pLocatedOnEdge,
370 G4ThreeVector& localPoint )
371{
372 G4VPhysicalVolume *motherPhysical, *pPhysical;
374 G4LogicalVolume *motherLogical;
375 G4ThreeVector localDir;
376 G4int replicaNo;
377
378 motherPhysical = history.GetTopVolume();
379 motherLogical = motherPhysical->GetLogicalVolume();
380
381 pPhysical = motherLogical->GetDaughter(0);
382 pParam = (G4PhantomParameterisation*)(pPhysical->GetParameterisation());
383
384 // Save parent history in touchable history
385 // ... for use as parent t-h in ComputeMaterial method of param
386 //
387 G4TouchableHistory parentTouchable( history );
388
389 // Get local direction
390 //
391 if( globalDirection != nullptr )
392 {
393 localDir = history.GetTopTransform().TransformAxis(*globalDirection);
394 }
395 else
396 {
397 localDir = G4ThreeVector(0.,0.,0.);
398 }
399
400 // Enter this daughter
401 //
402 replicaNo = pParam->GetReplicaNo( localPoint, localDir );
403
404 if( replicaNo < 0 || replicaNo >= G4int(pParam->GetNoVoxels()) )
405 {
406 return false;
407 }
408
409 // Set the correct copy number in physical
410 //
411 pPhysical->SetCopyNo(replicaNo);
412 pParam->ComputeTransformation(replicaNo,pPhysical);
413
414 history.NewLevel(pPhysical, kParameterised, replicaNo );
415 localPoint = history.GetTopTransform().TransformPoint(globalPoint);
416
417 // Set the correct solid and material in Logical Volume
418 //
419 G4LogicalVolume *pLogical = pPhysical->GetLogicalVolume();
420
421 pLogical->UpdateMaterial(pParam->ComputeMaterial(replicaNo,
422 pPhysical, &parentTouchable) );
423 return true;
424}
425
426//------------------------------------------------------------------
427void
429{
430 fnormalNav = fnormnav;
431}
@ JustWarning
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4ThreeVector InverseTransformAxis(const G4ThreeVector &axis) const
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
G4ThreeVector TransformAxis(const G4ThreeVector &axis) const
G4ThreeVector InverseTransformPoint(const G4ThreeVector &vec) const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4LogicalVolume represents a leaf node or unpositioned subtree in the geometry hierarchy....
G4VSolid * GetSolid() const
G4VPhysicalVolume * GetDaughter(const std::size_t i) const
void UpdateMaterial(G4Material *pMaterial)
G4NavigationHistory is a class responsible for the maintenance of the history of the path taken throu...
void NewLevel(G4VPhysicalVolume *pNewMother, EVolume vType=kNormal, G4int nReplica=-1)
const G4AffineTransform & GetTopTransform() const
std::size_t GetDepth() const
G4VPhysicalVolume * GetTopVolume() const
const G4AffineTransform & GetTransform(G4int n) const
G4NormalNavigation is a concrete utility class for navigation in volumes containing only G4PVPlacemen...
G4PhantomParameterisation describes regular parameterisations: a set of boxes of equal dimension in t...
virtual G4int GetReplicaNo(const G4ThreeVector &localPoint, const G4ThreeVector &localDir)
G4Material * ComputeMaterial(const G4int repNo, G4VPhysicalVolume *currentVol, const G4VTouchable *parentTouch=nullptr) override
void ComputeTransformation(const G4int, G4VPhysicalVolume *) const override
std::size_t GetNoVoxels() const
static G4RegularNavigationHelper * Instance()
void AddStepLength(G4int copyNo, G4double slen)
G4double ComputeStepSkippingEqualMaterials(G4ThreeVector &localPoint, const G4ThreeVector &localDirection, const G4double currentProposedStepLength, G4double &newSafety, G4NavigationHistory &history, G4bool &validExitNormal, G4ThreeVector &exitNormal, G4bool &exiting, G4bool &entering, G4VPhysicalVolume *(*pBlockedPhysical), G4int &blockedReplicaNo, G4VPhysicalVolume *pCurrentPhysical)
~G4RegularNavigation() override
G4bool LevelLocate(G4NavigationHistory &history, const G4VPhysicalVolume *blockedVol, const G4int blockedNum, const G4ThreeVector &globalPoint, const G4ThreeVector *globalDirection, const G4bool pLocatedOnEdge, G4ThreeVector &localPoint) 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
void SetNormalNavigation(G4NormalNavigation *fnormnav)
G4double ComputeSafety(const G4ThreeVector &localPoint, const G4NavigationHistory &history, const G4double pProposedMaxLength=DBL_MAX) final
G4TouchableHistory is an object representing a touchable detector element, and its history in the geo...
G4VPhysicalVolume is an abstract base class for the representation of a positioned volume....
virtual void SetCopyNo(G4int CopyNo)=0
G4LogicalVolume * GetLogicalVolume() const
virtual G4VPVParameterisation * GetParameterisation() const =0
G4VSolid is an abstract base class for solids, physical shapes that can be tracked through....
Definition G4VSolid.hh:80
virtual EInside Inside(const G4ThreeVector &p) const =0
virtual G4double DistanceToOut(const G4ThreeVector &p, const G4ThreeVector &v, const G4bool calcNorm=false, G4bool *validNorm=nullptr, G4ThreeVector *n=nullptr) const =0
@ kInside
Definition geomdefs.hh:70
@ kParameterised
Definition geomdefs.hh:86