Geant4-11
Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
G4ClassicalRK4 Class Reference

#include <G4ClassicalRK4.hh>

Inheritance diagram for G4ClassicalRK4:
G4MagErrorStepper G4MagIntegratorStepper

Public Member Functions

G4double DistChord () const
 
void DumbStepper (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[])
 
 G4ClassicalRK4 (const G4ClassicalRK4 &)=delete
 
 G4ClassicalRK4 (G4EquationOfMotion *EquationMotion, G4int numberOfVariables=6)
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
unsigned long GetfNoRHSCalls ()
 
G4int GetNumberOfStateVariables () const
 
G4int GetNumberOfVariables () const
 
G4int IntegrationOrder ()
 
G4int IntegratorOrder () const
 
G4bool IsFSAL () const
 
void NormalisePolarizationVector (G4double vec[12])
 
void NormaliseTangentVector (G4double vec[6])
 
G4ClassicalRK4operator= (const G4ClassicalRK4 &)=delete
 
void ResetfNORHSCalls ()
 
void RightHandSide (const G4double y[], G4double dydx[]) const
 
void RightHandSide (const G4double y[], G4double dydx[], G4double field[]) const
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
 
 ~G4ClassicalRK4 ()
 

Protected Member Functions

void SetFSAL (G4bool flag=true)
 
void SetIntegrationOrder (G4int order)
 

Private Member Functions

void StepWithEst (const G4double yIn[], const G4double dydx[], G4double h, G4double yOut[], G4double &alpha2, G4double &beta2, const G4double B1[], G4double B2[])
 

Private Attributes

G4doubledydxm
 
G4doubledydxMid
 
G4doubledydxt
 
G4EquationOfMotionfEquation_Rhs = nullptr
 
G4ThreeVector fFinalPoint
 
G4ThreeVector fInitialPoint
 
G4int fIntegrationOrder = -1
 
G4bool fIsFSAL = false
 
G4ThreeVector fMidPoint
 
const G4int fNoIntegrationVariables = 0
 
unsigned long fNoRHSCalls = 0UL
 
const G4int fNoStateVariables = 0
 
G4doubleyInitial
 
G4doubleyMiddle
 
G4doubleyOneStep
 
G4doubleyt
 

Detailed Description

Definition at line 40 of file G4ClassicalRK4.hh.

Constructor & Destructor Documentation

◆ G4ClassicalRK4() [1/2]

G4ClassicalRK4::G4ClassicalRK4 ( G4EquationOfMotion EquationMotion,
G4int  numberOfVariables = 6 
)

Definition at line 38 of file G4ClassicalRK4.cc.

