Geant4-11
Public Member Functions | Private Attributes | Friends
G4ErrorFreeTrajParam Class Reference

#include <G4ErrorFreeTrajParam.hh>

Public Member Functions

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

Private Attributes

G4Vector3D fDir
 
G4double fInvP
 
G4double fLambda
 
G4double fPhi
 
G4double fYPerp
 
G4double fZPerp
 

Friends

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

Detailed Description

Definition at line 48 of file G4ErrorFreeTrajParam.hh.

Constructor & Destructor Documentation

◆ G4ErrorFreeTrajParam() [1/4]

G4ErrorFreeTrajParam::G4ErrorFreeTrajParam ( )
inline

Definition at line 51 of file G4ErrorFreeTrajParam.hh.

◆ G4ErrorFreeTrajParam() [2/4]

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

Definition at line 39 of file G4ErrorFreeTrajParam.cc.

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

References pos, and SetParameters().

◆ G4ErrorFreeTrajParam() [3/4]

G4ErrorFreeTrajParam::G4ErrorFreeTrajParam ( const G4ErrorFreeTrajParam )
default

◆ G4ErrorFreeTrajParam() [4/4]

G4ErrorFreeTrajParam::G4ErrorFreeTrajParam ( G4ErrorFreeTrajParam &&  )
default

◆ ~G4ErrorFreeTrajParam()

virtual G4ErrorFreeTrajParam::~G4ErrorFreeTrajParam ( )
inlinevirtual

Definition at line 65 of file G4ErrorFreeTrajParam.hh.

65{}

Member Function Documentation

◆ GetDirection()

G4Vector3D G4ErrorFreeTrajParam::GetDirection ( ) const
inline

Definition at line 81 of file G4ErrorFreeTrajParam.hh.

81{ return fDir; }

References fDir.

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetInvP()

G4double G4ErrorFreeTrajParam::GetInvP ( ) const
inline

Definition at line 83 of file G4ErrorFreeTrajParam.hh.

83{ return fInvP; }

References fInvP.

◆ GetLambda()

G4double G4ErrorFreeTrajParam::GetLambda ( ) const
inline

◆ GetPhi()

G4double G4ErrorFreeTrajParam::GetPhi ( ) const
inline

Definition at line 85 of file G4ErrorFreeTrajParam.hh.

85{ return fPhi; }

References fPhi.

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

◆ GetYPerp()

G4double G4ErrorFreeTrajParam::GetYPerp ( ) const
inline

Definition at line 86 of file G4ErrorFreeTrajParam.hh.

86{ return fYPerp; }

References fYPerp.

◆ GetZPerp()

G4double G4ErrorFreeTrajParam::GetZPerp ( ) const
inline

Definition at line 87 of file G4ErrorFreeTrajParam.hh.

87{ return fZPerp; }

References fZPerp.

◆ operator=() [1/2]

G4ErrorFreeTrajParam & G4ErrorFreeTrajParam::operator= ( const G4ErrorFreeTrajParam )
default

◆ operator=() [2/2]

G4ErrorFreeTrajParam & G4ErrorFreeTrajParam::operator= ( G4ErrorFreeTrajParam &&  )
default

◆ SetParameters()

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

Definition at line 46 of file G4ErrorFreeTrajParam.cc.

48{
49 fInvP = 1. / mom.mag();
50 fDir = mom * fInvP;
51 fLambda = 90. * deg - mom.theta();
52 fPhi = mom.phi();
53 G4Vector3D vxPerp(0., 0., 0.);
54 if(mom.mag() > 0.)
55 {
56 vxPerp = mom / mom.mag();
57 }
58 G4Vector3D vyPerp = G4Vector3D(-vxPerp.y(), vxPerp.x(), 0.);
59 vyPerp /= vyPerp.mag();
60 G4Vector3D vzPerp = vxPerp.cross(vyPerp);
61 vzPerp /= vzPerp.mag();
62 // check if right handed
63 // fXPerp = pos.proj( mom );
64 G4ThreeVector posv(pos);
65 if(vyPerp.mag() != 0.)
66 {
67 // now all 2 scalar memeber variables retain the signs
68 // fYPerp = posv.project( vyPerp ).mag();
69 // fZPerp = posv.project( vzPerp ).mag();
70 fYPerp = posv.dot(vyPerp);
71 fZPerp = posv.dot(vzPerp);
72 }
73 else
74 {
75 fYPerp = 0.;
76 fZPerp = 0.;
77 }
78}
static constexpr double deg
Definition: G4SIunits.hh:132
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

References HepGeom::BasicVector3D< T >::cross(), deg, CLHEP::Hep3Vector::dot(), fDir, fInvP, fLambda, fPhi, fYPerp, fZPerp, HepGeom::BasicVector3D< T >::mag(), HepGeom::BasicVector3D< T >::phi(), pos, HepGeom::BasicVector3D< T >::theta(), HepGeom::BasicVector3D< T >::x(), and HepGeom::BasicVector3D< T >::y().

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

◆ Update()

void G4ErrorFreeTrajParam::Update ( const G4Track aTrack)

Definition at line 81 of file G4ErrorFreeTrajParam.cc.

82{
83 SetParameters(aTrack->GetPosition(), aTrack->GetMomentum());
84}
const G4ThreeVector & GetPosition() const
G4ThreeVector GetMomentum() const

References G4Track::GetMomentum(), G4Track::GetPosition(), and SetParameters().

Referenced by G4ErrorFreeTrajState::Update().

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 87 of file G4ErrorFreeTrajParam.cc.

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

Field Documentation

◆ fDir

G4Vector3D G4ErrorFreeTrajParam::fDir
private

Definition at line 90 of file G4ErrorFreeTrajParam.hh.

Referenced by GetDirection(), and SetParameters().

◆ fInvP

G4double G4ErrorFreeTrajParam::fInvP
private

Definition at line 91 of file G4ErrorFreeTrajParam.hh.

Referenced by GetInvP(), and SetParameters().

◆ fLambda

G4double G4ErrorFreeTrajParam::fLambda
private

Definition at line 92 of file G4ErrorFreeTrajParam.hh.

Referenced by GetLambda(), and SetParameters().

◆ fPhi

G4double G4ErrorFreeTrajParam::fPhi
private

Definition at line 93 of file G4ErrorFreeTrajParam.hh.

Referenced by GetPhi(), and SetParameters().

◆ fYPerp

G4double G4ErrorFreeTrajParam::fYPerp
private

Definition at line 94 of file G4ErrorFreeTrajParam.hh.

Referenced by GetYPerp(), and SetParameters().

◆ fZPerp

G4double G4ErrorFreeTrajParam::fZPerp
private

Definition at line 95 of file G4ErrorFreeTrajParam.hh.

Referenced by GetZPerp(), and SetParameters().


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