Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4BogackiShampine23.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// G4BogackiShampine23
27//
28// Class description:
29//
30// Bogacki-Shampine - 4 - 3(2) non-FSAL implementation
31//
32// An implementation of the embedded RK method from the paper
33// [1] P. Bogacki and L. F. Shampine,
34// "A 3(2) pair of Runge - Kutta formulas"
35// Appl. Math. Lett., vol. 2, no. 4, pp. 321-325, Jan. 1989.
36//
37// This version does not utilise the FSAL property of the method,
38// which would allow the reuse of the last derivative in the next step.
39// (Alternative FSAL implementation created with revised interface)
40
41// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 20.05.2015
42// Supervision: John Apostolakis (CERN)
43// --------------------------------------------------------------------
44#ifndef G4BOGACKI_SHAMPINE23_HH
45#define G4BOGACKI_SHAMPINE23_HH
46
48#include "G4FieldTrack.hh"
49
50/**
51 * @brief G4BogackiShampine23 is an integrator of particle's equation of
52 * motion based on the Bogacki-Shampine non-FSAL implementation.
53 */
54
56{
57 public:
58
59 /**
60 * Constructor for G4BogackiShampine23.
61 * @param[in] EqRhs Pointer to the provided equation of motion.
62 * @param[in] numberOfVariables The number of integration variables.
63 */
65 G4int numberOfVariables = 6);
66
67 /**
68 * Default Destructor.
69 */
70 ~G4BogackiShampine23() override = default;
71
72 /**
73 * Copy constructor and assignment operator not allowed.
74 */
77
78 /**
79 * The stepper for the Runge Kutta integration.
80 * The stepsize is fixed, with the step size given by 'hstep'.
81 * Integrates ODE starting values yInput[0 to 6].
82 * Outputs yOutput[] and its estimated error yError[].
83 * @param[in] yInput Starting values array of integration variables.
84 * @param[in] dydx Derivatives array.
85 * @param[in] hstep The given step size.
86 * @param[out] yOutput Integration output.
87 * @param[out] yError The estimated error.
88 */
89 void Stepper(const G4double yInput[],
90 const G4double dydx[],
91 G4double hstep,
92 G4double yOutput[],
93 G4double yError[]) override;
94
95 /**
96 * Same as the Stepper() function above, with dydx also in ouput.
97 * @param[in] yInput Starting values array of integration variables.
98 * @param[in] dydx Derivatives array.
99 * @param[in] hstep The given step size.
100 * @param[out] yOutput Integration output.
101 * @param[out] yError The estimated error.
102 * @param[out] dydxOutput dydx in output.
103 */
104 void Stepper(const G4double yInput[],
105 const G4double dydx[],
106 G4double hstep,
107 G4double yOutput[],
108 G4double yError[],
109 G4double dydxOutput[]);
110
111 /**
112 * Returns the distance from chord line.
113 */
114 G4double DistChord() const override;
115
116 /**
117 * Returns the order, 3, of integration.
118 */
119 inline G4int IntegratorOrder() const override { return 3; }
120
121 /**
122 * Returns the stepper type-ID, "kBogackiShampine23".
123 */
124 inline G4StepperType StepperType() const override { return kBogackiShampine23; }
125
126 private:
127
128 /**
129 * Utility method used in Stepper() for computing the actual step.
130 */
131 void makeStep(const G4double yInput[],
132 const G4double dydx[],
133 const G4double hstep,
134 G4double yOutput[],
135 G4double* dydxOutput = nullptr,
136 G4double* yError = nullptr) const;
137
138 private:
139
143 fdydxOut[G4FieldTrack::ncompSVEC];
144 G4double fhstep = -1.0;
145};
146
147#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kBogackiShampine23
G4BogackiShampine23.
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4double DistChord() const override
~G4BogackiShampine23() override=default
G4BogackiShampine23(const G4BogackiShampine23 &)=delete
G4BogackiShampine23 & operator=(const G4BogackiShampine23 &)=delete
G4StepperType StepperType() const override
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override
G4int IntegratorOrder() const override
G4BogackiShampine23(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6)
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)