Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RKG3_Stepper.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// G4RKG3_Stepper
27//
28// Class description:
29//
30// Runga-Kutta integrator stepper from Geant-3.
31
32// Authors: John Apostolakis & Vladimir Grichine (CERN), 30.01.1997
33// -------------------------------------------------------------------
34#ifndef G4RKG3_STEPPER_HH
35#define G4RKG3_STEPPER_HH
36
37#include "G4Types.hh"
39#include "G4ThreeVector.hh"
40
41class G4Mag_EqRhs;
42
43/**
44 * @brief G4RKG3_Stepper implements a Runga-Kutta integrator stepper
45 * used in Geant-3.
46 */
47
49{
50 public:
51
52 /**
53 * Constructor for G4RKG3_Stepper.
54 * @param[in] EqRhs Pointer to the provided equation of motion.
55 */
57
58 /**
59 * Default Destructor.
60 */
61 ~G4RKG3_Stepper() override = default;
62
63 /**
64 * Copy constructor and assignment operator not allowed.
65 */
68
69 /**
70 * The stepper for the Runge Kutta integration.
71 * Method provided, even if less efficient.
72 * The stepsize is fixed, with the step size given by 'h'.
73 * Integrates ODE starting values yInput[0 to 6].
74 * Outputs yOut[] and its estimated error yErr[].
75 * @param[in] yIn 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 */
81 void Stepper( const G4double yIn[],
82 const G4double dydx[],
83 G4double h,
84 G4double yOut[],
85 G4double yErr[] ) override;
86
87 /**
88 * Returns the distance from chord line.
89 */
90 G4double DistChord() const override ;
91
92 /**
93 * Integrator of Runge-Kutta Stepper from Geant-3 with only two field
94 * evaluation per Step. It is used in propagating the initial Step
95 * by small substeps after solution error and delta geometry
96 * considerations. B[3] is magnetic field which is passed from substep
97 * to substep.
98 */
99 void StepNoErr( const G4double tIn[8],
100 const G4double dydx[6],
101 G4double Step,
102 G4double tOut[8],
103 G4double B[3] );
104 /**
105 * Integrator for Runge-Kutta from Geant-3 with evaluation of error in
106 * solution and delta geometry based on naive similarity with the case
107 * of uniform magnetic field.
108 * B1[3] is in input and is the first magnetic field values
109 * B2[3] is the output and is the final magnetic field values.
110 */
111 void StepWithEst( const G4double tIn[8],
112 const G4double dydx[6],
113 G4double Step,
114 G4double tOut[8],
116 G4double& beta2,
117 const G4double B1[3],
118 G4double B2[3] );
119
120 /**
121 * Returns the order, 4, of integration.
122 */
123 inline G4int IntegratorOrder() const override { return 4; }
124
125 /**
126 * Returns the stepper type-ID, "kRKG3Stepper".
127 */
128 inline G4StepperType StepperType() const override { return kRKG3Stepper; }
129
130 private:
131
132 G4ThreeVector fyInitial,
133 fyMidPoint,
134 fyFinal;
135 G4ThreeVector fpInitial;
136 G4ThreeVector BfldIn;
137 G4double hStep = 0.0;
138};
139
140#endif
G4double B(G4double temperature)
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kRKG3Stepper
G4RKG3_Stepper.
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double alpha2
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4Mag_EqRhs is the "standard" equation of motion of a particle in a pure magnetic field.
G4double DistChord() const override
G4int IntegratorOrder() const override
~G4RKG3_Stepper() override=default
void StepWithEst(const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double &alpha2, G4double &beta2, const G4double B1[3], G4double B2[3])
G4RKG3_Stepper & operator=(const G4RKG3_Stepper &)=delete
G4StepperType StepperType() const override
G4RKG3_Stepper(const G4RKG3_Stepper &)=delete
void Stepper(const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double yErr[]) override
void StepNoErr(const G4double tIn[8], const G4double dydx[6], G4double Step, G4double tOut[8], G4double B[3])
G4RKG3_Stepper(G4Mag_EqRhs *EqRhs)