Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ErrorFreeTrajState Class Reference

#include <G4ErrorFreeTrajState.hh>

+ Inheritance diagram for G4ErrorFreeTrajState:

Public Member Functions

 G4ErrorFreeTrajState ()
 
 G4ErrorFreeTrajState (const G4String &partName, const G4Point3D &pos, const G4Vector3D &mom, const G4ErrorTrajErr &errmat=G4ErrorTrajErr(5, 0))
 
 G4ErrorFreeTrajState (const G4ErrorSurfaceTrajState &tpOS)
 
 ~G4ErrorFreeTrajState ()
 
virtual G4int Update (const G4Track *aTrack)
 
virtual G4int PropagateError (const G4Track *aTrack)
 
virtual void Dump (std::ostream &out=G4cout) const
 
virtual void SetPosition (const G4Point3D pos)
 
virtual void SetMomentum (const G4Vector3D &mom)
 
void SetParameters (const G4Point3D &pos, const G4Vector3D &mom)
 
G4ErrorFreeTrajParam GetParameters () const
 
G4ErrorMatrix GetTransfMat () const
 
- Public Member Functions inherited from G4ErrorTrajState
 G4ErrorTrajState ()
 
 G4ErrorTrajState (const G4String &partType, const G4Point3D &pos, const G4Vector3D &mom, const G4ErrorTrajErr &errmat=G4ErrorTrajErr(5, 0))
 
 G4ErrorTrajState (const G4ErrorTrajState &)
 
 G4ErrorTrajState (G4ErrorTrajState &&)
 
virtual ~G4ErrorTrajState ()
 
G4ErrorTrajStateoperator= (const G4ErrorTrajState &)
 
G4ErrorTrajStateoperator= (G4ErrorTrajState &&)
 
void SetData (const G4String &partType, const G4Point3D &pos, const G4Vector3D &mom)
 
void BuildCharge ()
 
virtual G4int PropagateError (const G4Track *)
 
virtual G4int Update (const G4Track *)
 
void UpdatePosMom (const G4Point3D &pos, const G4Vector3D &mom)
 
void DumpPosMomError (std::ostream &out=G4cout) const
 
virtual void Dump (std::ostream &out=G4cout) const =0
 
const G4StringGetParticleType () const
 
void SetParticleType (const G4String &partType)
 
G4Point3D GetPosition () const
 
virtual void SetPosition (const G4Point3D pos)
 
G4Vector3D GetMomentum () const
 
virtual void SetMomentum (const G4Vector3D &mom)
 
G4ErrorTrajErr GetError () const
 
virtual void SetError (G4ErrorTrajErr em)
 
G4TrackGetG4Track () const
 
void SetG4Track (G4Track *trk)
 
G4double GetCharge () const
 
void SetCharge (G4double ch)
 
virtual G4eTSType GetTSType () const
 

Friends

std::ostream & operator<< (std::ostream &, const G4ErrorFreeTrajState &ts)
 

Additional Inherited Members

- Protected Attributes inherited from G4ErrorTrajState
G4String fParticleType
 
G4Point3D fPosition
 
G4Vector3D fMomentum
 
G4double fCharge = 0.
 
G4ErrorTrajErr fError
 
G4eTSType theTSType
 
G4TracktheG4Track = nullptr
 
G4int iverbose = 0
 

Detailed Description

Definition at line 64 of file G4ErrorFreeTrajState.hh.

Constructor & Destructor Documentation

◆ G4ErrorFreeTrajState() [1/3]

G4ErrorFreeTrajState::G4ErrorFreeTrajState ( )
inline

Definition at line 68 of file G4ErrorFreeTrajState.hh.

68: theFirstStep(true) {}

◆ G4ErrorFreeTrajState() [2/3]

G4ErrorFreeTrajState::G4ErrorFreeTrajState ( const G4String partName,
const G4Point3D pos,
const G4Vector3D mom,
const G4ErrorTrajErr errmat = G4ErrorTrajErr(5,0) 
)

Definition at line 48 of file G4ErrorFreeTrajState.cc.

48 : G4ErrorTrajState( partName, pos, mom, errmat )
49{
50 fTrajParam = G4ErrorFreeTrajParam( pos, mom );
51 Init();
52}

◆ G4ErrorFreeTrajState() [3/3]

G4ErrorFreeTrajState::G4ErrorFreeTrajState ( const G4ErrorSurfaceTrajState tpOS)

Definition at line 56 of file G4ErrorFreeTrajState.cc.

