Geant4 11.1.1
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, Google Summer of Code 2017
39// Supervision: John Apostolakis, CERN
40// --------------------------------------------------------------------
41#ifndef G4INTEGRATIONDRIVER_HH
42#define G4INTEGRATIONDRIVER_HH
43
46
47template <class T>
49 public G4ChordFinderDelegate<G4IntegrationDriver<T>>
50{
51 public:
52
54 T* stepper,
55 G4int numberOfComponents = 6,
56 G4int statisticsVerbosity = 0 );
57
58 virtual ~G4IntegrationDriver() override;
59
62
64 G4double stepMax,
65 G4double epsStep,
66 G4double chordDistance) override;
67
68 virtual void OnStartTracking() override;
69 virtual void OnComputeStep() override {}
70 virtual G4bool DoesReIntegrate() const override { return true; }
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 G4bool QuickAdvance( G4FieldTrack& fieldTrack,
81 const G4double dydx[],
82 G4double hstep,
83 G4double& dchord_step,
84 G4double& dyerr) override;
85 // QuickAdvance just tries one Step - it does not ensure accuracy.
86
87 virtual void SetVerboseLevel(G4int newLevel) override;
88 virtual G4int GetVerboseLevel() const override;
89
90 virtual void StreamInfo( std::ostream& os ) const override;
91 // Write out the parameters / state of the driver
92
93 // Accessors
94 //
97
98 void OneGoodStep( G4double yVar[], // InOut
99 const G4double dydx[],
100 G4double& curveLength,
101 G4double htry,
102 G4double eps,
103 G4double& hdid,
104 G4double& hnext);
105 // This takes one Step that is of size htry, or as large
106 // as possible while satisfying the accuracy criterion of:
107 // yerr < eps * |y_end-y_start|
108
111
112 protected:
113
115
116 private:
117
118 void CheckStep(const G4ThreeVector& posIn,
119 const G4ThreeVector& posOut,
120 G4double hdid);
121
122 G4double fMinimumStep;
123 // Minimum Step allowed in a Step (in absolute units)
124
125 G4double fSmallestFraction;
126 // Smallest fraction of (existing) curve length - in relative units
127 // below this fraction the current step will be the last
128 // Expected range: smaller than 0.1 * epsilon and bigger than 5e-13
129 // Note: this range is not enforced.
130
131 G4int fVerboseLevel;
132 // Verbosity level for printing (debug, ..)
133 // Could be varied during tracking - to help identify issues
134
135 G4int fNoQuickAvanceCalls;
136 G4int fNoAccurateAdvanceCalls;
137 G4int fNoAccurateAdvanceBadSteps;
138 G4int fNoAccurateAdvanceGoodSteps;
139
142};
143
144#include "G4IntegrationDriver.icc"
145
146#endif
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void OneGoodStep(G4double yVar[], const G4double dydx[], G4double &curveLength, G4double htry, G4double eps, G4double &hdid, G4double &hnext)
virtual ~G4IntegrationDriver() override
G4IntegrationDriver(G4double hminimum, T *stepper, G4int numberOfComponents=6, G4int statisticsVerbosity=0)
virtual void StreamInfo(std::ostream &os) const override
G4double GetMinimumStep() const
void IncrementQuickAdvanceCalls()
virtual G4bool DoesReIntegrate() const override
const G4IntegrationDriver & operator=(const G4IntegrationDriver &)=delete
virtual void OnStartTracking() override
virtual void SetVerboseLevel(G4int newLevel) override
G4double GetSmallestFraction() const
G4IntegrationDriver(const G4IntegrationDriver &)=delete
void SetSmallestFraction(G4double val)
virtual G4bool QuickAdvance(G4FieldTrack &fieldTrack, const G4double dydx[], G4double hstep, G4double &dchord_step, G4double &dyerr) override
virtual G4int GetVerboseLevel() const override
void SetMinimumStep(G4double newval)
virtual G4bool AccurateAdvance(G4FieldTrack &track, G4double hstep, G4double eps, G4double hinitial=0) override
virtual void OnComputeStep() override
virtual G4double AdvanceChordLimited(G4FieldTrack &track, G4double stepMax, G4double epsStep, G4double chordDistance) override