40 : G4MagErrorStepper(EqRhs, numberOfVariables)
41{
42 unsigned int noVariables= std::max(numberOfVariables,8); // For Time .. 7+1
43
44 dydxm = new G4double[noVariables];
45 dydxt = new G4double[noVariables];
46 yt = new G4double[noVariables];
47}
double G4double
Definition: G4Types.hh:83
G4double * dydxm
G4double * dydxt
G4MagErrorStepper(G4EquationOfMotion *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References dydxm, dydxt, G4INCL::Math::max(), and yt.

◆ ~G4ClassicalRK4()

G4ClassicalRK4::~G4ClassicalRK4 ( )

Definition at line 53 of file G4ClassicalRK4.cc.

54{
55 delete [] dydxm;
56 delete [] dydxt;
57 delete [] yt;
58}

References dydxm, dydxt, and yt.

◆ G4ClassicalRK4() [2/2]

G4ClassicalRK4::G4ClassicalRK4 ( const G4ClassicalRK4 )
delete

Member Function Documentation

◆ DistChord()

G4double G4MagErrorStepper::DistChord ( ) const
virtualinherited

Implements G4MagIntegratorStepper.

Definition at line 100 of file G4MagErrorStepper.cc.

101{
102 // Estimate the maximum distance from the curve to the chord
103 //
104 // We estimate this using the distance of the midpoint to
105 // chord (the line between
106 //
107 // Method below is good only for angle deviations < 2 pi,
108 // This restriction should not a problem for the Runge cutta methods,
109 // which generally cannot integrate accurately for large angle deviations.
110
111 G4double distLine, distChord;
112
114 {
116 // This is a class method that gives distance of Mid
117 // from the Chord between the Initial and Final points.
118
119 distChord = distLine;
120 }
121 else
122 {
123 distChord = (fMidPoint-fInitialPoint).mag();
124 }
125
126 return distChord;
127}
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
G4ThreeVector fMidPoint
G4ThreeVector fFinalPoint
G4ThreeVector fInitialPoint

References G4LineSection::Distline(), G4MagErrorStepper::fFinalPoint, G4MagErrorStepper::fInitialPoint, and G4MagErrorStepper::fMidPoint.

◆ DumbStepper()

void G4ClassicalRK4::DumbStepper ( const G4double  yIn[],
const G4double  dydx[],
G4double  h,
G4double  yOut[] 
)
virtual

Implements G4MagErrorStepper.

Definition at line 71 of file G4ClassicalRK4.cc.

75{
76 const G4int nvar = GetNumberOfVariables(); // fNumberOfVariables();
77 G4int i;
78 G4double hh = h*0.5, h6 = h/6.0;
79
80 // Initialise time to t0, needed when it is not updated by the integration.
81 // [ Note: Only for time dependent fields (usually electric)
82 // is it neccessary to integrate the time.]
83 yt[7] = yIn[7];
84 yOut[7] = yIn[7];
85
86 for(i=0; i<nvar; ++i)
87 {
88 yt[i] = yIn[i] + hh*dydx[i] ; // 1st Step K1=h*dydx
89 }
90 RightHandSide(yt,dydxt) ; // 2nd Step K2=h*dydxt
91
92 for(i=0; i<nvar; ++i)
93 {
94 yt[i] = yIn[i] + hh*dydxt[i] ;
95 }
96 RightHandSide(yt,dydxm) ; // 3rd Step K3=h*dydxm
97
98 for(i=0; i<nvar; ++i)
99 {
100 yt[i] = yIn[i] + h*dydxm[i] ;
101 dydxm[i] += dydxt[i] ; // now dydxm=(K2+K3)/h
102 }
103 RightHandSide(yt,dydxt) ; // 4th Step K4=h*dydxt
104
105 for(i=0; i<nvar; ++i) // Final RK4 output
106 {
107 yOut[i] = yIn[i]+h6*(dydx[i]+dydxt[i]+2.0*dydxm[i]); //+K1/6+K4/6+(K2+K3)/3
108 }
109 if ( nvar == 12 ) { NormalisePolarizationVector ( yOut ); }
110
111} // end of DumbStepper ....................................................
int G4int
Definition: G4Types.hh:85
void NormalisePolarizationVector(G4double vec[12])
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const

References dydxm, dydxt, G4MagIntegratorStepper::GetNumberOfVariables(), G4MagIntegratorStepper::NormalisePolarizationVector(), G4MagIntegratorStepper::RightHandSide(), and yt.

◆ GetEquationOfMotion() [1/2]

G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( )
inlineinherited

◆ GetEquationOfMotion() [2/2]

const G4EquationOfMotion * G4MagIntegratorStepper::GetEquationOfMotion ( ) const
inlineinherited

◆ GetfNoRHSCalls()

unsigned long G4MagIntegratorStepper::GetfNoRHSCalls ( )
inlineinherited

◆ GetNumberOfStateVariables()

G4int G4MagIntegratorStepper::GetNumberOfStateVariables ( ) const
inlineinherited

◆ GetNumberOfVariables()

G4int G4MagIntegratorStepper::GetNumberOfVariables ( ) const
inlineinherited

◆ IntegrationOrder()

G4int G4MagIntegratorStepper::IntegrationOrder ( )
inlineinherited

◆ IntegratorOrder()

G4int G4ClassicalRK4::IntegratorOrder ( ) const
inlinevirtual

Implements G4MagIntegratorStepper.

Definition at line 70 of file G4ClassicalRK4.hh.

70{ return 4; }

◆ IsFSAL()

G4bool G4MagIntegratorStepper::IsFSAL ( ) const
inlineinherited

◆ NormalisePolarizationVector()

void G4MagIntegratorStepper::NormalisePolarizationVector ( G4double  vec[12])
inlineinherited

◆ NormaliseTangentVector()

void G4MagIntegratorStepper::NormaliseTangentVector ( G4double  vec[6])
inlineinherited

◆ operator=()

G4ClassicalRK4 & G4ClassicalRK4::operator= ( const G4ClassicalRK4 )
delete

◆ ResetfNORHSCalls()

void G4MagIntegratorStepper::ResetfNORHSCalls ( )
inlineinherited

◆ RightHandSide() [1/2]

void G4MagIntegratorStepper::RightHandSide ( const G4double  y[],
G4double  dydx[] 
) const
inlineinherited

◆ RightHandSide() [2/2]

void G4MagIntegratorStepper::RightHandSide ( const G4double  y[],
G4double  dydx[],
G4double  field[] 
) const
inlineinherited

◆ SetEquationOfMotion()

void G4MagIntegratorStepper::SetEquationOfMotion ( G4EquationOfMotion newEquation)
inlineinherited

◆ SetFSAL()

void G4MagIntegratorStepper::SetFSAL ( G4bool  flag = true)
inlineprotectedinherited

◆ SetIntegrationOrder()

void G4MagIntegratorStepper::SetIntegrationOrder ( G4int  order)
inlineprotectedinherited

◆ Stepper()

void G4MagErrorStepper::Stepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[],
G4double  yerr[] 
)
virtualinherited

