Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4InterpolationDriver.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// G4InterpolationDriver
27//
28// Class description:
29//
30// Driver class which uses Runge-Kutta stepper with interpolation property
31// to integrate track with error control
32
33// Author: Dmitry Sorokin (CERN), 26.09.2018
34// --------------------------------------------------------------------
35#ifndef G4INTERPOLATION_DRIVER_HH
36#define G4INTERPOLATION_DRIVER_HH
37
38#include "G4FieldUtils.hh"
40#include "globals.hh"
41
42#include <memory>
43#include <vector>
44
45/**
46 * @brief G4InterpolationDriver is a templated driver class which uses
47 * Runge-Kutta stepper with interpolation property to integrate track with
48 * error control.
49 */
50
51template <class T, G4bool StepperCachesDchord = true>
53{
54 public:
55
56 /**
57 * Constructor for G4IntegrationDriver.
58 * @param[in] hminimum Minimum allowed step.
59 * @param[in] stepper Pointer to the stepper algorithm.
60 * @param[in] numberOfComponents The number of integration variables,
61 * if not matching stepper's number of variables, issue exception.
62 * @param[in] statisticsVerbosity Verbosity level.
63 */
65 T* stepper,
66 G4int numberOfComponents = 6,
67 G4int statisticsVerbosity = 0);
68
69 /**
70 * Destructor. Provides statistics if verbosity level is greater than zero.
71 */
73
74 /**
75 * Copy constructor and assignment operator not allowed.
76 */
79
80 /**
81 * Computes the step to take, based on chord limits.
82 * @param[in,out] track The current track in field.
83 * @param[in] hstep Proposed step length.
84 * @param[in] eps Requested accuracy, y_err/hstep.
85 * @param[in] chordDistance Maximum sagitta distance.
86 * @returns The length of step taken.
87 */
89 G4double hstep,
90 G4double eps,
91 G4double chordDistance) override;
92
93 /**
94 * Dispatch interface method for initialisation/reset of driver.
95 */
96 void OnStartTracking() override;
97
98 /**
99 * Dispatch interface method for computing step. Does nothing here.
100 */
101 void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override;
102
103 /**
104 * The driver does not implement re-integration. Returns false.
105 */
106 G4bool DoesReIntegrate() const override { return false; }
107
108 /**
109 * Advances integration accurately by relative accuracy better than 'eps'.
110 * On output the track is replaced by the value at the end of interval.
111 * @param[in,out] track The current track in field.
112 * @param[in] hstep Proposed step length.
113 * @param[in] eps Requested accuracy, y_err/hstep.
114 * @param[in] hinitial Initial minimum integration step.
115 * @returns true if integration succeeds.
116 */
118 G4double hstep,
119 G4double eps, // Requested y_err/hstep
120 G4double hinitial = 0) override;
121
122 /**
123 * Setter and getter for verbosity.
124 */
125 void SetVerboseLevel(G4int level) override;
126 G4int GetVerboseLevel() const override;
127
128 /**
129 * Writes out to stream the parameters/state of the driver.
130 */
131 void StreamInfo(std::ostream& os) const override;
132
133 protected:
134
142
143 using StepperIterator = typename std::vector<InterpStepper>::iterator;
144 using ConstStepperIterator = typename std::vector<InterpStepper>::const_iterator;
145
146 /**
147 * Takes one Step that is as large as possible while satisfying the
148 * accuracy criterion.
149 * @param[in] it Stepper iterator.
150 * @param[in,out] y The current track state, y.
151 * @param[in] dydx dydx array.
152 * @param[in,out] hstep Step to attempt.
153 * @param[in] eps The relative accuracy.
154 * @param[in] curveLength Step start, x.
155 * @param[in,out] track Pointer to the Field track. Not used.
156 * @returns The step achieved.
157 */
160 field_utils::State& dydx,
161 G4double& hstep,
162 G4double eps,
163 G4double curveLength,
164 G4FieldTrack* track = nullptr);
165
166 /**
167 * Track interpolation.
168 * @param[in] curveLength Step start, x.
169 * @param[in,out] y The current track state, y.
170 */
171 void Interpolate(G4double curveLength, field_utils::State& y) const;
172
173 /**
174 * Wrapper method for interpolation.
175 */
176 void InterpolateImpl(G4double curveLength,
178 field_utils::State& y) const;
179
180 /**
181 * Methods for calculation of chord step and distance.
182 */
184 G4double curveLengthBegin,
185 const field_utils::State& yEnd,
186 G4double curveLengthEnd) const;
188 G4double curveLengthBegin,
189 field_utils::State& yEnd,
190 G4double curveLengthEnd,
191 G4double dChord,
192 G4double maxChordDistance);
194 G4double dChordStep,
195 G4double fDeltaChord);
196
197
198 /**
199 * Internal methods for printing/checking the state.
200 */
201 void PrintState() const;
202 void CheckState() const;
203
204 /**
205 * Increments number of trials and calls.
206 */
208
209 protected:
210
211 std::vector<InterpStepper> fSteppers;
214
215 /** Memory of last good step size for integration. */
217
218 /** Minimum Step allowed (in units of length). */
220
222 const G4double fFractionNextEstimate = 0.98; // Constant
223 const G4double fSmallestCurveFraction = 0.01; // Constant
224
225 G4int fVerboseLevel; // Parameter
226
229
230 G4int fMaxTrials = 100; // Constant
232
233 /** Statistics. */
237
239};
240
241#include "G4InterpolationDriver.icc"
242
243#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4FieldTrack defines a data structure bringing together a magnetic track's state (position,...
G4bool DoesReIntegrate() const override
G4double CalcChordStep(G4double stepTrialOld, G4double dChordStep, G4double fDeltaChord)
G4int GetVerboseLevel() const override
~G4InterpolationDriver() override
G4InterpolationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
virtual G4double OneGoodStep(StepperIterator it, field_utils::State &y, field_utils::State &dydx, G4double &hstep, G4double eps, G4double curveLength, G4FieldTrack *track=nullptr)
void AccumulateStatistics(G4int noTrials)
G4double AdvanceChordLimited(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
void CheckState() const
void PrintState() const
G4InterpolationDriver(const G4InterpolationDriver &)=delete
void StreamInfo(std::ostream &os) const override
G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0) override
const G4InterpolationDriver & operator=(const G4InterpolationDriver &)=delete
void InterpolateImpl(G4double curveLength, ConstStepperIterator it, field_utils::State &y) const
typename std::vector< InterpStepper >::iterator StepperIterator
void OnStartTracking() override
typename std::vector< InterpStepper >::const_iterator ConstStepperIterator
void Interpolate(G4double curveLength, field_utils::State &y) const
G4double FindNextChord(const field_utils::State &yBegin, G4double curveLengthBegin, field_utils::State &yEnd, G4double curveLengthEnd, G4double dChord, G4double maxChordDistance)
G4RKIntegrationDriver< T > Base
std::vector< InterpStepper > fSteppers
void SetVerboseLevel(G4int level) override
G4double DistChord(const field_utils::State &yBegin, G4double curveLengthBegin, const field_utils::State &yEnd, G4double curveLengthEnd) const
void OnComputeStep(const G4FieldTrack *=nullptr) override
G4RKIntegrationDriver(T *stepper)
G4double[G4FieldTrack::ncompSVEC] State
#define DBL_MAX
Definition templates.hh:62