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