Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FSALIntegrationDriver.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// G4FSALIntegrationDriver
27//
28// Class description:
29//
30// Driver class which controls the integration error of a Runge-Kutta stepper
31
32// Created: D.Sorokin, 2017
33// --------------------------------------------------------------------
34#ifndef G4FSALINTEGRATIONDRIVER_HH
35#define G4FSALINTEGRATIONDRIVER_HH
36
39
40template <class T>
42 : public G4RKIntegrationDriver<T> ,
43 public G4ChordFinderDelegate<G4FSALIntegrationDriver<T>>
44{
45 public:
46
48 T* stepper,
49 G4int numberOfComponents = 6,
50 G4int statisticsVerbosity = 1);
51
52 virtual ~G4FSALIntegrationDriver() override;
53
56
58 G4double hstep,
59 G4double eps,
60 G4double chordDistance) override;
61
62 virtual void OnStartTracking() override
63 {
65 }
66
67 virtual void OnComputeStep() override {}
68
69 virtual G4bool DoesReIntegrate() const override { return true; }
70
72 G4double hstep,
73 G4double eps, // Requested y_err/hstep
74 G4double hinitial = 0.0) override;
75 // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
76 // On output track is replaced by value at end of interval.
77 // The concept is similar to the odeint routine from NRC p.721-722.
78
79 virtual G4bool QuickAdvance( G4FieldTrack& fieldTrack,
80 const G4double dydx[],
81 G4double hstep,
82 G4double& dchord_step,
83 G4double& dyerr) override;
84 // QuickAdvance just tries one Step - it does not ensure accuracy.
85
86 virtual void SetVerboseLevel(G4int newLevel) override;
87 virtual G4int GetVerboseLevel() const override;
88
89 virtual void StreamInfo( std::ostream& os ) const override;
90 // Write out the parameters / state of the driver
91
92 // Accessors
93
96
97 void OneGoodStep(G4double y[], // InOut
98 G4double dydx[],
99 G4double& curveLength,
100 G4double htry,
101 G4double eps,
102 G4double& hdid,
103 G4double& hnext);
104 // This takes one Step that is of size htry, or as large
105 // as possible while satisfying the accuracy criterion of:
106 // yerr < eps * |y_end-y_start|
107
110
111 protected:
112
114
115 private:
116
117 void CheckStep(const G4ThreeVector& posIn,
118 const G4ThreeVector& posOut,
119 G4double hdid);
120
121 G4double fMinimumStep;
122 // Minimum Step allowed in a Step (in absolute units)
123
124 G4double fSmallestFraction;
125 // Smallest fraction of (existing) curve length - in relative units
126 // below this fraction the current step will be the last
127 // Expected range: smaller than 0.1 * epsilon and bigger than 5e-13
128 // ( Note: this range is not enforced. )
129
130 G4int fVerboseLevel;
131 // Verbosity level for printing (debug, ..)
132 // Could be varied during tracking - to help identify issues
133
134 G4int fNoQuickAvanceCalls;
135 G4int fNoAccurateAdvanceCalls;
136 G4int fNoAccurateAdvanceBadSteps;
137 G4int fNoAccurateAdvanceGoodSteps;
138
141};
142
143#include "G4FSALIntegrationDriver.icc"
144
145#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void StreamInfo(std::ostream &os) const override
G4FSALIntegrationDriver(const G4FSALIntegrationDriver &)=delete
G4double GetSmallestFraction() const
virtual G4bool DoesReIntegrate() const override
virtual G4bool QuickAdvance(G4FieldTrack &fieldTrack, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
virtual void SetVerboseLevel(G4int newLevel) override
void SetMinimumStep(G4double newval)
G4double GetMinimumStep() const
virtual void OnStartTracking() override
virtual G4int GetVerboseLevel() const override
G4FSALIntegrationDriver & operator=(const G4FSALIntegrationDriver &)=delete
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0.0) override
virtual void OnComputeStep() override
G4FSALIntegrationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=1)
void OneGoodStep(G4double y[], G4double dydx[], G4double &curveLength, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
virtual G4double AdvanceChordLimited(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
void SetSmallestFraction(G4double val)
virtual ~G4FSALIntegrationDriver() override