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
G4HelixMixedStepper.hh
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// class G4HelixMixedStepper
27//
28// Class description:
29//
30// G4HelixMixedStepper split the Method used for Integration in two:
31//
32// If Stepping Angle ( h / R_curve) < pi/3 : use Stepper for small step
33//
34// Else use HelixExplicitEuler Stepper
35//
36// Stepper for the small step is G4ClassicalRK4 by default, but
37// it possible to choose other stepper,like G4CashKarpRK45 or G4RKG3_Stepper,
38// by setting StepperNumber : new HelixMixedStepper(EqRhs,N)
39//
40// N=0 G4ExplicitEuler; N=1 G4ImplicitEuler;
41// N=2 G4SimpleRunge; N=3 G4SimpleHeum;
42// N=4 G4ClassicalRK4; N=5 G4HelixExplicitEuler;
43// N=6 G4HelixExplicitEuler; N=7 G4HelixSimpleRunge;
44// N=8 G4CashKarpRK45; N=9 G4ExactHelixStepper;
45// N=10 G4RKG3_Stepper;
46//
47// History:
48// Derived from ExactHelicalStepper 18/05/07
49//
50// -------------------------------------------------------------------
51
52#ifndef G4HELIXMIXEDSTEPPER_HH
53#define G4HELIXMIXEDSTEPPER_HH
54
56
57
59{
60
61 public:
62
63 G4HelixMixedStepper(G4Mag_EqRhs *EqRhs,G4int fStepperNumber=0);
65
66 void Stepper( const G4double y[],
67 const G4double dydx[],
68 G4double h,
69 G4double yout[],
70 G4double yerr[] );
71 // Step 'integration' for step size 'h'
72 // If SteppingAngle=h/R_curve<pi/3 uses RK4Stepper
73 // Else Helix Fast Method
74
75
76 void DumbStepper( const G4double y[],
77 G4ThreeVector Bfld,
78 G4double h,
79 G4double yout[]);
80 G4double DistChord() const;
81 // Estimate maximum distance of curved solution and chord ...
82
83
84 public: // with description
85
86 inline void SetVerbose (G4int newvalue){fVerbose=newvalue;}
87
88 public: // without description
89 void PrintCalls();
91
92
93 G4int IntegratorOrder() const { return 4; }
94 private:
95 // Mixed Integration RK4 for 'small' steps
96 G4MagIntegratorStepper* fRK4Stepper;
97
98 private:
99 // Used for statistic = how many calls to different steppers
100 G4int fVerbose;
101 G4int fNumCallsRK4;
102 G4int fNumCallsHelix;
103
104};
105
106#endif /* G4HELIXMIXEDSTEPPER_HH */
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
void SetVerbose(G4int newvalue)
G4double DistChord() const
void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])
G4int IntegratorOrder() const
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)