Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TMagErrorStepper.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// G4TMagErrorStepper
27//
28// Class description:
29//
30// Templated version of G4MagErrorStepper.
31// Adapted from G4G4TMagErrorStepper class.
32
33// Author: Josh Xie (CERN, Google Summer of Code 2014), June 2014
34// Supervisors: Sandro Wenzel, John Apostolakis (CERN)
35// --------------------------------------------------------------------
36#ifndef G4TMAG_ERROR_STEPPER_HH
37#define G4TMAG_ERROR_STEPPER_HH
38
39#include "G4Types.hh"
41#include "G4ThreeVector.hh"
42#include "G4LineSection.hh"
43
44/**
45 * @brief G4TMagErrorStepper is a templated version of G4MagErrorStepper.
46 */
47
48template <class T_Stepper, class T_Equation, unsigned int N>
50{
51 public:
52
53 G4TMagErrorStepper(T_Equation* EqRhs, G4int numberOfVariables,
54 G4int numStateVariables = 12)
55 : G4MagIntegratorStepper(EqRhs, numberOfVariables, numStateVariables)
56 , fEquation_Rhs(EqRhs) { ; }
57
58 virtual ~G4TMagErrorStepper() = default;
59
62
63 inline void RightHandSide(G4double y[], G4double dydx[])
64 {
65 fEquation_Rhs->T_Equation::RightHandSide(y, dydx);
66 }
67
68 inline void Stepper(const G4double yInput[], const G4double dydx[],
69 G4double hstep, G4double yOutput[], G4double yError[]) override final;
70
71 inline G4double DistChord() const override final;
72 G4StepperType StepperType() const override { return kTMagErrorStepper; }
73
74 private:
75
76 // STATE
77 G4ThreeVector fInitialPoint, fMidPoint, fFinalPoint;
78 // Data stored in order to find the chord
79
80 // Dependent Objects, owned --- part of the STATE
81 G4double yInitial[N < 8 ? 8 : N];
82 G4double yMiddle[N < 8 ? 8 : N];
83 G4double dydxMid[N < 8 ? 8 : N];
84 G4double yOneStep[N < 8 ? 8 : N];
85 // The following arrays are used only for temporary storage
86 // they are allocated at the class level only for efficiency -
87 // so that calls to new and delete are not made in Stepper().
88
89 T_Equation* fEquation_Rhs;
90};
91
92// ------------ Implementation -----------------------
93
94template <class T_Stepper, class T_Equation, unsigned int N >
96Stepper(const G4double yInput[],
97 const G4double dydx[],
98 G4double hstep,
99 G4double yOutput[],
100 G4double yError[])
101// The stepper for the Runge Kutta integration. The stepsize
102// is fixed, with the Step size given by hstep.
103// Integrates ODE starting values y[0 to N].
104// Outputs yout[] and its estimated error yerr[].
105{
106 const unsigned int maxvar = GetNumberOfStateVariables();
107
108 // Saving yInput because yInput and yOutput can be aliases for same array
109 for(unsigned int i = 0; i < N; ++i)
110 yInitial[i] = yInput[i];
111 yInitial[7] =
112 yInput[7]; // Copy the time in case ... even if not really needed
113 yMiddle[7] = yInput[7]; // Copy the time from initial value
114 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
115 // yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4
116 for(unsigned int i = N; i < maxvar; ++i)
117 yOutput[i] = yInput[i];
118
119 G4double halfStep = hstep * 0.5;
120
121 // Do two half steps
122
123 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, halfStep,
124 yMiddle);
125 this->RightHandSide(yMiddle, dydxMid);
126 static_cast<T_Stepper*>(this)->DumbStepper(yMiddle, dydxMid, halfStep,
127 yOutput);
128
129 // Store midpoint, chord calculation
130
131 fMidPoint = G4ThreeVector(yMiddle[0], yMiddle[1], yMiddle[2]);
132
133 // Do a full Step
134 static_cast<T_Stepper*>(this)->DumbStepper(yInitial, dydx, hstep, yOneStep);
135 for(unsigned int i = 0; i < N; ++i)
136 {
137 yError[i] = yOutput[i] - yOneStep[i];
138 yOutput[i] +=
139 yError[i] *
140 T_Stepper::IntegratorCorrection; // Provides accuracy increased
141 // by 1 order via the
142 // Richardson Extrapolation
143 }
144
145 fInitialPoint = G4ThreeVector(yInitial[0], yInitial[1], yInitial[2]);
146 fFinalPoint = G4ThreeVector(yOutput[0], yOutput[1], yOutput[2]);
147
148 return;
149}
150
151
152template <class T_Stepper, class T_Equation, unsigned int N >
153inline G4double
155{
156 // Estimate the maximum distance from the curve to the chord
157 //
158 // We estimate this using the distance of the midpoint to
159 // chord (the line between
160 //
161 // Method below is good only for angle deviations < 2 pi,
162 // This restriction should not a problem for the Runge cutta methods,
163 // which generally cannot integrate accurately for large angle deviations.
164 G4double distLine, distChord;
165
166 if(fInitialPoint != fFinalPoint)
167 {
168 distLine = G4LineSection::Distline(fMidPoint, fInitialPoint, fFinalPoint);
169 // This is a class method that gives distance of Mid
170 // from the Chord between the Initial and Final points.
171
172 distChord = distLine;
173 }
174 else
175 {
176 distChord = (fMidPoint - fInitialPoint).mag();
177 }
178
179 return distChord;
180}
181
182#endif
G4StepperType
G4StepperType defines the available integrator of particle's equation of motion in Geant4.
@ kTMagErrorStepper
G4TMagErrorStepper.
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)
G4int GetNumberOfStateVariables() const
G4TMagErrorStepper(const G4TMagErrorStepper &)=delete
G4double DistChord() const override final
G4TMagErrorStepper(T_Equation *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
G4TMagErrorStepper & operator=(const G4TMagErrorStepper &)=delete
virtual ~G4TMagErrorStepper()=default
void RightHandSide(G4double y[], G4double dydx[])
void Stepper(const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override final
G4StepperType StepperType() const override
#define N
Definition crc32.c:57