Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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
39#include "G4FieldUtils.hh"
40
41#include "globals.hh"
42
43#include <vector>
44#include <memory>
45
46template <class T>
48{
49 public:
50
52 T* stepper,
53 G4int numberOfComponents = 6,
54 G4int statisticsVerbosity = 0);
55
56 virtual ~G4InterpolationDriver() override;
57
60
62 G4double hstep,
63 G4double eps,
64 G4double chordDistance) override;
65
66 virtual void OnStartTracking() override;
67 virtual void OnComputeStep() override;
68 virtual G4bool DoesReIntegrate() const override { return false; }
69 // Interpolation driver does not recalculate when AccurateAdvance is called
70 // -- reintegration would require other calls
71
73 G4double hstep,
74 G4double eps, // Requested y_err/hstep
75 G4double hinitial = 0) override;
76 // Integrates ODE from current s (s=s0) to s=s0+h with accuracy eps.
77 // On output track is replaced by value at end of interval.
78 // The concept is similar to the odeint routine from NRC p.721-722.
79
80 virtual void SetVerboseLevel(G4int level) override;
81 virtual G4int GetVerboseLevel() const override;
82
83 virtual void StreamInfo( std::ostream& os ) const override;
84
85 private:
86
87 struct InterpStepper
88 {
89 std::unique_ptr<T> stepper;
90 G4double begin;
91 G4double end;
92 G4double inverseLength;
93 };
94
95 using StepperIterator = typename std::vector<InterpStepper>::iterator;
96 using ConstStepperIterator = typename std::vector<InterpStepper>::const_iterator;
97
98 G4double OneGoodStep(StepperIterator it,
100 field_utils::State& dydx,
101 G4double& hstep,
102 G4double eps,
103 G4double curveLength);
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,
112 ConstStepperIterator it,
113 field_utils::State& y) const;
114
115 G4double DistChord(const field_utils::State& yBegin,
116 G4double curveLengthBegin,
117 const field_utils::State& yEnd,
118 G4double curveLengthEnd) const;
119
120 G4double FindNextChord(const field_utils::State& yBegin,
121 G4double curveLengthBegin,
122 field_utils::State& yEnd,
123 G4double curveLengthEnd,
124 G4double dChord,
125 G4double maxChordDistance);
126
127 G4double CalcChordStep(G4double stepTrialOld,
128 G4double dChordStep,
129 G4double fDeltaChord);
130
131 void PrintState() const;
132
133 void CheckState() const;
134
135 void AccumulateStatistics(G4int noTrials);
136
137 private:
138
139 std::vector<InterpStepper> fSteppers;
140 StepperIterator fLastStepper;
141 G4bool fKeepLastStepper = false;
142
143 G4double fhnext = DBL_MAX; // Memory of last good step size for integration
144
145 // Minimum Step allowed in a Step (in units of length) // Parameter
146 G4double fMinimumStep;
147
148 G4double fChordStepEstimate = DBL_MAX;
149 const G4double fFractionNextEstimate = 0.98; // Constant
150 const G4double fSmallestCurveFraction = 0.01; // Constant
151
152 G4int fVerboseLevel; // Parameter
153
154 field_utils::State fdydx;
155 G4bool fFirstStep = true;
156
157 const G4int fMaxTrials = 100; // Constant
158 G4int fTotalStepsForTrack = 0;
159
160 // statistics
161 G4int fTotalNoTrials = 0;
162 G4int fNoCalls = 0;
163 G4int fmaxTrials = 0;
164
165 using Base = G4RKIntegrationDriver<T>;
166};
167
168#include "G4InterpolationDriver.icc"
169
170#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4InterpolationDriver(const G4InterpolationDriver &)=delete
virtual ~G4InterpolationDriver() override
const G4InterpolationDriver & operator=(const G4InterpolationDriver &)=delete
virtual void OnComputeStep() override
virtual void StreamInfo(std::ostream &os) const override
virtual void OnStartTracking() override
virtual G4bool DoesReIntegrate() const override
virtual G4int GetVerboseLevel() const override
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0) override
virtual G4double AdvanceChordLimited(G4FieldTrack &track, G4double hstep, G4double eps, G4double chordDistance) override
G4InterpolationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
virtual void SetVerboseLevel(G4int level) override
G4double[G4FieldTrack::ncompSVEC] State
Definition: G4FieldUtils.hh:45
#define DBL_MAX
Definition: templates.hh:62