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

#include <G4ErrorFreeTrajParam.hh>

Public Member Functions

 G4ErrorFreeTrajParam ()
 
 G4ErrorFreeTrajParam (const G4Point3D &pos, const G4Vector3D &mom)
 
virtual ~G4ErrorFreeTrajParam ()
 
void Update (const G4Track *aTrack)
 
void SetParameters (const G4Point3D &pos, const G4Vector3D &mom)
 
G4Vector3D GetDirection () const
 
G4double GetInvP () const
 
G4double GetLambda () const
 
G4double GetPhi () const
 
G4double GetYPerp () const
 
G4double GetZPerp () const
 

Friends

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

Detailed Description

Definition at line 49 of file G4ErrorFreeTrajParam.hh.

Constructor & Destructor Documentation

◆ G4ErrorFreeTrajParam() [1/2]

G4ErrorFreeTrajParam::G4ErrorFreeTrajParam ( )
inline

Definition at line 53 of file G4ErrorFreeTrajParam.hh.

54 : fInvP(0.), fLambda(0.), fPhi(0.), fYPerp(0.), fZPerp(0.){}

◆ G4ErrorFreeTrajParam() [2/2]

G4ErrorFreeTrajParam::G4ErrorFreeTrajParam ( const G4Point3D pos,
const G4Vector3D mom 
)

Definition at line 40 of file G4ErrorFreeTrajParam.cc.

42{
43 SetParameters( pos, mom );
44}
void SetParameters(const G4Point3D &pos, const G4Vector3D &mom)

◆ ~G4ErrorFreeTrajParam()

virtual G4ErrorFreeTrajParam::~G4ErrorFreeTrajParam ( )
inlinevirtual

Definition at line 58 of file G4ErrorFreeTrajParam.hh.

58{}

Member Function Documentation

◆ GetDirection()

G4Vector3D G4ErrorFreeTrajParam::GetDirection ( ) const
inline

Definition at line 70 of file G4ErrorFreeTrajParam.hh.

70{ return fDir;}

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetInvP()

G4double G4ErrorFreeTrajParam::GetInvP ( ) const
inline

Definition at line 72 of file G4ErrorFreeTrajParam.hh.

72{ return fInvP; }

◆ GetLambda()

G4double G4ErrorFreeTrajParam::GetLambda ( ) const
inline

Definition at line 73 of file G4ErrorFreeTrajParam.hh.

73{ return fLambda; }

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

◆ GetPhi()

G4double G4ErrorFreeTrajParam::GetPhi ( ) const
inline

◆ GetYPerp()

G4double G4ErrorFreeTrajParam::GetYPerp ( ) const
inline

Definition at line 75 of file G4ErrorFreeTrajParam.hh.

75{ return fYPerp; }

◆ GetZPerp()

G4double G4ErrorFreeTrajParam::GetZPerp ( ) const
inline

Definition at line 76 of file G4ErrorFreeTrajParam.hh.

76{ return fZPerp; }

◆ SetParameters()

void G4ErrorFreeTrajParam::SetParameters ( const G4Point3D pos,
const G4Vector3D mom 
)

Definition at line 48 of file G4ErrorFreeTrajParam.cc.

50{
51 fDir = mom;
52 fInvP = 1./mom.mag();
53 fLambda = 90.*deg - mom.theta();
54 fPhi = mom.phi();
55 G4Vector3D vxPerp(0.,0.,0.);
56 if( mom.mag() > 0.) {
57 vxPerp = mom/mom.mag();
58 }
59 G4Vector3D vyPerp = G4Vector3D( -vxPerp.y(), vxPerp.x(), 0.);
60 vyPerp /= vyPerp.mag();
61 G4Vector3D vzPerp = vxPerp.cross( vyPerp );
62 vzPerp /= vzPerp.mag();
63 // check if right handed
64 // fXPerp = pos.proj( mom );
65 G4ThreeVector posv(pos);
66 if( vyPerp.mag() != 0. ) {
67 fYPerp = posv.project( vyPerp ).mag();
68 fZPerp = posv.project( vzPerp ).mag();
69 } else {
70 fYPerp = 0.;
71 fZPerp = 0.;
72 }
73}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

Referenced by G4ErrorFreeTrajParam(), G4ErrorFreeTrajState::SetParameters(), and Update().

◆ Update()

void G4ErrorFreeTrajParam::Update ( const G4Track aTrack)

Definition at line 76 of file G4ErrorFreeTrajParam.cc.

77{
78 SetParameters( aTrack->GetPosition(), aTrack->GetMomentum() );
79
80}
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const

Referenced by G4ErrorFreeTrajState::Update().

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 84 of file G4ErrorFreeTrajParam.cc.

85{
86 G4int oldprc = out.precision(8);
87 out << " InvP= " << tp.fInvP << " Theta= "
88 << tp.fLambda << " Phi= " << tp.fPhi << " YPerp= " << tp.fYPerp
89 << " ZPerp= " << tp.fZPerp << G4endl;
90 out << " momentum direction= " << tp.fDir << G4endl;
91 out.precision(oldprc);
92
93 return out;
94}
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52

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