Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DoLoMcPriRK34.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// G4DoLoMcPriRK34
27//
28// Class description:
29//
30// Dormand-Lockyer-McGorrigan-Prince-6-3-4 non-FSAL method
31// ( 6 stage, 3rd & 4th order embedded RK method )
32
33// Author: Somnath Banerjee (CERN, Google Summer of Code 2015), 07.07.2015
34// Supervision: John Apostolakis (CERN)
35// --------------------------------------------------------------------
36#ifndef DOLO_MCPRI_RK34_HH
37#define DOLO_MCPRI_RK34_HH
38
40
41/**
42 * @brief G4DoLoMcPriRK34 implements the Dormand-Lockyer-McGorrigan-Prince-6-3-4
43 * non-FSAL method ( 6 stage, 3rd & 4th order embedded Runge-Kutta method ).
44 */
45
47{
48 public:
49
50 /**
51 * Constructor for G4DoLoMcPriRK34.
52 * @param[in] EqRhs Pointer to the provided equation of motion.
53 * @param[in] numberOfVariables The number of integration variables.
54 * @param[in] primary Flag for initialisation of the auxiliary stepper.
55 */
57 G4int numberOfVariables = 6,
58 G4bool primary = true );
59
60 /**
61 * Destructor.
62 */
63 ~G4DoLoMcPriRK34() override;
64
65 /**
66 * Copy constructor and assignment operator not allowed.
67 */
70
71 /**
72 * The stepper for the Runge Kutta integration.
73 * The stepsize is fixed, with the step size given by 'h'.
74 * Integrates ODE starting values y[0 to 6].
75 * Outputs yout[] and its estimated error yerr[].
76 * @param[in] y Starting values array of integration variables.
77 * @param[in] dydx Derivatives array.
78 * @param[in] h The given step size.
79 * @param[out] yout Integration output.
80 * @param[out] yerr The estimated error.
81 */
82 void Stepper( const G4double y[],
83 const G4double dydx[],
84 G4double h,
85 G4double yout[],
86 G4double yerr[] ) override ;
87
88 /**
89 * Interface method for interpolation setup. Does nothing here.
90 */
91 inline void SetupInterpolation() {}
92
93 /**
94 * Calculates the output at the tau fraction of Step.
95 * @param[in] yInput Starting values array of integration variables.
96 * @param[in] dydx Derivatives array.
97 * @param[in] Step The given step size.
98 * @param[out] yOut Interpolation output.
99 * @param[out] tau Fraction of step.
100 */
101 void Interpolate( const G4double yInput[],
102 const G4double dydx[],
103 const G4double Step,
104 G4double yOut[],
105 G4double tau );
106 void Interpolate( G4double tau,
107 G4double yOut[]);
108
109 /**
110 * Returns the distance from chord line.
111 */
112 G4double DistChord() const override;
113
114 /**
115 * Returns the order, 3, of integration.
116 */
117 inline G4int IntegratorOrder() const override { return 3; }
118
119
120 /**
121 * Returns the stepper type-ID, "kDoLoMcPriRK34".
122 */
123 inline G4StepperType StepperType() const override { return kDoLoMcPriRK34; }
124
125 private :
126
127 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *yTemp, *yIn;
128
129 G4double fLastStepLength = -1.0;
130
131 /** For DistChord calculations. */
132 G4double *fLastInitialVector, *fLastFinalVector,
133 *fLastDyDx, *fMidVector, *fMidError;
134
135 G4DoLoMcPriRK34* fAuxStepper = nullptr;
136};
137
138#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kDoLoMcPriRK34
G4DoLoMcPriRK34.
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4DoLoMcPriRK34 implements the Dormand-Lockyer-McGorrigan-Prince-6-3-4 non-FSAL method ( 6 stage,...
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
void Interpolate(const G4double yInput[], const G4double dydx[], const G4double Step, G4double yOut[], G4double tau)
G4DoLoMcPriRK34 & operator=(const G4DoLoMcPriRK34 &)=delete
G4int IntegratorOrder() const override
G4DoLoMcPriRK34(G4EquationOfMotion *EqRhs, G4int numberOfVariables=6, G4bool primary=true)
~G4DoLoMcPriRK34() override
G4DoLoMcPriRK34(const G4DoLoMcPriRK34 &)=delete
G4double DistChord() const override
G4StepperType StepperType() const override
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)