Geant4 11.2.2
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// Created: D.Sorokin, 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
45template <class T, G4bool StepperCachesDchord = true>
47{
48 public:
49
51 T* stepper,
52 G4int numberOfComponents = 6,
53 G4int statisticsVerbosity = 0);
54
56
59
61 G4double hstep,
62 G4double eps,
63 G4double chordDistance) override;
64
65 void OnStartTracking() override;
66 void OnComputeStep(const G4FieldTrack* /*track*/ = nullptr) override;
67 G4bool DoesReIntegrate() const override { return false; }
68 // Interpolation driver does not recalculate when AccurateAdvance
69 // is called -- reintegration would require other calls
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
79 void SetVerboseLevel(G4int level) override;
80 G4int GetVerboseLevel() const override;
81
82 void StreamInfo(std::ostream& os) const override;
83
84 protected:
85
93
94 using StepperIterator = typename std::vector<InterpStepper>::iterator;
95 using ConstStepperIterator = typename std::vector<InterpStepper>::const_iterator;
96
100 G4double& hstep,
101 G4double eps,
102 G4double curveLength,
103 G4FieldTrack* track = nullptr);
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 // return hdid
108
109 void Interpolate(G4double curveLength, field_utils::State& y) const;
110
111 void InterpolateImpl(G4double curveLength,
113 field_utils::State& y) const;
114
116 G4double curveLengthBegin,
117 const field_utils::State& yEnd,
118 G4double curveLengthEnd) const;
119
121 G4double curveLengthBegin,
122 field_utils::State& yEnd,
123 G4double curveLengthEnd,
124 G4double dChord,
125 G4double maxChordDistance);
126
128 G4double dChordStep,
129 G4double fDeltaChord);
130
131 void PrintState() const;
132
133 void CheckState() const;
134
136
137 protected:
138
139 std::vector<InterpStepper> fSteppers;
142
144 // Memory of last good step size for integration
145
147 // Minimum Step allowed in a Step (in units of length) // Parameter
148
150 const G4double fFractionNextEstimate = 0.98; // Constant
151 const G4double fSmallestCurveFraction = 0.01; // Constant
152
153 G4int fVerboseLevel; // Parameter
154
157
158 G4int fMaxTrials = 100; // Constant
160
161 // statistics
165
167};
168
169#include "G4InterpolationDriver.icc"
170
171#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
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)
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
G4double[G4FieldTrack::ncompSVEC] State
#define DBL_MAX
Definition templates.hh:62