source: trunk/source/geometry/magneticfield/include/G4CashKarpRKF45.hh@ 1315

Last change on this file since 1315 was 1231, checked in by garnier, 16 years ago

update from CVS

File size: 3.9 KB
Line 
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: G4CashKarpRKF45.hh,v 1.11 2008/01/11 15:23:54 japost Exp $
28// GEANT4 tag $Name: $
29//
30//
31// class G4CashKarpRKF45
32//
33// Class description:
34//
35// The Cash-Karp Runge-Kutta-Fehlberg 4/5 method is an embedded fourth
36// order method (giving fifth-order accuracy) for the solution of an ODE.
37// Two different fourth order estimates are calculated; their difference
38// gives an error estimate. [ref. Numerical Recipes in C, 2nd Edition]
39// It is used to integrate the equations of the motion of a particle
40// in a magnetic field.
41
42// History:
43// - Created. J.Apostolakis, V.Grichine - 30.1.97
44// -------------------------------------------------------------------
45
46#ifndef G4CashKARP_RKF45
47#define G4CashKARP_RKF45
48
49#include "G4MagIntegratorStepper.hh"
50
51class G4CashKarpRKF45 : public G4MagIntegratorStepper
52{
53
54 public: // with description
55
56 G4CashKarpRKF45( G4EquationOfMotion *EqRhs,
57 G4int numberOfVariables = 6,
58 G4bool primary= true ) ;
59 ~G4CashKarpRKF45() ;
60
61 void Stepper( const G4double y[],
62 const G4double dydx[],
63 G4double h,
64 G4double yout[],
65 G4double yerr[] ) ;
66
67 public: // without description
68
69 G4double DistChord() const;
70 G4int IntegratorOrder() const { return 4; }
71
72 private:
73
74 void StepWithEst( const G4double yIn[],
75 const G4double dydx[],
76 G4double Step,
77 G4double yOut[],
78 G4double& alpha2,
79 G4double& beta2,
80 const G4double B1[],
81 G4double B2[] );
82 // No longer used. Obsolete.
83
84 G4CashKarpRKF45(const G4CashKarpRKF45&);
85 G4CashKarpRKF45& operator=(const G4CashKarpRKF45&);
86 // Private copy constructor and assignment operator.
87
88 private:
89
90 // G4int fNumberOfVariables ; // Already kept in G4MagIntegratorStepper
91 G4double *ak2, *ak3, *ak4, *ak5, *ak6, *ak7, *yTemp, *yIn; // scratch space
92
93 // for DistChord calculations
94
95 G4double fLastStepLength;
96 G4double *fLastInitialVector, *fLastFinalVector,
97 *fLastDyDx, *fMidVector, *fMidError;
98 G4CashKarpRKF45* fAuxStepper;
99};
100
101#endif /* G4CashKARP_RKF45 */
Note: See TracBrowser for help on using the repository browser.