Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 G4ChordFinderDelegate<G4FSALIntegrationDriver<T>>
43{
44 public:
45
47 T* stepper,
48 G4int numberOfComponents = 6,
49 G4int statisticsVerbosity = 1);
50
52
55
57 G4double hstep,
58 G4double eps,
59 G4double chordDistance) override;
60
65
66 void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override {}
67
68 G4bool DoesReIntegrate() const override { return true; }
69
71 G4double hstep,
72 G4double eps, // Requested y_err/hstep
73 G4double hinitial = 0.0) override;
74 // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
75 // On output track is replaced by value at end of interval.
76 // The concept is similar to the odeint routine from NRC p.721-722.
77
79 const G4double dydx[],
80 G4double hstep,
81 G4double& dchord_step,
82 G4double& dyerr ) override;
83 // QuickAdvance just tries one Step - it does not ensure accuracy.
84
85 void SetVerboseLevel(G4int newLevel) override;
86 G4int GetVerboseLevel() const override;
87
88 void StreamInfo( std::ostream& os ) const override;
89 // Write out the parameters / state of the driver
90
91 // Accessors
92
95
96 void OneGoodStep(G4double y[], // InOut
97 G4double dydx[],
98 G4double& curveLength,
99 G4double htry,
100 G4double eps,
101 G4double& hdid,
102 G4double& hnext);
103 // This takes one Step that is of size htry, or as large
104 // as possible while satisfying the accuracy criterion of:
105 // yerr < eps * |y_end-y_start|
106
109
110 protected:
111
113
114 private:
115
116 void CheckStep(const G4ThreeVector& posIn,
117 const G4ThreeVector& posOut,
118 G4double hdid);
119
120 G4double fMinimumStep;
121 // Minimum Step allowed in a Step (in absolute units)
122
123 G4double fSmallestFraction{1e-12};
124 // Smallest fraction of (existing) curve length - in relative units
125 // below this fraction the current step will be the last
126 // Expected range: smaller than 0.1 * epsilon and bigger than 5e-13
127 // ( Note: this range is not enforced. )
128
129 G4int fVerboseLevel;
130 // Verbosity level for printing (debug, ..)
131 // Could be varied during tracking - to help identify issues
132
133 G4int fNoQuickAvanceCalls{0};
134 G4int fNoAccurateAdvanceCalls{0};
135 G4int fNoAccurateAdvanceBadSteps{0};
136 G4int fNoAccurateAdvanceGoodSteps{0};
137
138 using Base = G4RKIntegrationDriver<T>;
139 using ChordFinderDelegate = G4ChordFinderDelegate<G4FSALIntegrationDriver<T>>;
140};
141
142#include "G4FSALIntegrationDriver.icc"
143
144#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4int GetVerboseLevel() const override
void StreamInfo(std::ostream &os) const override
void OnComputeStep(const G4FieldTrack *=nullptr) override
G4FSALIntegrationDriver(const G4FSALIntegrationDriver &)=delete
G4double GetSmallestFraction() const
~G4FSALIntegrationDriver() override
G4double AdvanceChordLimited(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
G4bool DoesReIntegrate() const override
void SetMinimumStep(G4double newval)
G4double GetMinimumStep() const
G4FSALIntegrationDriver & operator=(const G4FSALIntegrationDriver &)=delete
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)
void SetSmallestFraction(G4double val)
G4bool QuickAdvance(G4FieldTrack &fieldTrack, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
void SetVerboseLevel(G4int newLevel) override
G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0.0) override