56 : G4ErrorTrajState( tpSD.GetParticleType(), tpSD.GetPosition(), tpSD.GetMomentum() )
57{
58 // G4ThreeVector planeNormal = tpSD.GetPlaneNormal();
59 // G4double fPt = tpSD.GetMomentum()*planeNormal;//mom projected on normal to plane
60 // G4ErrorSurfaceTrajParam tpSDparam = tpSD.GetParameters();
61 // G4ThreeVector Psc = fPt * planeNormal + tpSDparam.GetPU()*tpSDparam.GetVectorU() + tpSD.GetPV()*tpSD.GetVectorW();
62
64 Init();
65
66 //----- Get the error matrix in SC coordinates
67 G4ErrorSurfaceTrajParam tpSDparam = tpSD.GetParameters();
68 G4double mom = fMomentum.mag();
69 G4double mom2 = fMomentum.mag2();
70 G4double TVW1 = std::sqrt( mom2 / ( mom2 + tpSDparam.GetPV()*tpSDparam.GetPV() + tpSDparam.GetPW()*tpSDparam.GetPW()) );
71 G4ThreeVector vTVW( TVW1, tpSDparam.GetPV()/mom * TVW1, tpSDparam.GetPW()/mom * TVW1 );
72 G4Vector3D vectorU = tpSDparam.GetVectorV().cross( tpSDparam.GetVectorW() );
73 G4Vector3D vTN = vTVW.x()*vectorU + vTVW.y()*tpSDparam.GetVectorV() + vTVW.z()*tpSDparam.GetVectorW();
74
75#ifdef G4EVERBOSE
76 if( iverbose >= 5){
77 G4double pc2 = std::asin( vTN.z() );
78 G4double pc3 = std::atan (vTN.y()/vTN.x());
79
80 G4cout << " CHECK: pc2 " << pc2 << " = " << GetParameters().GetLambda() << " diff " << pc2-GetParameters().GetLambda() << G4endl;
81 G4cout << " CHECK: pc3 " << pc3 << " = " << GetParameters().GetPhi() << " diff " << pc3-GetParameters().GetPhi() << G4endl;
82 }
83#endif
84
85 //--- Get the unit vectors perp to P
86 G4double cosl = std::cos( GetParameters().GetLambda() );
87 if (cosl < 1.E-30) cosl = 1.E-30;
88 G4double cosl1 = 1./cosl;
89 G4Vector3D vUN(-vTN.y()*cosl1, vTN.x()*cosl1, 0. );
90 G4Vector3D vVN(-vTN.z()*vUN.y(), vTN.z()*vUN.x(), cosl );
91
92 G4Vector3D vUperp = G4Vector3D( -fMomentum.y(), fMomentum.x(), 0.);
93 G4Vector3D vVperp = vUperp.cross( fMomentum );
94 vUperp *= 1./vUperp.mag();
95 vVperp *= 1./vVperp.mag();
96
97#ifdef G4EVERBOSE
98 if( iverbose >= 5 ){
99 G4cout << " CHECK: vUN " << vUN << " = " << vUperp << " diff " << (vUN-vUperp).mag() << G4endl;
100 G4cout << " CHECK: vVN " << vVN << " = " << vVperp << " diff " << (vVN-vVperp).mag() << G4endl;
101 }
102#endif
103
104 //get the dot products of vectors perpendicular to direction and vector defining SD plane
105 G4double dUU = vUperp * tpSD.GetVectorV();
106 G4double dUV = vUperp * tpSD.GetVectorW();
107 G4double dVU = vVperp * tpSD.GetVectorV();
108 G4double dVV = vVperp * tpSD.GetVectorW();
109
110 //--- Get transformation first
111 G4ErrorMatrix transfM(5, 5, 1 );
112 //--- Get magnetic field
114 G4ThreeVector dir = fTrajParam.GetDirection();
115 G4double invCosTheta = 1./std::cos( dir.theta() );
116 G4cout << " dir="<<dir<<" invCosTheta "<<invCosTheta << G4endl;
117
118 if( fCharge != 0 && field ) {
119 G4double pos1[3]; pos1[0] = fPosition.x()*cm; pos1[1] = fPosition.y()*cm; pos1[2] = fPosition.z()*cm;
120 G4double h1[3];
121 field->GetFieldValue( pos1, h1 );
122 G4ThreeVector HPre = G4ThreeVector( h1[0], h1[1], h1[2] ) / tesla *10.;
123 G4double magHPre = HPre.mag();
124 G4double invP = 1./fMomentum.mag();
125 G4double magHPreM = magHPre * invP;
126 if( magHPre != 0. ) {
127 G4double magHPreM2 = fCharge / magHPre;
128
129 G4double Q = -magHPreM * c_light;
130 G4double sinz = -HPre*vUperp * magHPreM2;
131 G4double cosz = HPre*vVperp * magHPreM2;
132
133 transfM[1][3] = -Q*dir.y()*sinz;
134 transfM[1][4] = -Q*dir.z()*sinz;
135 transfM[2][3] = -Q*dir.y()*cosz*invCosTheta;
136 transfM[2][4] = -Q*dir.z()*cosz*invCosTheta;
137 }
138 }
139
140 transfM[0][0] = 1.;
141 transfM[1][1] = dir.x()*dVU;
142 transfM[1][2] = dir.x()*dVV;
143 transfM[2][1] = dir.x()*dUU*invCosTheta;
144 transfM[2][2] = dir.x()*dUV*invCosTheta;
145 transfM[3][3] = dUU;
146 transfM[3][4] = dUV;
147 transfM[4][3] = dVU;
148 transfM[4][4] = dVV;
149
150 fError = G4ErrorTrajErr( tpSD.GetError().similarity( transfM ) );
151
152#ifdef G4EVERBOSE
153 if( iverbose >= 1) G4cout << "error matrix SD2SC " << fError << G4endl;
154 if( iverbose >= 4) G4cout << "G4ErrorFreeTrajState from SD " << *this << G4endl;
155#endif
156}
G4ErrorSymMatrix G4ErrorTrajErr
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double theta() const
double x() const
double y() const
double mag() const
G4double GetLambda() const
G4Vector3D GetDirection() const
G4ErrorFreeTrajParam GetParameters() const
G4ErrorTrajErr fError
const G4Field * GetDetectorField() const
virtual void GetFieldValue(const G4double Point[4], G4double *fieldArr) const =0
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