Implements G4MagIntegratorStepper.

Definition at line 42 of file G4MagErrorStepper.cc.

47{
48 const G4int nvar = GetNumberOfVariables();
49 const G4int maxvar = GetNumberOfStateVariables();
50
51 // correction for Richardson Extrapolation.
52 //
53 G4double correction = 1. / ( (1 << IntegratorOrder()) -1 );
54
55 // Saving yInput because yInput and yOutput can be aliases for same array
56 //
57 for(G4int i=0; i<nvar; ++i)
58 {
59 yInitial[i]=yInput[i];
60 }
61 yInitial[7] = yInput[7]; // Copy the time in case...even if not really needed
62 yMiddle[7] = yInput[7]; // Copy the time from initial value
63 yOneStep[7] = yInput[7]; // As it contributes to final value of yOutput ?
64 // yOutput[7] = yInput[7]; // -> dumb stepper does it too for RK4
65
66 for(G4int i=nvar; i<maxvar; ++i)
67 {
68 yOutput[i]=yInput[i];
69 }
70 // yError[7] = 0.0;
71
72 G4double halfStep = hstep * 0.5;
73
74 // Do two half steps
75 //
76 DumbStepper (yInitial, dydx, halfStep, yMiddle);
78 DumbStepper (yMiddle, dydxMid, halfStep, yOutput);
79
80 // Store midpoint, chord calculation
81 //
83
84 // Do a full Step
85 //
86 DumbStepper(yInitial, dydx, hstep, yOneStep);
87 for(G4int i=0; i<nvar; ++i)
88 {
89 yError [i] = yOutput[i] - yOneStep[i] ;
90 yOutput[i] += yError[i]*correction ;
91 // Provides accuracy increased by 1 order via Richardson Extrapolation
92 }
93
95 fFinalPoint = G4ThreeVector( yOutput[0], yOutput[1], yOutput[2]);
96
97 return;
98}
CLHEP::Hep3Vector G4ThreeVector
virtual void DumbStepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[])=0
G4int GetNumberOfStateVariables() const
virtual G4int IntegratorOrder() const =0

