Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4MagErrorStepper Class Referenceabstract

#include <G4MagErrorStepper.hh>

Inheritance diagram for G4MagErrorStepper:
G4MagIntegratorStepper G4ClassicalRK4 G4ConstRK4 G4ExplicitEuler G4ImplicitEuler G4SimpleHeum G4SimpleRunge

Public Member Functions

 G4MagErrorStepper (G4EquationOfMotion *EqRhs, G4int numberOfVariables, G4int numStateVariables=12)
 
virtual ~G4MagErrorStepper ()
 
void Stepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
 
virtual void DumbStepper (const G4double y[], const G4double dydx[], G4double h, G4double yout[])=0
 
G4double DistChord () const
 
- Public Member Functions inherited from G4MagIntegratorStepper
 G4MagIntegratorStepper (G4EquationOfMotion *Equation, G4int numIntegrationVariables, G4int numStateVariables=12)
 
virtual ~G4MagIntegratorStepper ()
 
virtual void ComputeRightHandSide (const G4double y[], G4double dydx[])
 
void NormaliseTangentVector (G4double vec[6])
 
void NormalisePolarizationVector (G4double vec[12])
 
void RightHandSide (const double y[], double dydx[])
 
G4int GetNumberOfVariables () const
 
G4int GetNumberOfStateVariables () const
 
virtual G4int IntegratorOrder () const =0
 
G4EquationOfMotionGetEquationOfMotion ()
 
void SetEquationOfMotion (G4EquationOfMotion *newEquation)
 

Detailed Description

Definition at line 49 of file G4MagErrorStepper.hh.

Constructor & Destructor Documentation

G4MagErrorStepper::G4MagErrorStepper ( G4EquationOfMotion EqRhs,
G4int  numberOfVariables,
G4int  numStateVariables = 12 
)
G4MagErrorStepper::~G4MagErrorStepper ( )
virtual

Definition at line 34 of file G4MagErrorStepper.cc.

35 {
36  delete[] yMiddle;
37  delete[] dydxMid;
38  delete[] yInitial;
39  delete[] yOneStep;
40 }

Member Function Documentation

G4double G4MagErrorStepper::DistChord ( ) const
virtual

Implements G4MagIntegratorStepper.

Definition at line 96 of file G4MagErrorStepper.cc.

References G4LineSection::Distline().

97 {
98  // Estimate the maximum distance from the curve to the chord
99  //
100  // We estimate this using the distance of the midpoint to
101  // chord (the line between
102  //
103  // Method below is good only for angle deviations < 2 pi,
104  // This restriction should not a problem for the Runge cutta methods,
105  // which generally cannot integrate accurately for large angle deviations.
106  G4double distLine, distChord;
107 
108  if (fInitialPoint != fFinalPoint) {
109  distLine= G4LineSection::Distline( fMidPoint, fInitialPoint, fFinalPoint );
110  // This is a class method that gives distance of Mid
111  // from the Chord between the Initial and Final points.
112 
113  distChord = distLine;
114  }else{
115  distChord = (fMidPoint-fInitialPoint).mag();
116  }
117 
118  return distChord;
119 }
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
double G4double
Definition: G4Types.hh:76
virtual void G4MagErrorStepper::DumbStepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[] 
)
pure virtual
void G4MagErrorStepper::Stepper ( const G4double  y[],
const G4double  dydx[],
G4double  h,
G4double  yout[],
G4double  yerr[] 
)
virtual

Implements G4MagIntegratorStepper.

Definition at line 43 of file G4MagErrorStepper.cc.

References DumbStepper(), G4MagIntegratorStepper::GetNumberOfStateVariables(), G4MagIntegratorStepper::GetNumberOfVariables(), G4MagIntegratorStepper::IntegratorOrder(), and G4MagIntegratorStepper::RightHandSide().

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

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