Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorPropagationNavigator.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 G4ErrorPropagationNavigator implementation
27//
28// Author: Pedro Arce, CIEMAT
29// --------------------------------------------------------------------
30
32
33#include "globals.hh"
34#include "G4ThreeVector.hh"
37
40
42ComputeStep ( const G4ThreeVector &pGlobalPoint,
43 const G4ThreeVector &pDirection,
44 const G4double pCurrentProposedStepLength,
45 G4double &pNewSafety )
46{
47 G4double safetyGeom = DBL_MAX;
48
49 G4double Step = G4Navigator::ComputeStep(pGlobalPoint, pDirection,
50 pCurrentProposedStepLength,
51 safetyGeom);
52
55
56 if ( g4edata != nullptr )
57 {
58 const G4ErrorTarget* target = g4edata->GetTarget();
59 if( target != nullptr )
60 {
61 G4double StepPlane=target->GetDistanceFromPoint(pGlobalPoint,pDirection);
62
63 if( StepPlane < 0. ) // Negative means target is crossed,
64 { // will not be found
65 StepPlane = DBL_MAX;
66 }
67#ifdef G4VERBOSE
69 {
70 G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
71 << " Target step: " << StepPlane
72 << ", Transportation step: " << Step << G4endl;
73 target->Dump( "G4ErrorPropagationNavigator::ComputeStep Target " );
74 }
75#endif
76
77 if( StepPlane < Step )
78 {
79#ifdef G4VERBOSE
81 {
82 G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
83 << " TargetCloserThanBoundary: " << StepPlane << " < "
84 << Step << G4endl;
85 }
86#endif
87 Step = StepPlane;
89 }
90 else
91 {
93 }
94 }
95 }
96 G4double safetyTarget = TargetSafetyFromPoint(pGlobalPoint);
97
98 // Avoid call to G4Navigator::ComputeSafety - which could have side effects
99 //
100 pNewSafety = std::min(safetyGeom, safetyTarget);
101
102#ifdef G4VERBOSE
104 {
105 G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
106 << " Step: " << Step << ", ComputeSafety: " << pNewSafety
107 << G4endl;
108 }
109#endif
110
111 return Step;
112}
113
114//-------------------------------------------------------------------
115
117TargetSafetyFromPoint( const G4ThreeVector& pGlobalpoint )
118{
119 G4double safety = DBL_MAX;
120
121 G4ErrorPropagatorData* g4edata
123
124 if ( g4edata != nullptr )
125 {
126 const G4ErrorTarget* target = g4edata->GetTarget();
127 if( target != nullptr )
128 {
129 safety = target->GetDistanceFromPoint(pGlobalpoint);
130 }
131 }
132 return safety;
133}
134
135//-------------------------------------------------------------------
136
138ComputeSafety( const G4ThreeVector &pGlobalPoint,
139 const G4double pMaxLength,
140 const G4bool keepState )
141{
142 G4double safetyGeom = G4Navigator::ComputeSafety(pGlobalPoint,
143 pMaxLength, keepState);
144
145 G4double safetyTarget = TargetSafetyFromPoint( pGlobalPoint );
146
147 return std::min(safetyGeom, safetyTarget);
148}
149
150//-------------------------------------------------------------------
151
153GetGlobalExitNormal( const G4ThreeVector& point, G4bool* valid )
154{
155 G4ErrorPropagatorData* g4edata
157 const G4ErrorTarget* target = nullptr;
158
159 G4ThreeVector normal(0.0, 0.0, 0.0);
160 G4double distance= 0;
161
162 // Determine which 'geometry' limited the step
163 if ( g4edata != nullptr )
164 {
165 target = g4edata->GetTarget();
166 if( target != nullptr )
167 {
168 distance = target->GetDistanceFromPoint(point);
169 }
170 }
171
172 if( distance > kCarTolerance || (target == nullptr) )
173 // Not reached the target or if a target does not exist,
174 // this seems the best we can do
175 {
176 normal = G4Navigator::GetGlobalExitNormal(point, valid);
177 }
178 else
179 {
180 switch( target->GetType() )
181 {
183 // The volume is in the 'real' mass geometry
184 normal = G4Navigator::GetGlobalExitNormal(point, valid);
185 break;
187 normal = G4ThreeVector( 0.0, 0.0, 0.0);
188 *valid = false;
189 G4Exception("G4ErrorPropagationNavigator::GetGlobalExitNormal",
190 "Geometry1003",
191 JustWarning, "Unexpected value of Target type");
192 break;
195 const auto surfTarget= static_cast<const G4ErrorSurfaceTarget*>(target);
196 normal = surfTarget->GetTangentPlane(point).normal().unit();
197 *valid = true;
198 break;
199
200// default:
201// normal= G4ThreeVector( 0.0, 0.0, 0.0 );
202// *valid = false;
203// G4Exception("G4ErrorPropagationNavigator::GetGlobalExitNormal",
204// "Geometry:003",
205// FatalException, "Impossible value of Target type");
206// break;
207 }
208 }
209 return normal;
210}
211
@ G4ErrorState_TargetCloserThanBoundary
@ G4ErrorState_Propagating
@ G4ErrorTarget_PlaneSurface
@ G4ErrorTarget_CylindricalSurface
@ G4ErrorTarget_GeomVolume
@ G4ErrorTarget_TrkL
@ JustWarning
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
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector unit() const
G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety) override
G4double TargetSafetyFromPoint(const G4ThreeVector &pGlobalpoint)
G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid) override
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true) override
static G4ErrorPropagatorData * GetErrorPropagatorData()
const G4ErrorTarget * GetTarget(G4bool mustExist=false) const
void SetState(G4ErrorState sta)
virtual void Dump(const G4String &msg) const =0
G4ErrorTargetType GetType() const
virtual G4double GetDistanceFromPoint(const G4ThreeVector &, const G4ThreeVector &) const
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
G4double kCarTolerance
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
#define DBL_MAX
Definition templates.hh:62