Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 //
27 // $Id: G4ErrorPropagationNavigator.cc 87697 2014-12-17 09:40:21Z gcosmo $
28 //
29 //
30 // --------------------------------------------------------------------
31 // GEANT 4 class implementation file
32 // --------------------------------------------------------------------
33 
35 
36 #include "globals.hh"
37 #include "G4ThreeVector.hh"
38 #include "G4ErrorPropagatorData.hh"
39 #include "G4ErrorSurfaceTarget.hh"
40 
43 
44 
45 //-------------------------------------------------------------------
46 
48  : G4Navigator()
49 {
50 }
51 
52 //-------------------------------------------------------------------
53 
55 {
56 }
57 
58 //-------------------------------------------------------------------
59 
61 ComputeStep ( const G4ThreeVector &pGlobalPoint,
62  const G4ThreeVector &pDirection,
63  const G4double pCurrentProposedStepLength,
64  G4double &pNewSafety )
65 {
66  G4double safetyGeom= DBL_MAX;
67 
68  G4double Step = G4Navigator::ComputeStep(pGlobalPoint, pDirection,
69  pCurrentProposedStepLength,
70  safetyGeom);
71 
72  G4ErrorPropagatorData * g4edata
74 
75  if (g4edata !=0)
76  {
77  const G4ErrorTarget* target = g4edata->GetTarget();
78  if( target != 0 )
79  {
80  G4double StepPlane= target->GetDistanceFromPoint(pGlobalPoint,pDirection);
81 
82  if( StepPlane < 0. ) // Negative means target is crossed, will not be found
83  {
84  StepPlane = DBL_MAX;
85  }
86 #ifdef G4VERBOSE
88  {
89  G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
90  << " Target step: " << StepPlane
91  << ", Transportation step: " << Step << G4endl;
92  target->Dump( "G4ErrorPropagationNavigator::ComputeStep Target " );
93  }
94 #endif
95 
96  if(StepPlane<Step)
97  {
98 #ifdef G4VERBOSE
100  {
101  G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
102  << " TargetCloserThanBoundary: " << StepPlane << " < "
103  << Step << G4endl;
104  }
105 #endif
106  Step = StepPlane;
108  }
109  else
110  {
112  }
113  }
114  }
115  G4double safetyTarget = TargetSafetyFromPoint(pGlobalPoint);
116  // Avoid call to G4Navigator::ComputeSafety - which could have side effects
117  pNewSafety= std::min(safetyGeom, safetyTarget);
118 
119 #ifdef G4VERBOSE
120  if( G4ErrorPropagatorData::verbose() >= 3 )
121  {
122  G4cout << "G4ErrorPropagationNavigator::ComputeStep()" << G4endl
123  << " Step: " << Step << ", ComputeSafety: " << pNewSafety
124  << G4endl;
125  }
126 #endif
127 
128  return Step;
129 }
130 
131 //-------------------------------------------------------------------
132 
134 TargetSafetyFromPoint( const G4ThreeVector &pGlobalpoint )
135 {
136  G4double safety= DBL_MAX;
137 
138  G4ErrorPropagatorData *g4edata
140 
141  if (g4edata !=0)
142  {
143  const G4ErrorTarget* target = g4edata->GetTarget();
144  if( target != 0 )
145  {
146  safety = target->GetDistanceFromPoint(pGlobalpoint);
147  }
148  }
149  return safety;
150 }
151 
152 //-------------------------------------------------------------------
153 
155 ComputeSafety( const G4ThreeVector &pGlobalPoint,
156  const G4double pMaxLength,
157  const G4bool keepState )
158 {
159  G4double safetyGeom = G4Navigator::ComputeSafety(pGlobalPoint,
160  pMaxLength, keepState);
161 
162  G4double safetyTarget = TargetSafetyFromPoint( pGlobalPoint );
163 
164  return std::min(safetyGeom, safetyTarget);
165 }
166 
167 //-------------------------------------------------------------------
168 
171 {
172  G4ErrorPropagatorData *g4edata
174  const G4ErrorTarget* target = 0;
175 
176  G4ThreeVector normal(0.0, 0.0, 0.0);
177  G4double distance= 0;
178 
179  // Determine which 'geometry' limited the step
180  if (g4edata)
181  {
182  target = g4edata->GetTarget();
183  if(target)
184  {
185  distance = target->GetDistanceFromPoint(point);
186  }
187  }
188 
189  if( distance > kCarTolerance // Not reached the target.
190  || (!target) )
191  // If a target does not exist, this seems the best we can do
192  {
193  normal= G4Navigator::GetGlobalExitNormal(point, valid);
194  }
195  else
196  {
197  switch( target->GetType() )
198  {
200  // The volume is in the 'real' mass geometry
201  normal= G4Navigator::GetGlobalExitNormal(point, valid);
202  break;
203  case G4ErrorTarget_TrkL:
204  normal= G4ThreeVector( 0.0, 0.0, 0.0);
205  *valid= false;
206  G4Exception("G4ErrorPropagationNavigator::GetGlobalExitNormal",
207  "Geometry1003",
208  JustWarning, "Unexpected value of Target type");
209  break;
212  const G4ErrorSurfaceTarget* surfaceTarget=
213  static_cast<const G4ErrorSurfaceTarget*>(target);
214  normal= surfaceTarget->GetTangentPlane(point).normal().unit();
215  *valid= true;
216  break;
217 
218 // default:
219 // normal= G4ThreeVector( 0.0, 0.0, 0.0);
220 // *valid= false;
221 // G4Exception("G4ErrorPropagationNavigator::GetGlobalExitNormal",
222 // "Geometry:003",
223 // FatalException, "Impossible value of Target type");
224 // exit(1);
225 // break;
226  }
227  }
228  return normal;
229 }
230 
const XML_Char * target
Definition: expat.h:268
CLHEP::Hep3Vector G4ThreeVector
virtual G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
Definition: G4Navigator.cc:747
Normal3D< T > normal() const
Definition: Plane3D.h:90
G4ErrorTargetType GetType() const
static double normal(HepRandomEngine *eptr)
Definition: RandPoisson.cc:77
virtual void Dump(const G4String &msg) const =0
void SetState(G4ErrorState sta)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4double TargetSafetyFromPoint(const G4ThreeVector &pGlobalpoint)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double kCarTolerance
Definition: G4Navigator.hh:361
G4double ComputeStep(const G4ThreeVector &pGlobalPoint, const G4ThreeVector &pDirection, const G4double pCurrentProposedStepLength, G4double &pNewSafety)
virtual G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4double GetDistanceFromPoint(const G4ThreeVector &, const G4ThreeVector &) const
const G4ErrorTarget * GetTarget(G4bool mustExist=0) const
#define G4endl
Definition: G4ios.hh:61
G4ThreeVector GetGlobalExitNormal(const G4ThreeVector &point, G4bool *valid)
double G4double
Definition: G4Types.hh:76
static G4ErrorPropagatorData * GetErrorPropagatorData()
virtual G4double ComputeSafety(const G4ThreeVector &globalpoint, const G4double pProposedMaxLength=DBL_MAX, const G4bool keepState=true)
#define DBL_MAX
Definition: templates.hh:83
virtual G4Plane3D GetTangentPlane(const G4ThreeVector &point) const =0