Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4MagErrorStepper.cc
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//
27// $Id$
28//
29// --------------------------------------------------------------------
30
31#include "G4MagErrorStepper.hh"
32#include "G4LineSection.hh"
33
35{
36 delete[] yMiddle;
37 delete[] dydxMid;
38 delete[] yInitial;
39 delete[] yOneStep;
40}
41
42void
44 const G4double dydx[],
45 G4double hstep,
46 G4double yOutput[],
47 G4double yError [] )
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}
92
93
94
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}
120
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
static G4double Distline(const G4ThreeVector &OtherPnt, const G4ThreeVector &LinePntA, const G4ThreeVector &LinePntB)
virtual void DumbStepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[])=0
virtual ~G4MagErrorStepper()
G4double DistChord() const
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
G4int GetNumberOfVariables() const
void RightHandSide(const double y[], double dydx[])
G4int GetNumberOfStateVariables() const
virtual G4int IntegratorOrder() const =0