Geant4  10.01.p01
G4FieldTrack.cc
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 // $Id: G4FieldTrack.cc 88944 2015-03-16 16:11:10Z gcosmo $
28 //
29 // -------------------------------------------------------------------
30 
31 #include "G4FieldTrack.hh"
32 
33 std::ostream& operator<<( std::ostream& os, const G4FieldTrack& SixVec)
34 {
35  const G4double *SixV = SixVec.SixVector;
36  os << " ( ";
37  os << " X= " << SixV[0] << " " << SixV[1] << " "
38  << SixV[2] << " "; // Position
39  os << " P= " << SixV[3] << " " << SixV[4] << " "
40  << SixV[5] << " "; // Momentum
41  os << " Pmag= "
42  << G4ThreeVector(SixV[3], SixV[4], SixV[5]).mag(); // mom magnitude
43  os << " Ekin= " << SixVec.fKineticEnergy ;
44  os << " m0= " << SixVec.fRestMass_c2;
45  os << " Pdir= " << SixVec.fMomentumDir.mag();
46  os << " PolV= " << SixVec.GetPolarization();
47  os << " l= " << SixVec.GetCurveLength();
48  os << " t_lab= " << SixVec.fLabTimeOfFlight;
49  os << " t_proper= " << SixVec.fProperTimeOfFlight ;
50  os << " ) ";
51  return os;
52 }
53 
55  G4double LaboratoryTimeOfFlight,
56  const G4ThreeVector& pMomentumDirection,
57  G4double kineticEnergy,
58  G4double restMass_c2,
59  G4double charge,
60  const G4ThreeVector& vecPolarization,
61  G4double magnetic_dipole_moment,
62  G4double curve_length,
63  G4double pdgSpin )
64 : fDistanceAlongCurve(curve_length),
65  fKineticEnergy(kineticEnergy),
66  fRestMass_c2(restMass_c2),
67  fLabTimeOfFlight(LaboratoryTimeOfFlight),
68  fProperTimeOfFlight(0.),
69  // fMomentumDir(pMomentumDirection),
70  fChargeState( charge, magnetic_dipole_moment, pdgSpin )
71  // fChargeState( charge, magnetic_dipole_moment ) ,
72  // fPDGSpin( pdgSpin )
73 {
74  UpdateFourMomentum( kineticEnergy, pMomentumDirection );
75  // Sets momentum direction as well.
76 
77  SetPosition( pPosition );
78 
79  SetPolarization( vecPolarization );
80 }
81 
83  const G4ThreeVector& pMomentumDirection,
84  G4double curve_length,
85  G4double kineticEnergy,
86  const G4double restMass_c2,
87  G4double, // velocity
88  G4double pLaboratoryTimeOfFlight,
89  G4double pProperTimeOfFlight,
90  const G4ThreeVector* pPolarization,
91  G4double pdgSpin )
92  : fDistanceAlongCurve(curve_length),
93  fKineticEnergy(kineticEnergy),
94  fRestMass_c2(restMass_c2),
95  fLabTimeOfFlight(pLaboratoryTimeOfFlight),
96  fProperTimeOfFlight(pProperTimeOfFlight),
97  fChargeState( DBL_MAX, DBL_MAX, -1.0 ) // charge not set
98 {
99  UpdateFourMomentum( kineticEnergy, pMomentumDirection );
100  // Sets momentum direction as well.
101 
102  SetPosition( pPosition );
103  fChargeState.SetPDGSpin( pdgSpin );
104 
105  G4ThreeVector PolarVec(0.0, 0.0, 0.0);
106  if( pPolarization ) { PolarVec= *pPolarization; }
107  SetPolarization( PolarVec );
108 }
109 
110 G4FieldTrack::G4FieldTrack( char ) // Nothing is set !!
111  : fKineticEnergy(0.), fRestMass_c2(0.), fLabTimeOfFlight(0.),
112  fProperTimeOfFlight(0.), fChargeState( DBL_MAX , DBL_MAX, -1 )
113 {
114  G4ThreeVector Zero(0.0, 0.0, 0.0);
115  SetCurvePnt( Zero, Zero, 0.0 );
116  SetPolarization( Zero );
117  // fInitialMomentumMag= 0.00; // Invalid
118  // fLastMomentumMag= 0.0;
119 }
120 
121 void G4FieldTrack::
123  G4double magnetic_dipole_moment, // default= DBL_MAX - do not change
124  G4double electric_dipole_moment, // ditto
125  G4double magnetic_charge ) // ditto
126 {
128  magnetic_dipole_moment,
129  electric_dipole_moment,
130  magnetic_charge );
131 
132  // NOTE: Leaves Spin unchanged !
133  //
134  // G4double pdgSpin= fChargeState.GetSpin(); // New Property of ChargeState (not well documented! )
135 
136  // IDEA: Improve the implementation using handles
137  // -- and handle to the old one (which can be shared by other copies) and
138  // must not be left to hang loose
139  //
140  // fpChargeState= new G4ChargeState( charge, magnetic_dipole_moment,
141  // electric_dipole_moment, magnetic_charge );
142 }
143 
144 // Load values from array
145 //
146 // note that momentum direction must-be/is normalised
147 
148 void G4FieldTrack::LoadFromArray(const G4double valArrIn[ncompSVEC],
149  G4int noVarsIntegrated)
150 {
151  G4int i;
152 
153  // Fill the variables not integrated with zero -- so it's clear !!
154  G4double valArr[ncompSVEC];
155  for( i=0; i<noVarsIntegrated; i++){
156  valArr[i]= valArrIn[i];
157  }
158  for( i=noVarsIntegrated; i<ncompSVEC; i++) {
159  valArr[i]= 0.0;
160  }
161 
162  SixVector[0]=valArr[0];
163  SixVector[1]=valArr[1];
164  SixVector[2]=valArr[2];
165  SixVector[3]=valArr[3];
166  SixVector[4]=valArr[4];
167  SixVector[5]=valArr[5];
168 
169  G4ThreeVector Momentum(valArr[3],valArr[4],valArr[5]);
170 
171  G4double momentum_square= Momentum.mag2();
172  fMomentumDir= Momentum.unit();
173 
174  fKineticEnergy = momentum_square /
175  (std::sqrt(momentum_square+fRestMass_c2*fRestMass_c2)
176  + fRestMass_c2 );
177  // The above equation is stable for small and large momenta
178 
179  // The following components may or may not be
180  // integrated over -- integration is optional
181  // fKineticEnergy= valArr[6];
182 
183  fLabTimeOfFlight=valArr[7];
184  fProperTimeOfFlight=valArr[8];
185  G4ThreeVector vecPolarization= G4ThreeVector(valArr[9],valArr[10],valArr[11]);
186  SetPolarization( vecPolarization );
187 
188  // fMomentumDir=G4ThreeVector(valArr[13],valArr[14],valArr[15]);
189  // fDistanceAlongCurve= valArr[];
190 }
191 
void SetPosition(G4ThreeVector nPos)
G4double GetCurveLength() const
G4double SixVector[6]
CLHEP::Hep3Vector G4ThreeVector
void SetPolarization(const G4ThreeVector &vecPol)
G4double fLabTimeOfFlight
int G4int
Definition: G4Types.hh:78
void UpdateFourMomentum(G4double kineticEnergy, const G4ThreeVector &momentumDirection)
void SetChargesAndMoments(G4double charge, G4double magnetic_dipole_moment, G4double electric_dipole_moment, G4double magnetic_charge)
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:54
G4ThreeVector fMomentumDir
std::ostream & operator<<(std::ostream &os, const G4FieldTrack &SixVec)
Definition: G4FieldTrack.cc:33
G4ThreeVector GetPolarization() const
G4double fKineticEnergy
void LoadFromArray(const G4double valArr[ncompSVEC], G4int noVarsIntegrated)
G4double fRestMass_c2
G4double fProperTimeOfFlight
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
void SetChargeAndMoments(G4double charge, G4double magnetic_dipole_moment=DBL_MAX, G4double electric_dipole_moment=DBL_MAX, G4double magnetic_charge=DBL_MAX)
G4FieldTrack & SetCurvePnt(const G4ThreeVector &pPosition, const G4ThreeVector &pMomentum, G4double s_curve)
G4ChargeState fChargeState
void SetPDGSpin(G4double spin)