Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChordFinder.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//
27// $Id$
28//
29//
30// Class G4ChordFinder
31//
32// class description:
33//
34// A class that provides RK integration of motion ODE (as does g4magtr)
35// and also has a method that returns an Approximate point on the curve
36// near to a (chord) point.
37
38// History:
39// - 25.02.97 - John Apostolakis - Design and implementation
40// -------------------------------------------------------------------
41
42#ifndef G4CHORDFINDER_HH
43#define G4CHORDFINDER_HH
44
46#include "G4FieldTrack.hh"
47
48class G4MagneticField;
49
51{
52 public: // with description
53
54 G4ChordFinder( G4MagInt_Driver* pIntegrationDriver );
55
56 G4ChordFinder( G4MagneticField* itsMagField,
57 G4double stepMinimum = 1.0e-2, // * mm
58 G4MagIntegratorStepper* pItsStepper = 0 );
59 // A constructor that creates defaults for all "children" classes.
60
61 virtual ~G4ChordFinder();
62
64 G4double stepInitial,
65 G4double epsStep_Relative,
66 const G4ThreeVector latestSafetyOrigin,
67 G4double lasestSafetyRadius);
68 // Uses ODE solver's driver to find the endpoint that satisfies
69 // the chord criterion: that d_chord < delta_chord
70 // -> Returns Length of Step taken.
71
72 G4FieldTrack ApproxCurvePointS( const G4FieldTrack& curveAPointVelocity,
73 const G4FieldTrack& curveBPointVelocity,
74 const G4FieldTrack& ApproxCurveV,
75 const G4ThreeVector& currentEPoint,
76 const G4ThreeVector& currentFPoint,
77 const G4ThreeVector& PointG,
78 G4bool first, G4double epsStep);
79
80 G4FieldTrack ApproxCurvePointV( const G4FieldTrack& curveAPointVelocity,
81 const G4FieldTrack& curveBPointVelocity,
82 const G4ThreeVector& currentEPoint,
83 G4double epsStep);
84
85 inline G4double InvParabolic( const G4double xa, const G4double ya,
86 const G4double xb, const G4double yb,
87 const G4double xc, const G4double yc );
88
89 inline G4double GetDeltaChord() const;
90 inline void SetDeltaChord(G4double newval);
91
92 inline void SetChargeMomentumMass(G4double pCharge, // in e+ units
93 G4double pMomentum,
94 G4double pMass );
95 // Function to inform integration driver of charge, speed.
96
97 inline void SetIntegrationDriver(G4MagInt_Driver* IntegrationDriver);
99 // Access and set Driver.
100
101 inline void ResetStepEstimate();
102 // Clear internal state (last step estimate)
103
104 inline G4int GetNoCalls();
105 inline G4int GetNoTrials(); // Total number of trials
106 inline G4int GetNoMaxTrials(); // Maximum # of trials for one call
107 // Get statistics about number of calls & trials in FindNextChord
108
109 virtual void PrintStatistics();
110 // A report with the above -- and possibly other stats
111 inline G4int SetVerbose( G4int newvalue=1);
112 // Set verbosity and return old value
113
114 void SetFractions_Last_Next( G4double fractLast= 0.90,
115 G4double fractNext= 0.95 );
116 // Parameters for performance ... change with great care
117
118 inline void SetFirstFraction(G4double fractFirst);
119 // Parameter for performance ... change with great care
120
121 public: // without description
122
123 void TestChordPrint( G4int noTrials,
124 G4int lastStepTrial,
125 G4double dChordStep,
126 G4double nextStepTrial );
127
128 // Printing for monitoring ...
129
130 inline G4double GetFirstFraction(); // Originally 0.999
131 inline G4double GetFractionLast(); // Originally 1.000
132 inline G4double GetFractionNextEstimate(); // Originally 0.980
133 inline G4double GetMultipleRadius(); // No original value
134 // Parameters for adapting performance ... use with great care
135
136 protected: // .........................................................
137
138 inline void AccumulateStatistics( G4int noTrials );
139 // Accumulate the basic statistics
140 // - other specialised ones must be kept by derived classes
141
142 inline G4bool AcceptableMissDist(G4double dChordStep) const;
143
144 G4double NewStep( G4double stepTrialOld,
145 G4double dChordStep, // Current dchord estimate
146 G4double& stepEstimate_Unconstrained ) ;
147
148 virtual G4double FindNextChord( const G4FieldTrack& yStart,
149 G4double stepMax,
150 G4FieldTrack& yEnd,
151 G4double& dyErr, // Error of endpoint
152 G4double epsStep,
153 G4double* pNextStepForAccuracy, // = 0,
154 const G4ThreeVector latestSafetyOrigin,
155 G4double latestSafetyRadius
156 );
157
158 void PrintDchordTrial(G4int noTrials,
159 G4double stepTrial,
160 G4double oldStepTrial,
161 G4double dChordStep);
162
164 inline void SetLastStepEstimateUnc( G4double stepEst );
165
166 private: // ............................................................
167
169 G4ChordFinder& operator=(const G4ChordFinder&);
170 // Private copy constructor and assignment operator.
171
172 private: // ............................................................
173 // G4int nOK, nBAD;
174
175 // Constants
176 const G4double fDefaultDeltaChord; // SET in G4ChordFinder.cc = 0.25 mm
177
178 // PARAMETERS
179 // ---------------------
180 G4double fDeltaChord; // Maximum miss distance
181 // Internal parameters
182 G4double fFirstFraction, fFractionLast, fFractionNextEstimate;
183 G4double fMultipleRadius;
184 G4int fStatsVerbose; // if > 0, print Statistics in destructor
185
186 // DEPENDENT Objects
187 // ---------------------
188 G4MagInt_Driver* fIntgrDriver;
189 G4MagIntegratorStepper* fDriversStepper;
190 G4bool fAllocatedStepper; // Bookkeeping of dependent object
191 G4EquationOfMotion* fEquation;
192
193 // STATE information
194 // --------------------
195 G4double fLastStepEstimate_Unconstrained;
196 // State information for efficiency
197
198 // For Statistics
199 // -- G4int fNoTrials, fNoCalls;
200 G4int fTotalNoTrials_FNC, fNoCalls_FNC, fmaxTrials_FNC; // fnoTimesMaxTrFNC;
201};
202
203// Inline function implementation:
204
205#include "G4ChordFinder.icc"
206
207#endif // G4CHORDFINDER_HH
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void SetLastStepEstimateUnc(G4double stepEst)
void PrintDchordTrial(G4int noTrials, G4double stepTrial, G4double oldStepTrial, G4double dChordStep)
G4double InvParabolic(const G4double xa, const G4double ya, const G4double xb, const G4double yb, const G4double xc, const G4double yc)
G4int GetNoTrials()
G4bool AcceptableMissDist(G4double dChordStep) const
G4int GetNoMaxTrials()
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
G4int GetNoCalls()
G4double NewStep(G4double stepTrialOld, G4double dChordStep, G4double &stepEstimate_Unconstrained)
G4double GetLastStepEstimateUnc()
G4FieldTrack ApproxCurvePointS(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4FieldTrack &ApproxCurveV, const G4ThreeVector &currentEPoint, const G4ThreeVector &currentFPoint, const G4ThreeVector &PointG, G4bool first, G4double epsStep)
virtual ~G4ChordFinder()
G4double GetFractionNextEstimate()
G4double GetDeltaChord() const
G4int SetVerbose(G4int newvalue=1)
void AccumulateStatistics(G4int noTrials)
void SetChargeMomentumMass(G4double pCharge, G4double pMomentum, G4double pMass)
void SetIntegrationDriver(G4MagInt_Driver *IntegrationDriver)
void ResetStepEstimate()
void TestChordPrint(G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double nextStepTrial)
G4double GetFractionLast()
void SetDeltaChord(G4double newval)
virtual G4double FindNextChord(const G4FieldTrack &yStart, G4double stepMax, G4FieldTrack &yEnd, G4double &dyErr, G4double epsStep, G4double *pNextStepForAccuracy, const G4ThreeVector latestSafetyOrigin, G4double latestSafetyRadius)
G4double GetMultipleRadius()
void SetFirstFraction(G4double fractFirst)
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector latestSafetyOrigin, G4double lasestSafetyRadius)
virtual void PrintStatistics()
void SetFractions_Last_Next(G4double fractLast=0.90, G4double fractNext=0.95)
G4MagInt_Driver * GetIntegrationDriver()
G4double GetFirstFraction()