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