Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4MagIntegratorStepper.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// G4MagIntegratorStepper
27//
28// Class description:
29//
30// Abstract base class for integrator of particle's equation of motion,
31// used in tracking in space dependent magnetic field
32//
33// A Stepper must integrate over NumberOfVariables elements,
34// and also copy (from input to output) any of NoStateVariables
35// not included in the NumberOfVariables.
36//
37// So it is expected that NoStateVariables >= NumberOfVariables
38
39// Author: J.Apostolakis, CERN - 15.01.1997
40// --------------------------------------------------------------------
41#ifndef G4MAGINTEGRATORSTEPPER_HH
42#define G4MAGINTEGRATORSTEPPER_HH
43
44#include "G4Types.hh"
45#include "G4EquationOfMotion.hh"
48
50
52{
53 public: // with description
54
56 G4int numIntegrationVariables,
57 G4int numStateVariables = 12,
58 G4bool isFSAL = false );
59
60 virtual ~G4MagIntegratorStepper() = default;
61 // Constructor and destructor. No actions.
62
65
66 virtual void Stepper( const G4double y[],
67 const G4double dydx[],
68 G4double h,
69 G4double yout[],
70 G4double yerr[] ) = 0;
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
76 virtual G4double DistChord() const = 0;
77 // Estimate the maximum distance of a chord from the true path
78 // over the segment last integrated.
79
80 inline void NormaliseTangentVector( G4double vec[6] );
81 // Simple utility function to (re)normalise 'unit velocity' vector.
82
83 inline void NormalisePolarizationVector( G4double vec[12] );
84 // Simple utility function to (re)normalise 'unit spin' vector.
85
86 inline void RightHandSide( const G4double y[], G4double dydx[] ) const;
87 // Utility method to supply the standard Evaluation of the
88 // Right Hand side of the associated equation.
89
90 inline void RightHandSide( const G4double y[],
91 G4double dydx[],
92 G4double field[] ) const;
93 // Calculate dydx and field at point y.
94
96 // Get the number of variables that the stepper will integrate over.
97
99 // Get the number of variables of state variables (>= above, integration)
100
101 virtual G4int IntegratorOrder() const = 0;
102 // Returns the order of the integrator
103 // i.e. its error behaviour is of the order O(h^order).
104
106 // Replacement method - using new data member
107
110 // As some steppers (eg RKG3) require other methods of Eq_Rhs
111 // this function allows for access to them.
112
113 inline void SetEquationOfMotion(G4EquationOfMotion* newEquation);
114
115 inline unsigned long GetfNoRHSCalls();
116 inline void ResetfNORHSCalls();
117 // Count number of calls to RHS method(s)
118
119 inline G4bool IsFSAL() const;
120
121 // TODO - QSS
122 inline G4bool isQSS() const { return fIsQSS; }
123 void SetIsQSS(G4bool val){ fIsQSS= val;}
124
125protected:
126
127 inline void SetIntegrationOrder(G4int order);
128 inline void SetFSAL(G4bool flag = true);
129
130 private:
131
132 G4EquationOfMotion* fEquation_Rhs = nullptr;
133 const G4int fNoIntegrationVariables = 0; // Variables in integration
134 const G4int fNoStateVariables = 0; // Number required for FieldTrack
135
136 mutable unsigned long fNoRHSCalls = 0UL;
137 // Counter for calls to RHS method
138
139 // Parameters of a RK method -- must be shared by all steppers of a type
140 // -- Invariants for a class
141
142 G4int fIntegrationOrder = -1; // must be set by stepper !!!
143 // All ClassicalRK4 steppers are 4th order
144 G4bool fIsFSAL = false;
145 // Depends on RK method & implementation
146 G4bool fIsQSS = false;
147
148};
149
150#include "G4MagIntegratorStepper.icc"
151
152#endif /* G4MAGIntegratorSTEPPER */
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
void NormalisePolarizationVector(G4double vec[12])
G4bool IsFSAL() const
unsigned long GetfNoRHSCalls()
G4EquationOfMotion * GetEquationOfMotion()
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const
void NormaliseTangentVector(G4double vec[6])
virtual G4double DistChord() const =0
virtual ~G4MagIntegratorStepper()=default
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
G4MagIntegratorStepper & operator=(const G4MagIntegratorStepper &)=delete
void SetIntegrationOrder(G4int order)
void SetEquationOfMotion(G4EquationOfMotion *newEquation)
const G4EquationOfMotion * GetEquationOfMotion() const
G4MagIntegratorStepper(const G4MagIntegratorStepper &)=delete
void SetFSAL(G4bool flag=true)
void RightHandSide(const G4double y[], G4double dydx[], G4double field[]) const
G4int GetNumberOfStateVariables() const
virtual G4int IntegratorOrder() const =0