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

#include <G4RK547FEq3.hh>

Inheritance diagram for G4RK547FEq3:
G4MagIntegratorStepper

Public Member Functions

virtual G4double DistChord () const override
 
 G4RK547FEq3 (const G4RK547FEq3 &)=delete
 
 G4RK547FEq3 (G4EquationOfMotion *EqRhs, G4int integrationVariables=6)
 
G4EquationOfMotionGetEquationOfMotion ()
 
const G4EquationOfMotionGetEquationOfMotion () const
 
unsigned long GetfNoRHSCalls ()
 
G4int GetNumberOfStateVariables () const
 
G4int GetNumberOfVariables () const
 
G4int IntegrationOrder ()
 
virtual G4int IntegratorOrder () const override
 
G4bool IsFSAL () const
 
void NormalisePolarizationVector (G4double vec[12])
 
void NormaliseTangentVector (G4double vec[6])
 
G4RK547FEq3operator= (const G4RK547FEq3 &)=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)
 
virtual void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[]) override
 
void Stepper (const G4double yInput[], const G4double dydx[], G4double hstep, G4double yOutput[], G4double yError[], G4double dydxOutput[])
 

Protected Member Functions

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

Private Member Functions

void makeStep (const G4double yInput[], const G4double dydx[], const G4double hstep, G4double yOutput[], G4double *dydxOutput=nullptr, G4double *yError=nullptr) const
 

Private Attributes

G4double fdydx [G4FieldTrack::ncompSVEC]
 
G4double fdydxOut [G4FieldTrack::ncompSVEC]
 
G4EquationOfMotionfEquation_Rhs = nullptr
 
G4double fhstep = -1.0
 
G4int fIntegrationOrder = -1
 
G4bool fIsFSAL = false
 
const G4int fNoIntegrationVariables = 0
 
unsigned long fNoRHSCalls = 0UL
 
const G4int fNoStateVariables = 0
 
G4double fyIn [G4FieldTrack::ncompSVEC]
 
G4double fyOut [G4FieldTrack::ncompSVEC]
 

Detailed Description

Definition at line 44 of file G4RK547FEq3.hh.

Constructor & Destructor Documentation

◆ G4RK547FEq3() [1/2]

G4RK547FEq3::G4RK547FEq3 ( G4EquationOfMotion EqRhs,
G4int  integrationVariables = 6 
)

Definition at line 51 of file G4RK547FEq3.cc.

52 : G4MagIntegratorStepper(EqRhs, integrationVariables)
53{
54}
G4MagIntegratorStepper(G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12, G4bool isFSAL=false)

◆ G4RK547FEq3() [2/2]

G4RK547FEq3::G4RK547FEq3 ( const G4RK547FEq3 )
delete

Member Function Documentation

◆ DistChord()

G4double G4RK547FEq3::DistChord ( ) const
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 167 of file G4RK547FEq3.cc.

168{
170 makeStep(fyIn, fdydx, fhstep / 2., yMid);
171
172 const G4ThreeVector begin = makeVector(fyIn, Value3D::Position);
173 const G4ThreeVector mid = makeVector(yMid, Value3D::Position);
174 const G4ThreeVector end = makeVector(fyOut, Value3D::Position);
175
176 return G4LineSection::Distline(mid, begin, end);
177}
double G4double
Definition: G4Types.hh:83
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
void makeStep(const G4double yInput[], const G4double dydx[], const G4double hstep, G4double yOutput[], G4double *dydxOutput=nullptr, G4double *yError=nullptr) const
Definition: G4RK547FEq3.cc:56
G4double fhstep
Definition: G4RK547FEq3.hh:84
G4double fdydx[G4FieldTrack::ncompSVEC]
Definition: G4RK547FEq3.hh:80
G4double fyIn[G4FieldTrack::ncompSVEC]
Definition: G4RK547FEq3.hh:79
G4double fyOut[G4FieldTrack::ncompSVEC]
Definition: G4RK547FEq3.hh:81
G4ThreeVector makeVector(const ArrayType &array, Value3D value)

References G4LineSection::Distline(), fdydx, fhstep, fyIn, fyOut, makeStep(), field_utils::makeVector(), and G4FieldTrack::ncompSVEC.

◆ 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()

virtual G4int G4RK547FEq3::IntegratorOrder ( ) const
inlineoverridevirtual

Implements G4MagIntegratorStepper.

Definition at line 67 of file G4RK547FEq3.hh.

67{ return 4; }

◆ IsFSAL()

G4bool G4MagIntegratorStepper::IsFSAL ( ) const
inlineinherited

◆ makeStep()

void G4RK547FEq3::makeStep ( const G4double  yInput[],
const G4double  dydx[],
const G4double  hstep,
G4double  yOutput[],
G4double dydxOutput = nullptr,
G4double yError = nullptr 
) const
private

Definition at line 56 of file G4RK547FEq3.cc.