References G4MagErrorStepper::DumbStepper(), G4MagErrorStepper::dydxMid, G4MagErrorStepper::fFinalPoint, G4MagErrorStepper::fInitialPoint, G4MagErrorStepper::fMidPoint, G4MagIntegratorStepper::GetNumberOfStateVariables(), G4MagIntegratorStepper::GetNumberOfVariables(), G4MagIntegratorStepper::IntegratorOrder(), G4MagIntegratorStepper::RightHandSide(), G4MagErrorStepper::yInitial, G4MagErrorStepper::yMiddle, and G4MagErrorStepper::yOneStep.

◆ StepWithEst()

void G4ClassicalRK4::StepWithEst ( const G4double  yIn[],
const G4double  dydx[],
G4double  h,
G4double  yOut[],
G4double alpha2,
G4double beta2,
const G4double  B1[],
G4double  B2[] 
)
private

Definition at line 118 of file G4ClassicalRK4.cc.

126{
127 G4Exception("G4ClassicalRK4::StepWithEst()", "GeomField0001",
128 FatalException, "Method no longer used.");
129
130} // end of StepWithEst ......................................................
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35

References FatalException, and G4Exception().

Field Documentation

◆ dydxm

G4double* G4ClassicalRK4::dydxm
private

Definition at line 88 of file G4ClassicalRK4.hh.

Referenced by DumbStepper(), G4ClassicalRK4(), and ~G4ClassicalRK4().

◆ dydxMid

G4double * G4MagErrorStepper::dydxMid
privateinherited

◆ dydxt

G4double * G4ClassicalRK4::dydxt
private

Definition at line 88 of file G4ClassicalRK4.hh.

Referenced by DumbStepper(), G4ClassicalRK4(), and ~G4ClassicalRK4().

◆ fEquation_Rhs

G4EquationOfMotion* G4MagIntegratorStepper::fEquation_Rhs = nullptr
privateinherited

Definition at line 124 of file G4MagIntegratorStepper.hh.

◆ fFinalPoint

G4ThreeVector G4MagErrorStepper::fFinalPoint
privateinherited

◆ fInitialPoint

G4ThreeVector G4MagErrorStepper::fInitialPoint
privateinherited

◆ fIntegrationOrder

G4int G4MagIntegratorStepper::fIntegrationOrder = -1
privateinherited

Definition at line 134 of file G4MagIntegratorStepper.hh.

◆ fIsFSAL

G4bool G4MagIntegratorStepper::fIsFSAL = false
privateinherited

Definition at line 136 of file G4MagIntegratorStepper.hh.

◆ fMidPoint

G4ThreeVector G4MagErrorStepper::fMidPoint
privateinherited

◆ fNoIntegrationVariables

const G4int G4MagIntegratorStepper::fNoIntegrationVariables = 0
privateinherited

Definition at line 125 of file G4MagIntegratorStepper.hh.

◆ fNoRHSCalls

unsigned long G4MagIntegratorStepper::fNoRHSCalls = 0UL
mutableprivateinherited

Definition at line 128 of file G4MagIntegratorStepper.hh.

◆ fNoStateVariables

const G4int G4MagIntegratorStepper::fNoStateVariables = 0
privateinherited

Definition at line 126 of file G4MagIntegratorStepper.hh.

◆ yInitial

G4double* G4MagErrorStepper::yInitial
privateinherited

◆ yMiddle

G4double * G4MagErrorStepper::yMiddle
privateinherited

◆ yOneStep

G4double * G4MagErrorStepper::yOneStep
privateinherited

◆ yt

G4double * G4ClassicalRK4::yt
private

Definition at line 88 of file G4ClassicalRK4.hh.

Referenced by DumbStepper(), G4ClassicalRK4(), and ~G4ClassicalRK4().


The documentation for this class was generated from the following files: