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