62{
65 {
66 yOutput[i] = yTemp[i] = yInput[i];
67 }
68
74
75 const G4double b21 = 11./45.,
76 b31 = 11./120., b32 = 11./40.,
77 b41 = 106865./87808., b42 = -408375./87808.,
78 b43 = 193875./43904.,
79 b51 = 79503./121000., b52 = -1053./440.,
80 b53 = 147753./56870., b54 = 27048./710875.,
81 b61 = 89303./78045., b62 = -2025./473.,
82 b63 = 994650./244541., b64 = -2547216./28122215.,
83 b65 = 475./2967.,
84 b71 = 1247./10890., b72 = 0., b73 = 57375./108053.,
85 b74 = -1229312./1962015., b75 = 125./207.,
86 b76 = 43./114.;
87
88 const G4double dc1 = b71 - 21487./185130.,
89 dc2 = b72 - 0.,
90 dc3 = b73 - 963225./1836901.,
91 dc4 = b74 + 39864832./33354255.,
92 dc5 = b75 - 2575./3519.,
93 dc6 = b76 - 4472./4845.,
94 dc7 = 0. + 1./10.;
95
96 // RightHandSide(yInput, dydx);
97 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
98 yTemp[i] = yInput[i] + hstep * b21 * dydx[i];
99
100 RightHandSide(yTemp, ak2);
101 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
102 yTemp[i] = yInput[i] + hstep * (b31 * dydx[i] + b32 * ak2[i]);
103
104 RightHandSide(yTemp, ak3);
105 for(G4int i = 0;i < GetNumberOfVariables(); ++i)
106 yTemp[i] = yInput[i] + hstep * (b41 * dydx[i] + b42 * ak2[i] +
107 b43 * ak3[i]);
108
109 RightHandSide(yTemp, ak4);
110 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
111 yTemp[i] = yInput[i] + hstep * (b51 * dydx[i] + b52 * ak2[i] +
112 b53 * ak3[i] + b54 * ak4[i]);
113
114 RightHandSide(yTemp, ak5);
115 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
116 yTemp[i] = yInput[i] + hstep * (b61 * dydx[i] + b62 * ak2[i] +
117 b63 * ak3[i] + b64 * ak4[i] +
118 b65 * ak5[i]);
119
120 RightHandSide(yTemp, ak6);
121 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
122 yOutput[i] = yInput[i] + hstep * (b71 * dydx[i] + b72 * ak2[i] +
123 b73 * ak3[i] + b74 * ak4[i] +
124 b75 * ak5[i] + b76 * ak6[i]);
125 if (dydxOutput && yError)
126 {
127 RightHandSide(yOutput, dydxOutput);
128 for(G4int i = 0; i < GetNumberOfVariables(); ++i)
129 yError[i] = hstep * (dc1 * dydx[i] + dc2 * ak2[i] + dc3 * ak3[i] +
130 dc4 * ak4[i] + dc5 * ak5[i] + dc6 * ak6[i] +
131 dc7 * dydxOutput[i]);
132 }
133}
static const G4double ak2
int G4int
Definition: G4Types.hh:85
G4int GetNumberOfVariables() const
void RightHandSide(const G4double y[], G4double dydx[]) const
G4int GetNumberOfStateVariables() const

References ak2, G4MagIntegratorStepper::GetNumberOfStateVariables(), G4MagIntegratorStepper::GetNumberOfVariables(), G4FieldTrack::ncompSVEC, and G4MagIntegratorStepper::RightHandSide().

Referenced by DistChord(), and Stepper().

◆ NormalisePolarizationVector()

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

◆ NormaliseTangentVector()

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

◆ operator=()

G4RK547FEq3 & G4RK547FEq3::operator= ( const G4RK547FEq3 )
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() [1/2]

void G4RK547FEq3::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[] 
)
overridevirtual

Implements G4MagIntegratorStepper.

Definition at line 135 of file G4RK547FEq3.cc.

140{
141 copy(fyIn, yInput);
142 copy(fdydx, dydx);
143 fhstep = hstep;
144
145 makeStep(fyIn, fdydx, fhstep, fyOut, fdydxOut, yError);
146
147 copy(yOutput, fyOut);
148}
G4double fdydxOut[G4FieldTrack::ncompSVEC]
Definition: G4RK547FEq3.hh:82
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

References field_utils::copy(), fdydx, fdydxOut, fhstep, fyIn, fyOut, and makeStep().

◆ Stepper() [2/2]

void G4RK547FEq3::Stepper ( const G4double  yInput[],
const G4double  dydx[],
G4double  hstep,
G4double  yOutput[],
G4double  yError[],
G4double  dydxOutput[] 
)

Definition at line 150 of file G4RK547FEq3.cc.

156{
157 copy(fyIn, yInput);
158 copy(fdydx, dydx);
159 fhstep = hstep;
160
161 makeStep(fyIn, fdydx, fhstep, fyOut, fdydxOut, yError);
162
163 copy(yOutput, fyOut);
164 copy(dydxOutput, fdydxOut);
165}

References field_utils::copy(), fdydx, fdydxOut, fhstep, fyIn, fyOut, and makeStep().

Field Documentation

◆ fdydx

G4double G4RK547FEq3::fdydx[G4FieldTrack::ncompSVEC]
private

Definition at line 80 of file G4RK547FEq3.hh.

Referenced by DistChord(), and Stepper().

◆ fdydxOut

G4double G4RK547FEq3::fdydxOut[G4FieldTrack::ncompSVEC]
private

Definition at line 82 of file G4RK547FEq3.hh.

Referenced by Stepper().

◆ fEquation_Rhs

G4EquationOfMotion* G4MagIntegratorStepper::fEquation_Rhs = nullptr
privateinherited

Definition at line 124 of file G4MagIntegratorStepper.hh.

◆ fhstep

G4double G4RK547FEq3::fhstep = -1.0
private

Definition at line 84 of file G4RK547FEq3.hh.

Referenced by DistChord(), and Stepper().

◆ 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.

◆ 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.

◆ fyIn

G4double G4RK547FEq3::fyIn[G4FieldTrack::ncompSVEC]
private

Definition at line 79 of file G4RK547FEq3.hh.

Referenced by DistChord(), and Stepper().

◆ fyOut

G4double G4RK547FEq3::fyOut[G4FieldTrack::ncompSVEC]
private

Definition at line 81 of file G4RK547FEq3.hh.

Referenced by DistChord(), and Stepper().


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