Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4MagHelicalStepper.hh
이 파일의 문서화 페이지로 가기
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 93806 2015-11-02 11:21:01Z 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 // It is used for a set of steppers which use the helix as a sort of
40 // 'first order' solution.
41 // - Most obtain an error by breaking up the step in two
42 // - G4ExactHelicalStepper does not provide an error estimate
43 
44 // History:
45 // - 05.11.98 J.Apostolakis Creation of new ABC
46 // --------------------------------------------------------------------
47 
48 #ifndef G4MagHelicalStepper_hh
49 #define G4MagHelicalStepper_hh
50 
52 
53 #include "G4Types.hh"
55 #include "G4Mag_EqRhs.hh"
56 #include "G4ThreeVector.hh"
57 
59 {
60  public: // with description
61 
63  virtual ~G4MagHelicalStepper();
64 
65  virtual void Stepper( const G4double y[], // VIRTUAL for ExactHelix - temporary
66  const G4double dydx[],
67  G4double h,
68  G4double yout[],
69  G4double yerr[] );
70  // The stepper for the Runge Kutta integration.
71  // The stepsize is fixed, equal to h.
72  // Integrates ODE starting values y[0 to 6]
73  // Outputs yout[] and its estimated error yerr[].
74 
75  virtual void DumbStepper( const G4double y[],
76  G4ThreeVector Bfld,
77  G4double h,
78  G4double yout[] ) = 0;
79  // Performs a 'dump' Step without error calculation.
80 
81  G4double DistChord()const ;
82  // Estimate maximum distance of curved solution and chord ...
83 
84  protected: // with description
85 
86  inline void LinearStep( const G4double yIn[],
87  G4double h,
88  G4double yHelix[]) const;
89  // A linear Step in regions without magnetic field.
90 
91  void AdvanceHelix( const G4double yIn[],
92  G4ThreeVector Bfld,
93  G4double h,
94  G4double yHelix[],G4double yHelix2[]=0); // output
95  // A first order Step along a helix inside the field.
96 
97  inline void MagFieldEvaluate( const G4double y[], G4ThreeVector& Bfield );
98  // Evaluate the field at a certain point.
99 
100 
101  inline G4double GetInverseCurve( const G4double Momentum, const G4double Bmag );
102  // Evaluate Inverse of Curvature of Track
103 
104  // Store and use the parameters of track :
105  // Radius of curve, Stepping angle, Radius of projected helix
106  inline void SetAngCurve(const G4double Ang);
107  inline G4double GetAngCurve()const;
108 
109  inline void SetCurve(const G4double Curve);
110  inline G4double GetCurve()const;
111 
112  inline void SetRadHelix(const G4double Rad);
113  inline G4double GetRadHelix()const;
114 
115 
116  protected: // without description
117 
118  // void MagFieldEvaluate( const G4double y[], G4double B[] )
119  // { GetEquationOfMotion()-> GetFieldValue(y, B); }
120 
121  private:
122 
125  // Private copy constructor and assignment operator.
126 
127  static const G4double fUnitConstant; // As in G4Mag_EqRhs.hh/cc where it is not used.
128  private:
129 
131 
132  // Data stored in order to find the chord.
136  // Data stored in order to find the chord.
138 
139 
140 };
141 
142 #include "G4MagHelicalStepper.icc"
143 
144 #endif /* G4MagHelicalStepper_hh */
void SetCurve(const G4double Curve)
G4MagHelicalStepper(G4Mag_EqRhs *EqRhs)
G4double GetAngCurve() const
Float_t y
Definition: compare.C:6
void SetAngCurve(const G4double Ang)
void AdvanceHelix(const G4double yIn[], G4ThreeVector Bfld, G4double h, G4double yHelix[], G4double yHelix2[]=0)
double G4double
Definition: G4Types.hh:76
G4double GetRadHelix() const
G4double GetInverseCurve(const G4double Momentum, const G4double Bmag)
G4MagHelicalStepper & operator=(const G4MagHelicalStepper &)
void SetRadHelix(const G4double Rad)
static const G4double fUnitConstant
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
G4double GetCurve() const
void MagFieldEvaluate(const G4double y[], G4ThreeVector &Bfield)
void LinearStep(const G4double yIn[], G4double h, G4double yHelix[]) const
G4double DistChord() const