Geant4  10.03
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 // $Id: G4MagIntegratorStepper.hh 97387 2016-06-02 10:03:42Z gcosmo $
27 //
28 //
29 // class G4MagIntegratorStepper
30 //
31 // Class description:
32 //
33 // Abstract base class for integrator of particle's equation of motion,
34 // used in tracking in space dependent magnetic field
35 //
36 // A Stepper must integrate over NumberOfVariables elements,
37 // and also copy (from input to output) any of NoStateVariables
38 // not included in the NumberOfVariables.
39 //
40 // So it is expected that NoStateVariables >= NumberOfVariables
41 
42 // History:
43 // - 15.01.97 J. Apostolakis (J.Apostolakis@cern.ch)
44 // --------------------------------------------------------------------
45 
46 #ifndef G4MAGIntegratorSTEPPER
47 #define G4MAGIntegratorSTEPPER
48 
49 #include "G4Types.hh"
50 #include "G4EquationOfMotion.hh"
51 
53 {
54  public: // with description
55 
57  G4int numIntegrationVariables,
58  G4int numStateVariables=12,
59  bool isFSAL= false
60  // , G4int methodOrder
61  );
62  virtual ~G4MagIntegratorStepper();
63  // Constructor and destructor. No actions.
64 
65  virtual void Stepper( const G4double y[],
66  const G4double dydx[],
67  G4double h,
68  G4double yout[],
69  G4double yerr[] ) = 0 ;
70  // The stepper for the Runge Kutta integration.
71  // The stepsize is fixed, with the Step size given by h.
72  // Integrates ODE starting values y[0 to 6].
73  // Outputs yout[] and its estimated error yerr[].
74 
75  virtual G4double DistChord() const = 0;
76  // Estimate the maximum distance of a chord from the true path
77  // over the segment last integrated.
78 
79  virtual void ComputeRightHandSide( const G4double y[], G4double dydx[] );
80  // Must compute the RightHandSide as in the method below
81  // Optionally can cache the input y[] and the dydx[] values computed.
82 
83  inline void NormaliseTangentVector( G4double vec[6] );
84  // Simple utility function to (re)normalise 'unit velocity' vector.
85 
86  inline void NormalisePolarizationVector( G4double vec[12] );
87  // Simple utility function to (re)normalise 'unit spin' vector.
88 
89  inline void RightHandSide( const double y[], double dydx[] );
90  // Utility method to supply the standard Evaluation of the
91  // Right Hand side of the associated equation.
92 
93  inline G4int GetNumberOfVariables() const;
94  // Get the number of variables that the stepper will integrate over.
95 
96  inline G4int GetNumberOfStateVariables() const;
97  // Get the number of variables of state variables (>= above, integration)
98 
99  virtual G4int IntegratorOrder() const = 0;
100  // Returns the order of the integrator
101  // i.e. its error behaviour is of the order O(h^order).
102 
104  // Replacement method - using new data member
105 
107  // As some steppers (eg RKG3) require other methods of Eq_Rhs
108  // this function allows for access to them.
109  inline void SetEquationOfMotion(G4EquationOfMotion* newEquation);
110 
111  inline unsigned long GetfNoRHSCalls(){ return fNoRHSCalls; }
112  // void IncrementRHSCalls() { fNoRHSCalls++; }
113  inline void ResetfNORHSCalls(){ fNoRHSCalls = 0; }
114  // Count number of calls to RHS method(s)
115 
116  bool IsFSAL() { return fIsFSAL; }
117 
118  protected:
119  void SetIntegrationOrder(int order) { fIntegrationOrder= order; }
120  void SetFSAL( bool flag= true) { fIsFSAL= flag; }
121 
122  private:
123 
126  // Private copy constructor and assignment operator.
127 
128  private:
129 
131  const G4int fNoIntegrationVariables; // Number of Variables in integration
132  const G4int fNoStateVariables; // Number required for FieldTrack
133  // const G4int fNumberOfVariables;
134 
135  // Counter for calls to RHS method
136  mutable unsigned long fNoRHSCalls;
137 
138  // Parameters of a RK method -- must be shared by all steppers of a type
139  // -- Invariants for a class
140  /* const */ int fIntegrationOrder; // all ClassicalRK4 steppers are 4th order
141  /* const */ bool fIsFSAL; // Depends on RK method & implementation
142 };
143 
145 
146 #endif /* G4MAGIntegratorSTEPPER */
void SetEquationOfMotion(G4EquationOfMotion *newEquation)
virtual void ComputeRightHandSide(const G4double y[], G4double dydx[])
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])=0
void NormaliseTangentVector(G4double vec[6])
int G4int
Definition: G4Types.hh:78
G4int GetNumberOfVariables() const
G4EquationOfMotion * GetEquationOfMotion()
G4EquationOfMotion * fEquation_Rhs
virtual G4int IntegratorOrder() const =0
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, bool isFSAL=false)
G4MagIntegratorStepper & operator=(const G4MagIntegratorStepper &)
virtual G4double DistChord() const =0
void SetFSAL(bool flag=true)
G4int GetNumberOfStateVariables() const
void RightHandSide(const double y[], double dydx[])
void NormalisePolarizationVector(G4double vec[12])
double G4double
Definition: G4Types.hh:76