◆ ~G4ErrorFreeTrajState()

G4ErrorFreeTrajState::~G4ErrorFreeTrajState ( )
inline

Definition at line 78 of file G4ErrorFreeTrajState.hh.

78{}

Member Function Documentation

◆ Dump()

void G4ErrorFreeTrajState::Dump ( std::ostream &  out = G4cout) const
virtual

Implements G4ErrorTrajState.

Definition at line 169 of file G4ErrorFreeTrajState.cc.

170{
171 out << *this;
172}

◆ GetParameters()

G4ErrorFreeTrajParam G4ErrorFreeTrajState::GetParameters ( ) const
inline

Definition at line 107 of file G4ErrorFreeTrajState.hh.

108 { return fTrajParam; }

Referenced by G4ErrorSurfaceTrajState::BuildErrorMatrix(), and G4ErrorFreeTrajState().

◆ GetTransfMat()

G4ErrorMatrix G4ErrorFreeTrajState::GetTransfMat ( ) const
inline

Definition at line 110 of file G4ErrorFreeTrajState.hh.

111 { return theTransfMat; }

◆ PropagateError()

G4int G4ErrorFreeTrajState::PropagateError ( const G4Track aTrack)
virtual

Reimplemented from G4ErrorTrajState.

Definition at line 202 of file G4ErrorFreeTrajState.cc.

