Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RKIntegrationDriver.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// G4RKIntegrationDriver
27//
28// Class description:
29//
30// Driver class which controls the integration error of a
31// Runge-Kutta stepper
32
33// Created: D.Sorokin, 2017
34// --------------------------------------------------------------------
35#ifndef G4RKINTEGRATIONDRIVER_HH
36#define G4RKINTEGRATIONDRIVER_HH
37
39
40template <class T>
42{
43 public:
44
46
49
50 virtual void GetDerivatives(const G4FieldTrack& track,
51 G4double dydx[]) const override;
52
53 virtual void GetDerivatives(const G4FieldTrack& track,
54 G4double dydx[],
55 G4double field[]) const override;
56
57 virtual G4double ComputeNewStepSize(G4double errMaxNorm, // normalised error
58 G4double hstepCurrent) override final;
59 // Taking the last step's normalised error, calculate
60 // a step size for the next step.
61 // - Limits the next step's size within a factor of the current one.
62
64 virtual void SetEquationOfMotion(G4EquationOfMotion* equation) override;
65
66 virtual const T* GetStepper() const override;
67 virtual T* GetStepper() override;
68
69 virtual void StreamInfo( std::ostream& os ) const override;
70
71 // Accessors.
75
76 virtual void RenewStepperAndAdjust(G4MagIntegratorStepper* stepper) override;
77
78 void ReSetParameters(G4double safety = 0.9);
79 void SetSafety(G4double valS);
80 // i) sets the exponents (pgrow & pshrnk),
81 // using the current Stepper's order,
82 // ii) sets the safety
83
86 // Modify and Get the Maximum number of Steps that can be
87 // taken for the integration of a single segment -
88 // (ie a single call to AccurateAdvance).
89
92
93 protected:
94
97
100
102
103 private:
104
105 inline void RenewStepperAndAdjustImpl(T* stepper);
106
107 G4int fMaxNoSteps;
108
109 G4int fMaxStepBase;
110 // The (default) maximum number of steps is Base
111 // divided by the order of Stepper
112
113 // Parameters used to grow and shrink trial stepsize.
114 //
115 G4double safety;
116 G4double pshrnk; // exponent for shrinking
117 G4double pgrow; // exponent for growth
118
119 // muximum error values for shrinking / growing (optimisation).
120 //
121 G4double errorConstraintShrink;
122 G4double errorConstraintGrow;
123
124 T* pIntStepper = nullptr;
125};
126
127#include "G4RKIntegrationDriver.icc"
128
129#endif
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double GetSmallestFraction() const
G4RKIntegrationDriver(const G4RKIntegrationDriver &)=delete
virtual void GetDerivatives(const G4FieldTrack &track, G4double dydx[], G4double field[]) const override
void SetSafety(G4double valS)
void ReSetParameters(G4double safety=0.9)
G4double GetPgrow() const
virtual void RenewStepperAndAdjust(G4MagIntegratorStepper *stepper) override
G4RKIntegrationDriver(T *stepper)
void SetMaxNoSteps(G4int val)
G4double GetPshrnk() const
G4int GetMaxNoSteps() const
G4double GrowStepSize(G4double h, G4double error) const
virtual const T * GetStepper() const override
G4double GrowStepSize2(G4double h, G4double error2) const
void SetSmallestFraction(G4double val)
G4double ShrinkStepSize(G4double h, G4double error) const
G4double ShrinkStepSize2(G4double h, G4double error2) const
G4double GetSafety() const
virtual void GetDerivatives(const G4FieldTrack &track, G4double dydx[]) const override
virtual void StreamInfo(std::ostream &os) const override
virtual G4double ComputeNewStepSize(G4double errMaxNorm, G4double hstepCurrent) override final
virtual T * GetStepper() override
G4RKIntegrationDriver & operator=(const G4RKIntegrationDriver &)=delete
virtual G4EquationOfMotion * GetEquationOfMotion() override
virtual void SetEquationOfMotion(G4EquationOfMotion *equation) override