Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FSALDormandPrince745.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// G4FSALDormandPrince745
27//
28// Class description:
29//
30// DormandPrince7 - 5(4) FSAL stepper
31
32// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 25.05.2015
33// Supervision: John Apostolakis (CERN)
34// --------------------------------------------------------------------
35#ifndef G4FSALDORMANDPRINCE745_HH
36#define G4FSALDORMANDPRINCE745_HH
37
39
40/**
41 * @brief G4FSALDormandPrince745 is an integrator of particle's equation of
42 * motion based on the DormandPrince7 - 5(4) FSAL implementation.
43 */
44
46{
47 public:
48
49 /**
50 * Constructor for G4FSALDormandPrince745.
51 * @param[in] EqRhs Pointer to the provided equation of motion.
52 * @param[in] numberOfVariables The number of integration variables.
53 * @param[in] primary Flag for initialisation of the auxiliary stepper.
54 */
56 G4int numberOfVariables = 6,
57 G4bool primary = true);
58
59 /**
60 * Destructor.
61 */
62 ~G4FSALDormandPrince745() override;
63
64 /**
65 * Copy constructor and assignment operator not allowed.
66 */
69
70 /**
71 * The stepper for the Runge Kutta integration.
72 * The stepsize is fixed, with the step size given by 'h'.
73 * Integrates ODE starting values y[0 to 6].
74 * Outputs yout[] and its estimated error yerr[].
75 * @param[in] y Starting values array of integration variables.
76 * @param[in] dydx Derivatives array.
77 * @param[in] h The given step size.
78 * @param[out] yout Integration output.
79 * @param[out] yerr The estimated error.
80 * @param[out] nextDydx Last derivatives array for the next step.
81 */
82 void Stepper( const G4double y[],
83 const G4double dydx[],
84 G4double h,
85 G4double yout[],
86 G4double yerr[],
87 G4double nextDydx[]) override ;
88
89 /**
90 * Calculates the output at the tau fraction of step.
91 * @param[in] yInput Starting values array of integration variables.
92 * @param[in] dydx Derivatives array.
93 * @param[out] yOut Interpolation output.
94 * @param[in] Step The given step size.
95 * @param[in] tau The tau fraction of the step.
96 */
97 void interpolate( const G4double yInput[],
98 const G4double dydx[],
99 G4double yOut[],
100 G4double Step,
101 G4double tau ) ;
102
103 /**
104 * Calculates the output at the tau fraction of step. Same as above
105 * for higher order interpolant.
106 */
107 void Interpolate( const G4double yInput[],
108 const G4double dydx[],
109 const G4double Step,
110 G4double yOut[],
111 G4double tau );
112
113 /**
114 * Setup method for higher order interpolant.
115 * @param[in] yInput Starting values array of integration variables.
116 * @param[in] dydx Derivatives array.
117 * @param[in] Step The given step size.
118 */
119 void SetupInterpolate( const G4double yInput[],
120 const G4double dydx[],
121 const G4double Step );
122
123 /**
124 * Returns the distance from chord line.
125 */
126 G4double DistChord() const override;
127
128 /**
129 * Returns the order, 4, of integration.
130 */
131 inline G4int IntegratorOrder() const override { return 4; }
132
133 /**
134 * Returns true as this is a FSAL integrator.
135 */
136 inline G4bool isFSAL() const { return true; }
137
138 private:
139
140 /** Working arrays -- used during stepping. */
141 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7,
142 *ak8, *ak9, // For additional stages in the interpolant
143 *yTemp, *yIn;
144
145 /** Only for use with DistChord(). */
146 G4double* pseudoDydx_for_DistChord;
147
148 G4double fLastStepLength = -1.0;
149 G4double *fLastInitialVector, *fLastFinalVector,
150 *fInitialDyDx, *fLastDyDx,
151 *fMidVector, *fMidError; // For DistChord() calculations
152
153 G4FSALDormandPrince745* fAuxStepper = nullptr;
154};
155
156#endif
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4EquationOfMotion is the abstract base class for the right hand size of the equation of motion of a ...
G4FSALDormandPrince745 is an integrator of particle's equation of motion based on the DormandPrince7 ...
G4FSALDormandPrince745(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[], G4double nextDydx[]) override
void interpolate(const G4double yInput[], const G4double dydx[], G4double yOut[], G4double Step, G4double tau)
G4FSALDormandPrince745 & operator=(const G4FSALDormandPrince745 &)=delete
G4FSALDormandPrince745(const G4FSALDormandPrince745 &)=delete
G4double DistChord() const override
void SetupInterpolate(const G4double yInput[], const G4double dydx[], const G4double Step)
G4int IntegratorOrder() const override
G4VFSALIntegrationStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12)