Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 (CERN, Google Summer of Code 2017), 20.10.2017
39// Supervision: John Apostolakis (CERN)
40// --------------------------------------------------------------------
41#ifndef G4INTEGRATIONDRIVER_HH
42#define G4INTEGRATIONDRIVER_HH
43
46
47/**
48 * @brief G4IntegrationDriver is a templated driver class which controls the
49 * integration error of a Runge-Kutta stepper.
50 * It serves as the driver of choice for steppers which do not have extra
51 * capabilities, in particular First Same As Last (FSAL) and/or interpolation.
52 */
53
54template <class T>
56 public G4ChordFinderDelegate<G4IntegrationDriver<T>>
57{
58 public:
59
60 /**
61 * Constructor for G4IntegrationDriver.
62 * @param[in] hminimum Minimum allowed step.
63 * @param[in] stepper Pointer to the stepper algorithm.
64 * @param[in] numberOfComponents The number of integration variables,
65 * if not matching stepper's number of variables, issue exception.
66 * @param[in] statisticsVerbosity Verbosity level.
67 */
68 inline G4IntegrationDriver( G4double hminimum,
69 T* stepper,
70 G4int numberOfComponents = 6,
71 G4int statisticsVerbosity = 0 );
72
73 /**
74 * Destructor. Provides statistics if verbosity level is greater than zero.
75 */
76 inline ~G4IntegrationDriver() override;
77
78 /**
79 * Copy constructor and assignment operator not allowed.
80 */
83
84 /**
85 * Computes the step to take, based on chord limits.
86 * @param[in,out] track The current track in field.
87 * @param[in] stepMax Proposed step length.
88 * @param[in] epsStep Requested accuracy, y_err/hstep.
89 * @param[in] chordDistance Maximum sagitta distance.
90 * @returns The length of step taken.
91 */
93 G4double stepMax,
94 G4double epsStep,
95 G4double chordDistance) override;
96
97 /**
98 * Dispatch interface method for initialisation/reset of driver.
99 */
100 inline void OnStartTracking() override;
101
102 /**
103 * Dispatch interface method for computing step. Does nothing here.
104 */
105 inline void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override;
106
107 /**
108 * The driver does implement re-integration. Returns true.
109 */
110 inline G4bool DoesReIntegrate() const override;
111
112 /**
113 * Advances integration accurately by relative accuracy better than 'eps'.
114 * On output the track is replaced by the value at the end of interval.
115 * @param[in,out] track The current track in field.
116 * @param[in] hstep Proposed step length.
117 * @param[in] eps Requested accuracy, y_err/hstep.
118 * @param[in] hinitial Initial minimum integration step.
119 * @returns true if integration succeeds.
120 */
122 G4double hstep,
123 G4double eps, // Requested y_err/hstep
124 G4double hinitial = 0 ) override;
125
126 /**
127 * Attempts one integration step, and returns estimated error 'dyerr'.
128 * It does not ensure accuracy.
129 * @param[in,out] fieldTrack The current track in field.
130 * @param[in] dydx dydx array.
131 * @param[in] hstep Proposed step length.
132 * @param[out] dchord_step Estimated sagitta distance.
133 * @param[out] dyerr Estimated error.
134 * @returns true if integration succeeds.
135 */
136 inline G4bool QuickAdvance(G4FieldTrack& fieldTrack,
137 const G4double dydx[],
138 G4double hstep,
139 G4double& dchord_step,
140 G4double& dyerr) override;
141
142 /**
143 * Takes one Step that is as large as possible while satisfying the
144 * accuracy criterion.
145 * @param[in,out] yVar The current track state, y.
146 * @param[in] dydx dydx array.
147 * @param[in,out] curveLength Step start, x.
148 * @param[in] htry Step to attempt.
149 * @param[in] eps The relative accuracy.
150 * @param[out] hdid Step achieved.
151 * @param[out] hnext Proposed next step.
152 */
153 inline void OneGoodStep(G4double yVar[], // InOut
154 const G4double dydx[],
155 G4double& curveLength,
156 G4double htry,
157 G4double eps,
158 G4double& hdid,
159 G4double& hnext);
160
161 /**
162 * Setter and getter for verbosity.
163 */
164 inline void SetVerboseLevel(G4int newLevel) override;
165 inline G4int GetVerboseLevel() const override;
166
167 /**
168 * Writes out to stream the parameters/state of the driver.
169 */
170 inline void StreamInfo( std::ostream& os ) const override;
171
172 /**
173 * Getter and Setter for minimum allowed step.
174 */
175 inline G4double GetMinimumStep() const;
176 inline void SetMinimumStep(G4double newval);
177
178 /**
179 * Getter and Setter for smallest fraction.
180 */
183
184 protected:
185
186 /**
187 * Increments the counter for the number of calls to QuickAdvance().
188 */
190
191 private:
192
193 /**
194 * Checks accuracy of step distance on the end point.
195 */
196 inline void CheckStep(const G4ThreeVector& posIn,
197 const G4ThreeVector& posOut, G4double hdid);
198
199 private:
200
201 /** Minimum Step allowed in a Step (in absolute units). */
202 G4double fMinimumStep;
203
204 /** Smallest fraction of (existing) curve length in relative units.
205 * Below this fraction the current step will be the last.
206 * The expected range: smaller than 0.1 * epsilon and bigger than 5e-13
207 * (range not enforced). */
208 G4double fSmallestFraction{1e-12};
209
210 /** Verbosity level for printing (debug, etc..)
211 * Could be varied during tracking to help identifying issues. */
212 G4int fVerboseLevel;
213
214 G4int fNoQuickAvanceCalls{0};
215 G4int fNoAccurateAdvanceCalls{0};
216 G4int fNoAccurateAdvanceBadSteps{0};
217 G4int fNoAccurateAdvanceGoodSteps{0};
218
219 using Base = G4RKIntegrationDriver<T>;
220 using ChordFinderDelegate = G4ChordFinderDelegate<G4IntegrationDriver<T>>;
221};
222
223#include "G4IntegrationDriver.icc"
224
225#endif
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4ChordFinderDelegate is a templated class for a common algorithm of finding step size with distance ...
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
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
G4RKIntegrationDriver(T *stepper)