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