203{
204 G4double stepLengthCm = aTrack->GetStep()->GetStepLength()/cm;
205 if( G4ErrorPropagatorData::GetErrorPropagatorData()->GetStage() == G4ErrorStage_Deflation ) stepLengthCm *= -1.;
206
208
209 if( std::fabs(stepLengthCm) <= kCarTolerance/cm ) return 0;
210
211#ifdef G4EVERBOSE
212 if( iverbose >= 2 )G4cout << " G4ErrorFreeTrajState::PropagateError " << G4endl;
213 G4cout << "G4EP: iverbose="<< iverbose << G4endl;
214#endif
215
216 // * *** ERROR PROPAGATION ON A HELIX ASSUMING SC VARIABLES
217 G4Point3D vposPost = aTrack->GetPosition()/cm;
218 G4Vector3D vpPost = aTrack->GetMomentum()/GeV;
219 // G4Point3D vposPre = fPosition/cm;
220 // G4Vector3D vpPre = fMomentum/GeV;
221 G4Point3D vposPre = aTrack->GetStep()->GetPreStepPoint()->GetPosition()/cm;
222 G4Vector3D vpPre = aTrack->GetStep()->GetPreStepPoint()->GetMomentum()/GeV;
223 //correct to avoid propagation along Z
224 if( vpPre.mag() == vpPre.z() ) vpPre.setX( 1.E-6*MeV );
225 if( vpPost.mag() == vpPost.z() ) vpPost.setX( 1.E-6*MeV );
226
227 G4double pPre = vpPre.mag();
228 G4double pPost = vpPost.mag();
229#ifdef G4EVERBOSE
230 if( iverbose >= 2 ) {
231 G4cout << "G4EP: vposPre " << vposPre << G4endl
232 << "G4EP: vposPost " << vposPost << G4endl;
233 G4cout << "G4EP: vpPre " << vpPre << G4endl
234 << "G4EP: vpPost " << vpPost << G4endl;
235 G4cout << " err start step " << fError << G4endl;
236 G4cout << "G4EP: stepLengthCm " << stepLengthCm << G4endl;
237 }
238#endif
239
240 if( pPre == 0. || pPost == 0 ) return 2;
241 G4double pInvPre = 1./pPre;
242 G4double pInvPost = 1./pPost;
243 G4double deltaPInv = pInvPost - pInvPre;
244 if( iverbose >= 2 ) G4cout << "G4EP: pInvPre" << pInvPre<< " pInvPost:" << pInvPost<<" deltaPInv:" << deltaPInv<< G4endl;
245
246
247 G4Vector3D vpPreNorm = vpPre * pInvPre;
248 G4Vector3D vpPostNorm = vpPost * pInvPost;
249 if( iverbose >= 2 ) G4cout << "G4EP: vpPreNorm " << vpPreNorm << " vpPostNorm " << vpPostNorm << G4endl;
250 //return if propagation along Z??
251 if( 1. - std::fabs(vpPreNorm.z()) < kCarTolerance ) return 4;
252 if( 1. - std::fabs(vpPostNorm.z()) < kCarTolerance ) return 4;
253 G4double sinpPre = std::sin( vpPreNorm.theta() ); //cosine perpendicular to pPre = sine pPre
254 G4double sinpPost = std::sin( vpPostNorm.theta() ); //cosine perpendicular to pPost = sine pPost
255 G4double sinpPostInv = 1./std::sin( vpPostNorm.theta() );
256
257#ifdef G4EVERBOSE
258 if( iverbose >= 2 ) G4cout << "G4EP: cosl " << sinpPre << " cosl0 " << sinpPost << G4endl;
259#endif
260 //* *** DEFINE TRANSFORMATION MATRIX BETWEEN X1 AND X2 FOR
261 //* *** NEUTRAL PARTICLE OR FIELDFREE REGION
262 G4ErrorMatrix transf(5, 5, 0 );
263
264 transf[3][2] = stepLengthCm * sinpPost;
265 transf[4][1] = stepLengthCm;
266 for( size_t ii=0;ii < 5; ii++ ){
267 transf[ii][ii] = 1.;
268 }
269#ifdef G4EVERBOSE
270 if( iverbose >= 2 ) {
271 G4cout << "G4EP: transf matrix neutral " << transf;
272 }
273#endif
274
275 // charge X propagation direction
276 G4double charge = aTrack->GetDynamicParticle()->GetCharge();
278 charge *= -1.;
279 }
280 // G4cout << " charge " << charge << G4endl;
281 //t check if particle has charge
282 //t if( charge == 0 ) goto 45;
283 // check if the magnetic field is = 0.
284
285 //position is from geant4, it is assumed to be in mm (for debugging, eventually it will not be transformed)
286 //it is assumed vposPre[] is in cm and pos1[] is in mm.
287 G4double pos1[3]; pos1[0] = vposPre.x()*cm; pos1[1] = vposPre.y()*cm; pos1[2] = vposPre.z()*cm;
288 G4double pos2[3]; pos2[0] = vposPost.x()*cm; pos2[1] = vposPost.y()*cm; pos2[2] = vposPost.z()*cm;
289 G4double h1[3], h2[3];
290
292 if( !field ) return 0; //goto 45
293
294
295
296 // calculate transformation except it NEUTRAL PARTICLE OR FIELDFREE REGION
297 if( charge != 0. && field ) {
298
299 field->GetFieldValue( pos1, h1 ); //here pos1[], pos2[] are in mm, not changed
300 field->GetFieldValue( pos2, h2 );
301 G4ThreeVector HPre = G4ThreeVector( h1[0], h1[1], h1[2] ) / tesla *10.; //10. is to get same dimensions as GEANT3 (kilogauss)
302 G4ThreeVector HPost= G4ThreeVector( h2[0], h2[1], h2[2] ) / tesla *10.;
303 G4double magHPre = HPre.mag();
304 G4double magHPost = HPost.mag();
305#ifdef G4EVERBOSE
306 if( iverbose >= 2 ) {
307 G4cout << "G4EP: h1 = "
308 << h1[0] << ", " << h1[1] << ", " << h1[2] << G4endl;
309 G4cout << "G4EP: pos1/mm = "
310 << pos1[0] << ", " << pos1[1] << ", " << pos1[2] << G4endl;
311 G4cout << "G4EP: pos2/mm = "
312 << pos2[0] << ", " << pos2[1] << ", " << pos2[2] << G4endl;
313 G4cout << "G4EP: B-filed in KGauss HPre " << HPre << G4endl
314 << "G4EP: in KGauss HPost " << HPost << G4endl;
315 }
316#endif
317
318 if( magHPre + magHPost != 0. ) {
319
320 //* *** CHECK WHETHER H*ALFA/P IS TOO DIFFERENT AT X1 AND X2
322 if( magHPost != 0. ){
323 gam = HPost * vpPostNorm / magHPost;
324 }else {
325 gam = HPre * vpPreNorm / magHPre;
326 }
327
328 // G4eMagneticLimitsProcess will limit the step, but based on an straight line trajectory
329 G4double alphaSqr = 1. - gam * gam;
330 G4double diffHSqr = ( HPre * pInvPre - HPost * pInvPost ).mag2();
331 G4double delhp6Sqr = 300.*300.;
332#ifdef G4EVERBOSE
333 if( iverbose >= 2 ) {
334 G4cout << " G4EP: gam " << gam << " alphaSqr " << alphaSqr
335 << " diffHSqr " << diffHSqr << G4endl;
336 G4cout << " alpha= " << std::sqrt(alphaSqr) << G4endl;
337 }
338#endif
339 if( diffHSqr * alphaSqr > delhp6Sqr ) return 3;
340
341
342 //* *** DEFINE AVERAGE MAGNETIC FIELD AND GRADIENT
343 G4double pInvAver = 1./(pInvPre + pInvPost );
344 G4double CFACT8 = 2.997925E-4;
345 //G4double HAver
346 G4ThreeVector vHAverNorm( (HPre*pInvPre + HPost*pInvPost ) * pInvAver * charge * CFACT8 );
347 G4double HAver = vHAverNorm.mag();
348 G4double invHAver = 1./HAver;
349 vHAverNorm *= invHAver;
350#ifdef G4EVERBOSE
351 if( iverbose >= 2 ) G4cout << " G4EP: HaverNorm " << vHAverNorm << " magHAver " << HAver << " charge " << charge<< G4endl;
352#endif
353
354 G4double pAver = (pPre+pPost)*0.5;
355 G4double QAver = -HAver/pAver;
356 G4double thetaAver = QAver * stepLengthCm;
357 G4double sinThetaAver = std::sin(thetaAver);
358 G4double cosThetaAver = std::cos(thetaAver);
359 G4double gamma = vHAverNorm * vpPostNorm;
360 G4ThreeVector AN2 = vHAverNorm.cross( vpPostNorm );
361
362#ifdef G4EVERBOSE
363 if( iverbose >= 2 ) G4cout << " G4EP: AN2 " << AN2 << " gamma:"<<gamma<< " theta="<< thetaAver<<G4endl;
364#endif
365 G4double AU = 1./vpPreNorm.perp();
366 //t G4ThreeVector vU( vpPreNorm.cross( G4ThreeVector(0.,0.,1.) ) * AU );
367 G4ThreeVector vUPre( -AU*vpPreNorm.y(),
368 AU*vpPreNorm.x(),
369 0. );
370 G4ThreeVector vVPre( -vpPreNorm.z()*vUPre.y(),
371 vpPreNorm.z()*vUPre.x(),
372 vpPreNorm.x()*vUPre.y() - vpPreNorm.y()*vUPre.x() );
373
374 //
375 AU = 1./vpPostNorm.perp();
376 //t G4ThreeVector vU( vpPostNorm.cross( G4ThreeVector(0.,0.,1.) ) * AU );
377 G4ThreeVector vUPost( -AU*vpPostNorm.y(),
378 AU*vpPostNorm.x(),
379 0. );
380 G4ThreeVector vVPost( -vpPostNorm.z()*vUPost.y(),
381 vpPostNorm.z()*vUPost.x(),
382 vpPostNorm.x()*vUPost.y() - vpPostNorm.y()*vUPost.x() );
383#ifdef G4EVERBOSE
384 G4cout << " vpPostNorm " << vpPostNorm << G4endl;
385 if( iverbose >= 2 ) G4cout << " G4EP: AU " << AU << " vUPre " << vUPre << " vVPre " << vVPre << " vUPost " << vUPost << " vVPost " << vVPost << G4endl;
386#endif
387 G4Point3D deltaPos( vposPre - vposPost );
388
389 // * *** COMPLETE TRANSFORMATION MATRIX BETWEEN ERRORS AT X1 AND X2
390 // * *** FIELD GRADIENT PERPENDICULAR TO TRACK IS PRESENTLY NOT
391 // * *** TAKEN INTO ACCOUNT
392
393 G4double QP = QAver * pAver; // = -HAver
394#ifdef G4EVERBOSE
395 if( iverbose >= 2) G4cout << " G4EP: QP " << QP << " QAver " << QAver << " pAver " << pAver << G4endl;
396#endif
397 G4double ANV = -( vHAverNorm.x()*vUPost.x() + vHAverNorm.y()*vUPost.y() );
398 G4double ANU = ( vHAverNorm.x()*vVPost.x() + vHAverNorm.y()*vVPost.y() + vHAverNorm.z()*vVPost.z() );
399 G4double OMcosThetaAver = 1. - cosThetaAver;
400#ifdef G4EVERBOSE
401 if( iverbose >= 2) G4cout << "G4EP: OMcosThetaAver " << OMcosThetaAver << " cosThetaAver " << cosThetaAver << " thetaAver " << thetaAver << " QAver " << QAver << " stepLengthCm " << stepLengthCm << G4endl;
402#endif
403 G4double TMSINT = thetaAver - sinThetaAver;
404#ifdef G4EVERBOSE
405 if( iverbose >= 2 ) G4cout << " G4EP: ANV " << ANV << " ANU " << ANU << G4endl;
406#endif
407
408 G4ThreeVector vHUPre( -vHAverNorm.z() * vUPre.y(),
409 vHAverNorm.z() * vUPre.x(),
410 vHAverNorm.x() * vUPre.y() - vHAverNorm.y() * vUPre.x() );
411#ifdef G4EVERBOSE
412 // if( iverbose >= 2 ) G4cout << "G4EP: HUPre(1) " << vHUPre.x() << " " << vHAverNorm.z() << " " << vUPre.y() << G4endl;
413#endif
414 G4ThreeVector vHVPre( vHAverNorm.y() * vVPre.z() - vHAverNorm.z() * vVPre.y(),
415 vHAverNorm.z() * vVPre.x() - vHAverNorm.x() * vVPre.z(),
416 vHAverNorm.x() * vVPre.y() - vHAverNorm.y() * vVPre.x() );
417#ifdef G4EVERBOSE
418 if( iverbose >= 2 ) G4cout << " G4EP: HUPre " << vHUPre << " HVPre " << vHVPre << G4endl;
419#endif
420
421 //------------------- COMPUTE MATRIX
422 //---------- 1/P
423
424 transf[0][0] = 1.-deltaPInv*pAver*(1.+(vpPostNorm.x()*deltaPos.x()+vpPostNorm.y()*deltaPos.y()+vpPostNorm.z()*deltaPos.z())/stepLengthCm)
425 +2.*deltaPInv*pAver;
426
427 transf[0][1] = -deltaPInv/thetaAver*
428 ( TMSINT*gamma*(vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z()) +
429 sinThetaAver*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z()) +
430 OMcosThetaAver*(vHVPre.x()*vpPostNorm.x()+vHVPre.y()*vpPostNorm.y()+vHVPre.z()*vpPostNorm.z()) );
431
432 transf[0][2] = -sinpPre*deltaPInv/thetaAver*
433 ( TMSINT*gamma*(vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() ) +
434 sinThetaAver*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() ) +
435 OMcosThetaAver*(vHUPre.x()*vpPostNorm.x()+vHUPre.y()*vpPostNorm.y()+vHUPre.z()*vpPostNorm.z()) );
436
437 transf[0][3] = -deltaPInv/stepLengthCm*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() );
438
439 transf[0][4] = -deltaPInv/stepLengthCm*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z());
440
441 // *** Lambda
442 transf[1][0] = -QP*ANV*(vpPostNorm.x()*deltaPos.x()+vpPostNorm.y()*deltaPos.y()+vpPostNorm.z()*deltaPos.z())
443 *(1.+deltaPInv*pAver);
444#ifdef G4EVERBOSE
445 if(iverbose >= 3) G4cout << "ctransf10= " << transf[1][0] << " " << -QP<< " " << ANV<< " " << vpPostNorm.x()<< " " << deltaPos.x()<< " " << vpPostNorm.y()<< " " << deltaPos.y()<< " " << vpPostNorm.z()<< " " << deltaPos.z()
446 << " " << deltaPInv<< " " << pAver << G4endl;
447#endif
448
449 transf[1][1] = cosThetaAver*(vVPre.x()*vVPost.x()+vVPre.y()*vVPost.y()+vVPre.z()*vVPost.z()) +
450 sinThetaAver*(vHVPre.x()*vVPost.x()+vHVPre.y()*vVPost.y()+vHVPre.z()*vVPost.z()) +
451 OMcosThetaAver*(vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z())*
452 (vHAverNorm.x()*vVPost.x()+vHAverNorm.y()*vVPost.y()+vHAverNorm.z()*vVPost.z()) +
453 ANV*( -sinThetaAver*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z()) +
454 OMcosThetaAver*(vVPre.x()*AN2.x()+vVPre.y()*AN2.y()+vVPre.z()*AN2.z()) -
455 TMSINT*gamma*(vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z()) );
456
457 transf[1][2] = cosThetaAver*(vUPre.x()*vVPost.x()+vUPre.y()*vVPost.y() ) +
458 sinThetaAver*(vHUPre.x()*vVPost.x()+vHUPre.y()*vVPost.y()+vHUPre.z()*vVPost.z()) +
459 OMcosThetaAver*(vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() )*
460 (vHAverNorm.x()*vVPost.x()+vHAverNorm.y()*vVPost.y()+vHAverNorm.z()*vVPost.z()) +
461 ANV*( -sinThetaAver*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() ) +
462 OMcosThetaAver*(vUPre.x()*AN2.x()+vUPre.y()*AN2.y() ) -
463 TMSINT*gamma*(vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() ) );
464 transf[1][2] = sinpPre*transf[1][2];
465
466 transf[1][3] = -QAver*ANV*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() );
467
468 transf[1][4] = -QAver*ANV*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z());
469
470 // *** Phi
471
472 transf[2][0] = -QP*ANU*(vpPostNorm.x()*deltaPos.x()+vpPostNorm.y()*deltaPos.y()+vpPostNorm.z()*deltaPos.z())*sinpPostInv
473 *(1.+deltaPInv*pAver);
474#ifdef G4EVERBOSE
475 if(iverbose >= 3)G4cout <<"ctransf20= " << transf[2][0] <<" "<< -QP<<" "<<ANU<<" "<<vpPostNorm.x()<<" "<<deltaPos.x()<<" "<<vpPostNorm.y()<<" "<<deltaPos.y()<<" "<<vpPostNorm.z()<<" "<<deltaPos.z()<<" "<<sinpPostInv
476 <<" "<<deltaPInv<<" "<<pAver<< G4endl;
477#endif
478 transf[2][1] = cosThetaAver*(vVPre.x()*vUPost.x()+vVPre.y()*vUPost.y() ) +
479 sinThetaAver*(vHVPre.x()*vUPost.x()+vHVPre.y()*vUPost.y() ) +
480 OMcosThetaAver*(vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z())*
481 (vHAverNorm.x()*vUPost.x()+vHAverNorm.y()*vUPost.y() ) +
482 ANU*( -sinThetaAver*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z()) +
483 OMcosThetaAver*(vVPre.x()*AN2.x()+vVPre.y()*AN2.y()+vVPre.z()*AN2.z()) -
484 TMSINT*gamma*(vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z()) );
485 transf[2][1] = sinpPostInv*transf[2][1];
486
487 transf[2][2] = cosThetaAver*(vUPre.x()*vUPost.x()+vUPre.y()*vUPost.y() ) +
488 sinThetaAver*(vHUPre.x()*vUPost.x()+vHUPre.y()*vUPost.y() ) +
489 OMcosThetaAver*(vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() )*
490 (vHAverNorm.x()*vUPost.x()+vHAverNorm.y()*vUPost.y() ) +
491 ANU*( -sinThetaAver*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() ) +
492 OMcosThetaAver*(vUPre.x()*AN2.x()+vUPre.y()*AN2.y() ) -
493 TMSINT*gamma*(vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() ) );
494 transf[2][2] = sinpPostInv*sinpPre*transf[2][2];
495
496 transf[2][3] = -QAver*ANU*(vUPre.x()*vpPostNorm.x()+vUPre.y()*vpPostNorm.y() )*sinpPostInv;
497#ifdef G4EVERBOSE
498 if(iverbose >= 3)G4cout <<"ctransf23= " << transf[2][3] <<" "<< -QAver<<" "<<ANU<<" "<<vUPre.x()<<" "<<vpPostNorm.x()<<" "<< vUPre.y()<<" "<<vpPostNorm.y()<<" "<<sinpPostInv<<G4endl;
499#endif
500
501 transf[2][4] = -QAver*ANU*(vVPre.x()*vpPostNorm.x()+vVPre.y()*vpPostNorm.y()+vVPre.z()*vpPostNorm.z())*sinpPostInv;
502
503 // *** Yt
504
505 transf[3][0] = pAver*(vUPost.x()*deltaPos.x()+vUPost.y()*deltaPos.y() )
506 *(1.+deltaPInv*pAver);
507#ifdef G4EVERBOSE
508 if(iverbose >= 3) G4cout <<"ctransf30= " << transf[3][0] <<" "<< pAver<<" "<<vUPost.x()<<" "<<deltaPos.x()<<" "<<vUPost.y()<<" "<<deltaPos.y()
509 <<" "<<deltaPInv<<" "<<pAver<<G4endl;
510#endif
511
512 transf[3][1] = ( sinThetaAver*(vVPre.x()*vUPost.x()+vVPre.y()*vUPost.y() ) +
513 OMcosThetaAver*(vHVPre.x()*vUPost.x()+vHVPre.y()*vUPost.y() ) +
514 TMSINT*(vHAverNorm.x()*vUPost.x()+vHAverNorm.y()*vUPost.y() )*
515 (vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z()) )/QAver;
516
517 transf[3][2] = ( sinThetaAver*(vUPre.x()*vUPost.x()+vUPre.y()*vUPost.y() ) +
518 OMcosThetaAver*(vHUPre.x()*vUPost.x()+vHUPre.y()*vUPost.y() ) +
519 TMSINT*(vHAverNorm.x()*vUPost.x()+vHAverNorm.y()*vUPost.y() )*
520 (vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() ) )*sinpPre/QAver;
521#ifdef G4EVERBOSE
522 if(iverbose >= 3) G4cout <<"ctransf32= " << transf[3][2] <<" "<< sinThetaAver<<" "<<vUPre.x()<<" "<<vUPost.x()<<" "<<vUPre.y()<<" "<<vUPost.y() <<" "<<
523 OMcosThetaAver<<" "<<vHUPre.x()<<" "<<vUPost.x()<<" "<<vHUPre.y()<<" "<<vUPost.y() <<" "<<
524 TMSINT<<" "<<vHAverNorm.x()<<" "<<vUPost.x()<<" "<<vHAverNorm.y()<<" "<<vUPost.y() <<" "<<
525 vHAverNorm.x()<<" "<<vUPre.x()<<" "<<vHAverNorm.y()<<" "<<vUPre.y() <<" "<<sinpPre<<" "<<QAver<<G4endl;
526#endif
527
528 transf[3][3] = (vUPre.x()*vUPost.x()+vUPre.y()*vUPost.y() );
529
530 transf[3][4] = (vVPre.x()*vUPost.x()+vVPre.y()*vUPost.y() );
531
532 // *** Zt
533 transf[4][0] = pAver*(vVPost.x()*deltaPos.x()+vVPost.y()*deltaPos.y()+vVPost.z()*deltaPos.z())
534 *(1.+deltaPInv*pAver);
535
536 transf[4][1] = ( sinThetaAver*(vVPre.x()*vVPost.x()+vVPre.y()*vVPost.y()+vVPre.z()*vVPost.z()) +
537 OMcosThetaAver*(vHVPre.x()*vVPost.x()+vHVPre.y()*vVPost.y()+vHVPre.z()*vVPost.z()) +
538 TMSINT*(vHAverNorm.x()*vVPost.x()+vHAverNorm.y()*vVPost.y()+vHAverNorm.z()*vVPost.z())*
539 (vHAverNorm.x()*vVPre.x()+vHAverNorm.y()*vVPre.y()+vHAverNorm.z()*vVPre.z()) )/QAver;
540#ifdef G4EVERBOSE
541 if(iverbose >= 3)G4cout <<"ctransf41= " << transf[4][1] <<" "<< sinThetaAver<<" "<< OMcosThetaAver <<" "<<TMSINT<<" "<< vVPre <<" "<<vVPost <<" "<<vHVPre<<" "<<vHAverNorm <<" "<< QAver<<G4endl;
542#endif
543
544 transf[4][2] = ( sinThetaAver*(vUPre.x()*vVPost.x()+vUPre.y()*vVPost.y() ) +
545 OMcosThetaAver*(vHUPre.x()*vVPost.x()+vHUPre.y()*vVPost.y()+vHUPre.z()*vVPost.z()) +
546 TMSINT*(vHAverNorm.x()*vVPost.x()+vHAverNorm.y()*vVPost.y()+vHAverNorm.z()*vVPost.z())*
547 (vHAverNorm.x()*vUPre.x()+vHAverNorm.y()*vUPre.y() ) )*sinpPre/QAver;
548
549 transf[4][3] = (vUPre.x()*vVPost.x()+vUPre.y()*vVPost.y() );
550
551 transf[4][4] = (vVPre.x()*vVPost.x()+vVPre.y()*vVPost.y()+vVPre.z()*vVPost.z());
552 // if(iverbose >= 3) G4cout <<"ctransf44= " << transf[4][4] <<" "<< vVPre.x() <<" "<<vVPost.x() <<" "<< vVPre.y() <<" "<< vVPost.y() <<" "<< vVPre.z() <<" "<< vVPost.z() << G4endl;
553
554
555#ifdef G4EVERBOSE
556 if( iverbose >= 1 ) G4cout << "G4EP: transf matrix computed " << transf << G4endl;
557#endif
558 /* for( G4int ii=0;ii<5;ii++){
559 for( G4int jj=0;jj<5;jj++){
560 G4cout << transf[ii][jj] << " ";
561 }
562 G4cout << G4endl;
563 } */
564 }
565 }
566 // end of calculate transformation except it NEUTRAL PARTICLE OR FIELDFREE REGION
567 /* if( iverbose >= 1 ) G4cout << "G4EP: transf not updated but initialized " << theFirstStep << G4endl;
568 if( theFirstStep ) {
569 theTransfMat = transf;
570 theFirstStep = false;
571 }else{
572 theTransfMat = theTransfMat * transf;
573 if( iverbose >= 1 ) G4cout << "G4EP: transf matrix accumulated" << theTransfMat << G4endl;
574 }
575 */
576 theTransfMat = transf;
577#ifdef G4EVERBOSE
578 if( iverbose >= 1 ) G4cout << "G4EP: error matrix before transformation " << fError << G4endl;
579 if( iverbose >= 2 ) G4cout << " tf * err " << theTransfMat * fError << G4endl
580 << " transf matrix " << theTransfMat.T() << G4endl;
581#endif
582
583 fError = fError.similarity(theTransfMat).T();
584 //- fError = transf * fError * transf.T();
585#ifdef G4EVERBOSE
586 if( iverbose >= 1 ) G4cout << "G4EP: error matrix propagated " << fError << G4endl;
587#endif
588
589 //? S = B*S*BT S.similarity(B)
590 //? R = S
591 // not needed * *** TRANSFORM ERROR MATRIX FROM INTERNAL TO EXTERNAL VARIABLES;
592
593 PropagateErrorMSC( aTrack );
594
595 PropagateErrorIoni( aTrack );
596
597 return 0;
598}
const G4double kCarTolerance
@ G4ErrorMode_PropBackwards
@ G4ErrorStage_Deflation
Hep3Vector cross(const Hep3Vector &) const
G4double GetCharge() const
G4ErrorMatrix T() const
static G4ErrorPropagatorData * GetErrorPropagatorData()
G4ErrorSymMatrix similarity(const G4ErrorMatrix &m1) const
G4ErrorSymMatrix T() const
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
G4ThreeVector GetMomentum() const
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPreStepPoint() const
G4double GetStepLength() const
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const
const G4DynamicParticle * GetDynamicParticle() const
const G4Step * GetStep() const

