Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RegularNavigation.hh
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
27//
28// Class description:
29//
30// Utility for fast navigation in volumes containing a regular
31// parameterisation. If two contiguous voxels have the same material,
32// navigation does not stop at the surface
33
34// History:
35// - Created. P. Arce, May 2007
36// --------------------------------------------------------------------
37#ifndef G4RegularNavigation_HH
38#define G4RegularNavigation_HH
39
40#include <vector>
41
42#include "G4Types.hh"
43#include "G4ThreeVector.hh"
44
47class G4Navigator;
49
51{
52 public: // with description
53
56
58 const G4VPhysicalVolume* blockedVol,
59 const G4int blockedNum,
60 const G4ThreeVector& globalPoint,
61 const G4ThreeVector* globalDirection,
62 const G4bool pLocatedOnEdge,
63 G4ThreeVector& localPoint );
64 // Locate point using its position with respect to regular
65 // parameterisation container volume.
66
67 G4double ComputeStep( const G4ThreeVector& globalPoint,
68 const G4ThreeVector& globalDirection,
69 const G4double currentProposedStepLength,
70 G4double& newSafety,
71 G4NavigationHistory& history,
72 G4bool& validExitNormal,
73 G4ThreeVector& exitNormal,
74 G4bool& exiting,
75 G4bool& entering,
76 G4VPhysicalVolume *(*pBlockedPhysical),
77 G4int& blockedReplicaNo );
78 // Method never called because to be called the daughter has to be a
79 // 'regular' volume. This would only happen if the track is in the
80 // mother of voxels volume. But the voxels fill completely their mother,
81 // so when a track enters the mother it automatically enters a voxel.
82
84 G4ThreeVector& localPoint,
85 const G4ThreeVector& globalDirection,
86 const G4double currentProposedStepLength,
87 G4double& newSafety,
88 G4NavigationHistory& history,
89 G4bool& validExitNormal,
90 G4ThreeVector& exitNormal,
91 G4bool& exiting,
92 G4bool& entering,
93 G4VPhysicalVolume *(*pBlockedPhysical),
94 G4int& blockedReplicaNo,
95 G4VPhysicalVolume* pCurrentPhysical);
96 // Compute the step skipping surfaces when they separate voxels with
97 // equal materials. Loop to voxels until a different material is found:
98 // invokes G4NormalNavigation::ComputeStep() in each voxel and move the
99 // point to the next voxel.
100
101 G4double ComputeSafety( const G4ThreeVector& localPoint,
102 const G4NavigationHistory& history,
103 const G4double pProposedMaxLength = DBL_MAX );
104 // Method never called because to be called the daughter has to be a
105 // 'regular' volume. This would only happen if the track is in the
106 // mother of voxels volume. But the voxels fill completely their mother,
107 // so when a track enters the mother it automatically enters a voxel.
108
109 public: // without description
110
111 // Set and Get methods
112
113 void SetVerboseLevel(G4int level) { fverbose = level; }
114 void CheckMode(G4bool mode) { fcheck = mode; }
116 { fnormalNav = fnormnav; }
117
118 private:
119
120 G4int fverbose = false;
121 G4bool fcheck = false;
122
123 G4NormalNavigation* fnormalNav = nullptr;
124 G4double kCarTolerance;
125 G4double fMinStep;
126
127 G4bool fLastStepWasZero = false;
128 // Whether the last ComputeStep moved Zero. Used to check for edges.
129 G4int fNumberZeroSteps = 0;
130 // Number of preceding moves that were Zero. Reset to 0 after finite step
131 G4int fActionThreshold_NoZeroSteps = 2;
132 // After this many failed/zero steps, act (push etc)
133 G4int fAbandonThreshold_NoZeroSteps = 25;
134 // After this many failed/zero steps, abandon track
135 G4int fNoStepsAllowed = 10000;
136 // Maximum number of steps a track can travel skipping voxels (if there are more, track is assumed to be stuck and it is killed)
137};
138
139#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetVerboseLevel(G4int level)
void CheckMode(G4bool mode)
void SetNormalNavigation(G4NormalNavigation *fnormnav)
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)
#define DBL_MAX
Definition: templates.hh:62