Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BogackiShampine45.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// G4BogackiShampine45
27//
28// Class description:
29//
30// An implementation of the embedded RK method from the following paper
31// by P. Bogacki and L. F. Shampine:
32// "An efficient Runge-Kutta (4,5) pair"
33// Comput. Math. with Appl., vol. 32, no. 6, pp. 15-28, Sep. 1996.
34//
35// An interpolation method provides the value of an intermediate
36// point in a step -- if a step was sucessful.
37//
38// This version can provide the FSAL property of the method,
39// which allows the reuse of the last derivative in the next step,
40// but only by using the additional method GetLastDyDx() (an alternative
41// interface for simpler use of FSAL is under development).
42
43// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 25.05.2015
44// Supervision: John Apostolakis (CERN)
45// --------------------------------------------------------------------
46#ifndef BOGACKI_SHAMPINE_45_HH
47#define BOGACKI_SHAMPINE_45_HH
48
50
51/**
52 * @brief G4BogackiShampine45 is an integrator of particle's equation of
53 * motion based on the Bogacki-Shampine method with FSAL property, allowing
54 * the reuse of the last derivative in the next step.
55 * This Stepper provides 'dense output'. After a successful step, it is
56 * possible to obtain an estimate of the value of the function at an
57 * intermediate point of the interval. This requires only two additional
58 * evaluations of the derivative (and thus the field).
59 */
60
62{
63 public:
64
65 /**
66 * Constructor for G4BogackiShampine45.
67 * @param[in] EqRhs Pointer to the provided equation of motion.
68 * @param[in] numberOfVariables The number of integration variables.
69 * @param[in] primary Flag for initialisation of the auxiliary stepper.
70 */
72 G4int numberOfVariables = 6,
73 G4bool primary = true);
74
75 /**
76 * Destructor.
77 */
78 ~G4BogackiShampine45() override;
79
80 /**
81 * Copy constructor and assignment operator not allowed.
82 */
85
86 /**
87 * The stepper for the Runge Kutta integration.
88 * The stepsize is fixed, with the step size given by 'h'.
89 * Integrates ODE starting values y[0 to 6].
90 * Outputs yout[] and its estimated error yerr[].
91 * @param[in] y Starting values array of integration variables.
92 * @param[in] dydx Derivatives array.
93 * @param[in] h The given step size.
94 * @param[out] yout Integration output.
95 * @param[out] yerr The estimated error.
96 */
97 void Stepper( const G4double y[],
98 const G4double dydx[],
99 G4double h,
100 G4double yout[],
101 G4double yerr[] ) override;
102
103 /**
104 * Setup all coefficients for interpolation.
105 */
106 void SetupInterpolationHigh(); // ( yInput, dydx, Step);
108
109 /**
110 * Calculates the output at the tau fraction of step.
111 * @param[in] tau The tau fraction of the step.
112 * @param[out] yOut Interpolation output.
113 */
114 void InterpolateHigh( G4double tau, G4double yOut[] ) const;
115 inline void Interpolate( G4double tau,
116 G4double yOut[] ) { InterpolateHigh( tau, yOut); }
117
118 /**
119 * Returns the distance from chord line.
120 */
121 G4double DistChord() const override;
122
123 /**
124 * Returns the order, 4, of integration.
125 */
126 inline G4int IntegratorOrder() const override { return 4; }
127
128 /**
129 * Returns the stepper type-ID, "kBogackiShampine45".
130 */
131 inline G4StepperType StepperType() const override { return kBogackiShampine45; }
132
133 /**
134 * Acccessor for dydx array.
135 */
136 void GetLastDydx( G4double dyDxLast[] );
137
138 /**
139 * Initialises the values of the bi[][] array.
140 */
141 void PrepareConstants();
142
143 private:
144
145 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *ak8,
146 *ak9, *ak10, *ak11, *yTemp, *yIn;
147
148 G4double *p[6];
149
150 G4double fLastStepLength = -1.0;
151
152 /** For DistChord() calculations. */
153 G4double *fLastInitialVector, *fLastFinalVector, *fLastDyDx,
154 *fMidVector, *fMidError;
155
156 G4BogackiShampine45* fAuxStepper = nullptr;
157
158 /** For chord - until interpolation is proven. */
159 G4bool fPreparedInterpolation = false;
160
161 /** Class constants. */
162 static G4bool fPreparedConstants;
163 static G4double bi[12][7];
164};
165
166#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kBogackiShampine45
G4BogackiShampine45.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4BogackiShampine45 is an integrator of particle's equation of motion based on the Bogacki-Shampine m...
G4BogackiShampine45(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
G4double DistChord() const override
G4BogackiShampine45(const G4BogackiShampine45 &)=delete
void GetLastDydx(G4double dyDxLast[])
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
void Interpolate(G4double tau, G4double yOut[])
G4StepperType StepperType() const override
G4int IntegratorOrder() const override
void InterpolateHigh(G4double tau, G4double yOut[]) const
G4BogackiShampine45 & operator=(const G4BogackiShampine45 &)=delete
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)