Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4IntegrationDriver.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// G4IntegrationDriver
27//
28// Class description:
29//
30// Templated driver class which controls the integration error of a
31// Runge-Kutta stepper.
32// It's purpose is to provide a replacement of G4MagIntegratorDriver
33// and work for all types of steppers.
34// It will serve as the driver of choice for steppers which do not
35// have extra capabilities, in particular First Same As Last (FSAL)
36// and/or interpolation.
37
38// Author: Dmitry Sorokin, Google Summer of Code 2017
39// Supervision: John Apostolakis, CERN
40// --------------------------------------------------------------------
41#ifndef G4INTEGRATIONDRIVER_HH
42#define G4INTEGRATIONDRIVER_HH
43
46
47template <class T>
49 public G4ChordFinderDelegate<G4IntegrationDriver<T>>
50{
51 public:
52
54 T* stepper,
55 G4int numberOfComponents = 6,
56 G4int statisticsVerbosity = 0 );
58
61
63 G4double stepMax,
64 G4double epsStep,
65 G4double chordDistance) override;
66
67 void OnStartTracking() override;
68 void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override {}
69 G4bool DoesReIntegrate() const override { return true; }
70
72 G4double hstep,
73 G4double eps, // Requested y_err/hstep
74 G4double hinitial = 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
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 void SetVerboseLevel(G4int newLevel) override;
87 G4int GetVerboseLevel() const override;
88
89 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 yVar[], // InOut
98 const 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{1e-12};
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{0};
135 G4int fNoAccurateAdvanceCalls{0};
136 G4int fNoAccurateAdvanceBadSteps{0};
137 G4int fNoAccurateAdvanceGoodSteps{0};
138
139 using Base = G4RKIntegrationDriver<T>;
140 using ChordFinderDelegate = G4ChordFinderDelegate<G4IntegrationDriver<T>>;
141};
142
143#include "G4IntegrationDriver.icc"
144
145#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void OneGoodStep(G4double yVar[], const G4double dydx[], G4double &curveLength, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
G4IntegrationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
G4bool QuickAdvance(G4FieldTrack &fieldTrack, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
G4bool DoesReIntegrate() const override
G4double GetMinimumStep() const
void OnStartTracking() override
void StreamInfo(std::ostream &os) const override
G4int GetVerboseLevel() const override
void IncrementQuickAdvanceCalls()
const G4IntegrationDriver & operator=(const G4IntegrationDriver &)=delete
~G4IntegrationDriver() override
G4double GetSmallestFraction() const
G4IntegrationDriver(const G4IntegrationDriver &)=delete
void SetSmallestFraction(G4double val)
G4double AdvanceChordLimited(G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override
void SetMinimumStep(G4double newval)
void OnComputeStep(const G4FieldTrack *=nullptr) override
void SetVerboseLevel(G4int newLevel) override
G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0) override