Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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// 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=2 G4SimpleRunge; N=3 G4SimpleHeum;
41// N=4 G4ClassicalRK4;
42// N=6 G4HelixImplicitEuler; N=7 G4HelixSimpleRunge;
43// N=8 G4CashKarpRK45; N=9 G4ExactHelixStepper;
44// N=10 G4RKG3_Stepper; N=13 G4NystromRK4
45// N=23 BogackiShampine23 N=145 TsitourasRK45
46// N=45 BogackiShampine45 N=745 DormandPrince745 (ie DoPri5)
47//
48// For completeness also available are:
49// N=11 G4ExplicitEuler N=12 G4ImplicitEuler; -- Likely poor
50// N=5 G4HelixExplicitEuler (testing only)
51// For recommendations see comments in 'SetupStepper' method.
52//
53// Note: Like other helix steppers, only applicable in pure magnetic field
54
55// Created: T.Nikitina, CERN - 18.05.2007, derived from G4ExactHelicalStepper
56// -------------------------------------------------------------------
57#ifndef G4HELIXMIXEDSTEPPER_HH
58#define G4HELIXMIXEDSTEPPER_HH
59
61
63{
64 public:
65
67 G4int StepperNumber = -1,
68 G4double Angle_threshold = -1.0);
69 ~G4HelixMixedStepper() override;
70
71 void Stepper( const G4double y[],
72 const G4double dydx[],
73 G4double h,
74 G4double yout[],
75 G4double yerr[] ) override;
76 // Step 'integration' for step size 'h'
77 // If SteppingAngle= h/R_curve < pi/3 uses default RK stepper
78 // else use Helix Fast Method
79
80 void DumbStepper( const G4double y[],
81 G4ThreeVector Bfld,
82 G4double h,
83 G4double yout[]) override;
84
85 G4double DistChord() const override;
86 // Estimate maximum distance of curved solution and chord ...
87
88 inline void SetVerbose (G4int newvalue) { fVerbose = newvalue; }
89
90 void PrintCalls();
92
93 inline void SetAngleThreshold( G4double val ) { fAngle_threshold = val; }
94 inline G4double GetAngleThreshold() { return fAngle_threshold; }
95 inline G4int IntegratorOrder() const override { return 4; }
96
97 private:
98
99 G4MagIntegratorStepper* fRK4Stepper = nullptr;
100 // Mixed Integration RK4 for 'small' steps
101 G4int fStepperNumber = -1;
102 // Int ID of RK stepper
103 G4double fAngle_threshold = -1.0;
104 // Threshold angle (in radians ) - above it Helical stepper is used
105
106 private:
107
108 G4int fVerbose = 0;
109
110 G4int fNumCallsRK4 = 0;
111 G4int fNumCallsHelix = 0;
112 // Used for statistic = how many calls to different steppers
113};
114
115#endif
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
G4HelixMixedStepper(G4Mag_EqRhs *EqRhs, G4int StepperNumber=-1, G4double Angle_threshold=-1.0)
void SetAngleThreshold(G4double val)
void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[]) override
void SetVerbose(G4int newvalue)
G4int IntegratorOrder() const override
void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[]) override
G4double DistChord() const override
G4MagIntegratorStepper * SetupStepper(G4Mag_EqRhs *EqRhs, G4int StepperName)