Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FieldTrack.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 // $Id: G4FieldTrack.hh 87868 2015-01-16 08:22:22Z gcosmo $
28 //
29 //
30 // class G4FieldTrack
31 //
32 // Class description:
33 //
34 // Data structure bringing together a magnetic track's state.
35 // (position, momentum direction & modulus, energy, spin, ... )
36 // Uses/abilities:
37 // - does not maintain any relationship between its data (eg energy/momentum).
38 // - for use in Runge-Kutta solver (in passing it the values right now).
39 
40 // History
41 // - First version: Oct 14, 1996 John Apostolakis
42 // - Modified: Oct 24, 1996 JA: Added dist_on_curve, deleted constructor
43 // Nov 5, 1998 JA: Added energy, momentum, TOF, spin &
44 // several constructor, access, set methods
45 // May 10, 2006 JA: Added charge, "default" constructor
46 // -------------------------------------------------------------------
47 
48 #ifndef G4FieldTrack_HH
49 #define G4FieldTrack_HH
50 
51 #include "G4ThreeVector.hh"
52 #include "G4ChargeState.hh"
53 
55 {
56  public: // with description
57 
58  G4FieldTrack( const G4ThreeVector& pPosition,
59  G4double LaboratoryTimeOfFlight,
60  const G4ThreeVector& pMomentumDirection,
62  G4double restMass_c2,
63  G4double charge,
64  const G4ThreeVector& polarization,
65  G4double magnetic_dipole_moment= 0.0,
66  G4double curve_length= 0.0,
67  G4double PDGspin = -1.0 );
68 
69  G4FieldTrack( const G4FieldTrack& pFieldTrack );
70  G4FieldTrack( char ); // Almost default constructor
71 
72  ~G4FieldTrack();
73  // End of preferred Constructors / Destructor
74 
75  inline void
76  UpdateState( const G4ThreeVector& pPosition,
77  G4double LaboratoryTimeOfFlight,
78  const G4ThreeVector& pMomentumDirection,
79  G4double kineticEnergy);
80  // Update four-vectors for space/time and momentum/energy
81  // Also resets curve length.
82  inline
83  void UpdateFourMomentum( G4double kineticEnergy,
84  const G4ThreeVector& momentumDirection );
85  // Update momentum (and direction), and kinetic energy
86 
87  void SetChargeAndMoments(G4double charge,
88  G4double magnetic_dipole_moment= DBL_MAX,
89  G4double electric_dipole_moment= DBL_MAX,
90  G4double magnetic_charge=DBL_MAX );
91  // Sets the charges and moments that are not given as DBL_MAX
92 
93  inline void SetPDGSpin(G4double pdgSpin){ fChargeState.SetPDGSpin(pdgSpin); }
94  inline G4double GetPDGSpin(){ return fChargeState.GetPDGSpin(); }
95 
96  G4FieldTrack( const G4ThreeVector& pPosition,
97  const G4ThreeVector& pMomentumDirection,
98  G4double curve_length,
99  G4double kineticEnergy,
100  const G4double restMass_c2,
101  G4double velocity,
102  G4double LaboratoryTimeOfFlight=0.0,
103  G4double ProperTimeOfFlight=0.0,
104  const G4ThreeVector* pPolarization=0,
105  G4double PDGspin = -1.0 );
106  // Older constructor
107  // ---> Misses charge !!!
108 
109  inline G4FieldTrack& operator = ( const G4FieldTrack & rStVec );
110  // Assignment operator
111 
112  inline G4ThreeVector GetMomentum() const;
113  inline G4ThreeVector GetPosition() const;
114  inline const G4ThreeVector& GetMomentumDir() const;
115  inline G4ThreeVector GetMomentumDirection() const;
116  inline G4double GetCurveLength() const;
117  // Distance along curve of point.
118 
119  inline G4ThreeVector GetPolarization() const;
120  inline void SetPolarization( const G4ThreeVector& vecPol );
121 
122  inline G4double GetLabTimeOfFlight() const;
123  inline G4double GetProperTimeOfFlight() const;
124  inline G4double GetKineticEnergy() const;
125  inline G4double GetCharge() const;
126  inline G4double GetRestMass() const { return fRestMass_c2; }
127  // Accessors.
128 
129  inline void SetPosition(G4ThreeVector nPos);
130  inline void SetMomentum(G4ThreeVector nMomDir);
131  // Does change mom-dir too.
132 
133  inline void SetMomentumDir(G4ThreeVector nMomDir);
134  // Does NOT change Momentum or Velocity Vector.
135 
136  inline void SetRestMass(G4double Mass_c2) { fRestMass_c2= Mass_c2; }
137 
138  inline void SetCurveLength(G4double nCurve_s);
139  // Distance along curve.
140  inline void SetKineticEnergy(G4double nEnergy);
141  // Does not modify momentum.
142 
143  inline void SetLabTimeOfFlight(G4double tofLab);
144  inline void SetProperTimeOfFlight(G4double tofProper);
145  // Modifiers
146 
147  public: // without description
148 
149  enum { ncompSVEC = 12 };
150  // Needed and should be used only for RK integration driver
151 
152  inline void DumpToArray(G4double valArr[ncompSVEC]) const;
153  void LoadFromArray(const G4double valArr[ncompSVEC],
154  G4int noVarsIntegrated);
155  friend std::ostream&
156  operator<<( std::ostream& os, const G4FieldTrack& SixVec);
157 
158  public: // Obsolete methods -- due to potential confusion with PDG spin
159  inline void InitialiseSpin( const G4ThreeVector& vecPolarization )
160  { SetPolarization( vecPolarization ); }
161  inline G4ThreeVector GetSpin() const { return GetPolarization(); }
162  inline void SetSpin(G4ThreeVector vSpin){ SetPolarization(vSpin); }
163 
164  private: // Implementation method -- Obsolete
165  inline G4FieldTrack& SetCurvePnt(const G4ThreeVector& pPosition,
166  const G4ThreeVector& pMomentum,
167  G4double s_curve );
168  private:
169 
171  G4double fDistanceAlongCurve; // distance along curve of point
178  // G4double fInitialMomentumMag; // At 'track' creation.
179  // G4double fLastMomentumMag; // From last Update (for checking.)
180 
182 
183  public: // Access
184 
185  const G4ChargeState* GetChargeState() const { return &fChargeState; }
186 };
187 
188 #include "G4FieldTrack.icc"
189 
190 #endif /* End of ifndef G4FieldTrack_HH */
void SetPDGSpin(G4double pdgSpin)
Definition: G4FieldTrack.hh:93
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
void SetKineticEnergy(G4double nEnergy)
G4double GetProperTimeOfFlight() const
void SetProperTimeOfFlight(G4double tofProper)
G4double fKineticEnergy
G4double GetRestMass() const
void SetChargeAndMoments(G4double charge, G4double magnetic_dipole_moment=DBL_MAX, G4double electric_dipole_moment=DBL_MAX, G4double magnetic_charge=DBL_MAX)
void DumpToArray(G4double valArr[ncompSVEC]) const
G4ChargeState fChargeState
void UpdateFourMomentum(G4double kineticEnergy, const G4ThreeVector &momentumDirection)
const G4ChargeState * GetChargeState() const
void SetMomentum(G4ThreeVector nMomDir)
G4FieldTrack & SetCurvePnt(const G4ThreeVector &pPosition, const G4ThreeVector &pMomentum, G4double s_curve)
void SetPosition(G4ThreeVector nPos)
void SetPDGSpin(G4double spin)
void SetLabTimeOfFlight(G4double tofLab)
G4ThreeVector GetSpin() const
G4ThreeVector GetPolarization() const
G4ThreeVector GetPosition() const
double G4double
Definition: G4Types.hh:76
G4FieldTrack(const G4ThreeVector &pPosition, G4double LaboratoryTimeOfFlight, const G4ThreeVector &pMomentumDirection, G4double kineticEnergy, G4double restMass_c2, G4double charge, const G4ThreeVector &polarization, G4double magnetic_dipole_moment=0.0, G4double curve_length=0.0, G4double PDGspin=-1.0)
Definition: G4FieldTrack.cc:72
G4ThreeVector GetMomentumDirection() const
G4double GetKineticEnergy() const
G4double GetPDGSpin() const
void SetCurveLength(G4double nCurve_s)
G4FieldTrack & operator=(const G4FieldTrack &rStVec)
G4double fRestMass_c2
G4ThreeVector GetMomentum() const
G4double GetPDGSpin()
Definition: G4FieldTrack.hh:94
void SetSpin(G4ThreeVector vSpin)
G4double fProperTimeOfFlight
G4double fLabTimeOfFlight
void SetRestMass(G4double Mass_c2)
void SetPolarization(const G4ThreeVector &vecPol)
G4ThreeVector fMomentumDir
int G4int
Definition: G4Types.hh:78
void SetMomentumDir(G4ThreeVector nMomDir)
G4double fDistanceAlongCurve
void InitialiseSpin(const G4ThreeVector &vecPolarization)
const G4ThreeVector & GetMomentumDir() const
G4ThreeVector fPolarization
friend std::ostream & operator<<(std::ostream &os, const G4FieldTrack &SixVec)
Definition: G4FieldTrack.cc:33
void UpdateState(const G4ThreeVector &pPosition, G4double LaboratoryTimeOfFlight, const G4ThreeVector &pMomentumDirection, G4double kineticEnergy)
G4double GetCharge() const
#define DBL_MAX
Definition: templates.hh:83
G4double GetLabTimeOfFlight() const
G4double GetCurveLength() const
G4double SixVector[6]