Geant4  10.02
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: G4ChordFinder.hh 69699 2013-05-13 08:50:30Z gcosmo $
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 
45 #include "G4MagIntegratorDriver.hh"
46 #include "G4FieldTrack.hh"
47 
48 class 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 SetIntegrationDriver(G4MagInt_Driver* IntegrationDriver);
94  // Access and set Driver.
95 
96  inline void ResetStepEstimate();
97  // Clear internal state (last step estimate)
98 
99  inline G4int GetNoCalls();
100  inline G4int GetNoTrials(); // Total number of trials
101  inline G4int GetNoMaxTrials(); // Maximum # of trials for one call
102  // Get statistics about number of calls & trials in FindNextChord
103 
104  virtual void PrintStatistics();
105  // A report with the above -- and possibly other stats
106  inline G4int SetVerbose( G4int newvalue=1);
107  // Set verbosity and return old value
108 
109  void SetFractions_Last_Next( G4double fractLast= 0.90,
110  G4double fractNext= 0.95 );
111  // Parameters for performance ... change with great care
112 
113  inline void SetFirstFraction(G4double fractFirst);
114  // Parameter for performance ... change with great care
115 
116  public: // without description
117 
118  void TestChordPrint( G4int noTrials,
119  G4int lastStepTrial,
120  G4double dChordStep,
121  G4double nextStepTrial );
122 
123  // Printing for monitoring ...
124 
125  inline G4double GetFirstFraction(); // Originally 0.999
126  inline G4double GetFractionLast(); // Originally 1.000
127  inline G4double GetFractionNextEstimate(); // Originally 0.980
128  inline G4double GetMultipleRadius(); // No original value
129  // Parameters for adapting performance ... use with great care
130 
131  protected: // .........................................................
132 
133  inline void AccumulateStatistics( G4int noTrials );
134  // Accumulate the basic statistics
135  // - other specialised ones must be kept by derived classes
136 
137  inline G4bool AcceptableMissDist(G4double dChordStep) const;
138 
139  G4double NewStep( G4double stepTrialOld,
140  G4double dChordStep, // Current dchord estimate
141  G4double& stepEstimate_Unconstrained ) ;
142 
143  virtual G4double FindNextChord( const G4FieldTrack& yStart,
144  G4double stepMax,
145  G4FieldTrack& yEnd,
146  G4double& dyErr, // Error of endpoint
147  G4double epsStep,
148  G4double* pNextStepForAccuracy, // = 0,
149  const G4ThreeVector latestSafetyOrigin,
150  G4double latestSafetyRadius
151  );
152 
153  void PrintDchordTrial(G4int noTrials,
154  G4double stepTrial,
155  G4double oldStepTrial,
156  G4double dChordStep);
157 
159  inline void SetLastStepEstimateUnc( G4double stepEst );
160 
161  private: // ............................................................
162 
165  // Private copy constructor and assignment operator.
166 
167  private: // ............................................................
168  // G4int nOK, nBAD;
169 
170  // Constants
171  const G4double fDefaultDeltaChord; // SET in G4ChordFinder.cc = 0.25 mm
172 
173  // PARAMETERS
174  // ---------------------
175  G4double fDeltaChord; // Maximum miss distance
176  // Internal parameters
179  G4int fStatsVerbose; // if > 0, print Statistics in destructor
180 
181  // DEPENDENT Objects
182  // ---------------------
185  G4bool fAllocatedStepper; // Bookkeeping of dependent object
187 
188  // STATE information
189  // --------------------
191  // State information for efficiency
192 
193  // For Statistics
194  // -- G4int fNoTrials, fNoCalls;
196 };
197 
198 // Inline function implementation:
199 
200 #include "G4ChordFinder.icc"
201 
202 #endif // G4CHORDFINDER_HH
G4double InvParabolic(const G4double xa, const G4double ya, const G4double xb, const G4double yb, const G4double xc, const G4double yc)
G4MagIntegratorStepper * fDriversStepper
G4double GetFirstFraction()
CLHEP::Hep3Vector G4ThreeVector
G4double fMultipleRadius
G4ChordFinder(G4MagInt_Driver *pIntegrationDriver)
G4double GetFractionNextEstimate()
void PrintDchordTrial(G4int noTrials, G4double stepTrial, G4double oldStepTrial, G4double dChordStep)
G4int fTotalNoTrials_FNC
G4FieldTrack ApproxCurvePointV(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4ThreeVector &currentEPoint, G4double epsStep)
void SetLastStepEstimateUnc(G4double stepEst)
G4double fDeltaChord
void SetIntegrationDriver(G4MagInt_Driver *IntegrationDriver)
const G4double fDefaultDeltaChord
int G4int
Definition: G4Types.hh:78
G4double fFractionLast
virtual void PrintStatistics()
G4int SetVerbose(G4int newvalue=1)
virtual G4double FindNextChord(const G4FieldTrack &yStart, G4double stepMax, G4FieldTrack &yEnd, G4double &dyErr, G4double epsStep, G4double *pNextStepForAccuracy, const G4ThreeVector latestSafetyOrigin, G4double latestSafetyRadius)
G4double fLastStepEstimate_Unconstrained
void TestChordPrint(G4int noTrials, G4int lastStepTrial, G4double dChordStep, G4double nextStepTrial)
G4bool fAllocatedStepper
bool G4bool
Definition: G4Types.hh:79
G4int GetNoTrials()
G4bool AcceptableMissDist(G4double dChordStep) const
G4double fFirstFraction
G4double fFractionNextEstimate
void AccumulateStatistics(G4int noTrials)
G4double AdvanceChordLimited(G4FieldTrack &yCurrent, G4double stepInitial, G4double epsStep_Relative, const G4ThreeVector latestSafetyOrigin, G4double lasestSafetyRadius)
G4double GetMultipleRadius()
G4int GetNoMaxTrials()
void SetFractions_Last_Next(G4double fractLast=0.90, G4double fractNext=0.95)
G4double NewStep(G4double stepTrialOld, G4double dChordStep, G4double &stepEstimate_Unconstrained)
G4EquationOfMotion * fEquation
G4MagInt_Driver * fIntgrDriver
G4double GetLastStepEstimateUnc()
G4int GetNoCalls()
G4double GetDeltaChord() const
virtual ~G4ChordFinder()
G4ChordFinder & operator=(const G4ChordFinder &)
double G4double
Definition: G4Types.hh:76
G4MagInt_Driver * GetIntegrationDriver()
G4double GetFractionLast()
void SetDeltaChord(G4double newval)
G4FieldTrack ApproxCurvePointS(const G4FieldTrack &curveAPointVelocity, const G4FieldTrack &curveBPointVelocity, const G4FieldTrack &ApproxCurveV, const G4ThreeVector &currentEPoint, const G4ThreeVector &currentFPoint, const G4ThreeVector &PointG, G4bool first, G4double epsStep)
void SetFirstFraction(G4double fractFirst)
void ResetStepEstimate()