Referenced by G4ErrorPropagator::MakeOneStep().

◆ SetMomentum()

virtual void G4ErrorFreeTrajState::SetMomentum ( const G4Vector3D mom)
inlinevirtual

Reimplemented from G4ErrorTrajState.

Definition at line 97 of file G4ErrorFreeTrajState.hh.

98 { SetParameters( fPosition, mom ); }
void SetParameters(const G4Point3D &pos, const G4Vector3D &mom)

◆ SetParameters()

void G4ErrorFreeTrajState::SetParameters ( const G4Point3D pos,
const G4Vector3D mom 
)
inline

Definition at line 100 of file G4ErrorFreeTrajState.hh.

101 {
102 fPosition = pos;
103 fMomentum = mom;
104 fTrajParam.SetParameters( pos, mom );
105 }
void SetParameters(const G4Point3D &pos, const G4Vector3D &mom)

Referenced by SetMomentum(), and SetPosition().

◆ SetPosition()

virtual void G4ErrorFreeTrajState::SetPosition ( const G4Point3D  pos)
inlinevirtual

Reimplemented from G4ErrorTrajState.

Definition at line 94 of file G4ErrorFreeTrajState.hh.

95 { SetParameters( pos, fMomentum ); }

◆ Update()

G4int G4ErrorFreeTrajState::Update ( const G4Track aTrack)
virtual

Reimplemented from G4ErrorTrajState.

Definition at line 175 of file G4ErrorFreeTrajState.cc.

176{
177 G4int ierr = 0;
178 fTrajParam.Update( aTrack );
179 UpdatePosMom( aTrack->GetPosition(), aTrack->GetMomentum() );
180 return ierr;
181}
int G4int
Definition: G4Types.hh:85
void Update(const G4Track *aTrack)
void UpdatePosMom(const G4Point3D &pos, const G4Vector3D &mom)

Referenced by G4ErrorPropagator::MakeOneStep().

Friends And Related Function Documentation

◆ operator<<

std::ostream & operator<< ( std::ostream &  out,
const G4ErrorFreeTrajState ts 
)
friend

Definition at line 185 of file G4ErrorFreeTrajState.cc.

186{
187 std::ios::fmtflags orig_flags = out.flags();
188
189 out.setf(std::ios::fixed,std::ios::floatfield);
190
191 ts.DumpPosMomError( out );
192
193 out << " G4ErrorFreeTrajState: Params: " << ts.fTrajParam << G4endl;
194
195 out.flags(orig_flags);
196
197 return out;
198}
void DumpPosMomError(std::ostream &out=G4cout) const

The documentation for this class was generated from the following files: