Geant4  10.00.p03
G4MagHelicalStepper.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 //
27 //
28 // $Id: G4MagHelicalStepper.hh 66356 2012-12-18 09:02:32Z gcosmo $
29 //
30 //
31 //
32 // class G4MagHelicalStepper
33 //
34 // Class description:
35 //
36 // Abstract base class for integrator of particle's equation of motion,
37 // used in tracking in space dependent magnetic field
38 
39 // History:
40 // - 05.11.98 J.Apostolakis Creation of new ABC
41 // --------------------------------------------------------------------
42 
43 #ifndef G4MagHelicalStepper_hh
44 #define G4MagHelicalStepper_hh
45 
46 #include <CLHEP/Units/PhysicalConstants.h>
47 
48 #include "G4Types.hh"
50 #include "G4Mag_EqRhs.hh"
51 #include "G4ThreeVector.hh"
52 
54 {
55  public: // with description
56 
58  virtual ~G4MagHelicalStepper();
59 
60  virtual void Stepper( const G4double y[], // VIRTUAL for ExactHelix - temporary
61  const G4double dydx[],
62  G4double h,
63  G4double yout[],
64  G4double yerr[] );
65  // The stepper for the Runge Kutta integration.
66  // The stepsize is fixed, equal to h.
67  // Integrates ODE starting values y[0 to 6]
68  // Outputs yout[] and its estimated error yerr[].
69 
70  virtual void DumbStepper( const G4double y[],
71  G4ThreeVector Bfld,
72  G4double h,
73  G4double yout[] ) = 0;
74  // Performs a 'dump' Step without error calculation.
75 
76  G4double DistChord()const ;
77  // Estimate maximum distance of curved solution and chord ...
78 
79  protected: // with description
80 
81  inline void LinearStep( const G4double yIn[],
82  G4double h,
83  G4double yHelix[]) const;
84  // A linear Step in regions without magnetic field.
85 
86  void AdvanceHelix( const G4double yIn[],
87  G4ThreeVector Bfld,
88  G4double h,
89  G4double yHelix[],G4double yHelix2[]=0); // output
90  // A first order Step along a helix inside the field.
91 
92  inline void MagFieldEvaluate( const G4double y[], G4ThreeVector& Bfield );
93  // Evaluate the field at a certain point.
94 
95 
96  inline G4double GetInverseCurve( const G4double Momentum, const G4double Bmag );
97  // Evaluate Inverse of Curvature of Track
98 
99  // Store and use the parameters of track :
100  // Radius of curve, Stepping angle, Radius of projected helix
101  inline void SetAngCurve(const G4double Ang);
102  inline G4double GetAngCurve()const;
103 
104  inline void SetCurve(const G4double Curve);
105  inline G4double GetCurve()const;
106 
107  inline void SetRadHelix(const G4double Rad);
108  inline G4double GetRadHelix()const;
109 
110 
111  protected: // without description
112 
113  // void MagFieldEvaluate( const G4double y[], G4double B[] )
114  // { GetEquationOfMotion()-> GetFieldValue(y, B); }
115 
116  private:
117 
120  // Private copy constructor and assignment operator.
121 
122  static const G4double fUnitConstant; // As in G4Mag_EqRhs.hh/cc where it is not used.
123  private:
124 
126 
127  // Data stored in order to find the chord.
131  // Data stored in order to find the chord.
133 
134 
135 };
136 
137 #include "G4MagHelicalStepper.icc"
138 
139 #endif /* G4MagHelicalStepper_hh */
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
virtual void Stepper(const G4double y[], const G4double dydx[], G4double h, G4double yout[], G4double yerr[])
virtual void DumbStepper(const G4double y[], G4ThreeVector Bfld, G4double h, G4double yout[])=0
CLHEP::Hep3Vector G4ThreeVector
G4double GetRadHelix() const
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
G4double GetCurve() const
void LinearStep(const G4double yIn[], G4double h, G4double yHelix[]) const
G4MagHelicalStepper(G4Mag_EqRhs *EqRhs)
G4double GetInverseCurve(const G4double Momentum, const G4double Bmag)
void SetRadHelix(const G4double Rad)
void SetAngCurve(const G4double Ang)
G4MagHelicalStepper & operator=(const G4MagHelicalStepper &)
void SetCurve(const G4double Curve)
double G4double
Definition: G4Types.hh:76
G4double GetAngCurve() const
G4double DistChord() const
static const G4double fUnitConstant