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

#include <G4ErrorSurfaceTrajParam.hh>

Public Member Functions

 G4ErrorSurfaceTrajParam ()
 
 G4ErrorSurfaceTrajParam (const G4Point3D &pos, const G4Vector3D &mom, const G4Vector3D &vecV, const G4Vector3D &vecW)
 
 G4ErrorSurfaceTrajParam (const G4Point3D &pos, const G4Vector3D &mom, const G4Plane3D &plane)
 
 G4ErrorSurfaceTrajParam (const G4ErrorSurfaceTrajParam &)=default
 
 G4ErrorSurfaceTrajParam (G4ErrorSurfaceTrajParam &&)=default
 
virtual ~G4ErrorSurfaceTrajParam ()
 
G4ErrorSurfaceTrajParamoperator= (const G4ErrorSurfaceTrajParam &)=default
 
G4ErrorSurfaceTrajParamoperator= (G4ErrorSurfaceTrajParam &&)=default
 
void SetParameters (const G4Point3D &pos, const G4Vector3D &mom, const G4Vector3D &vecV, const G4Vector3D &vecW)
 
void SetParameters (const G4Point3D &pos, const G4Vector3D &mom, const G4Plane3D &plane)
 
G4Vector3D GetDirection () const
 
G4Vector3D GetPlaneNormal () const
 
G4Vector3D GetVectorV () const
 
G4Vector3D GetVectorW () const
 
G4double GetPV () const
 
G4double GetPW () const
 
G4double GetV () const
 
G4double GetW () const
 
G4double GetInvP () const
 

Friends

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

Detailed Description

Definition at line 51 of file G4ErrorSurfaceTrajParam.hh.

Constructor & Destructor Documentation

◆ G4ErrorSurfaceTrajParam() [1/5]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( )
inline

Definition at line 54 of file G4ErrorSurfaceTrajParam.hh.

55 : fInvP(0.)
56 , fPV(0.)
57 , fPW(0.)
58 , fV(0.)
59 , fW(0.)
60 {}

◆ G4ErrorSurfaceTrajParam() [2/5]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( const G4Point3D pos,
const G4Vector3D mom,
const G4Vector3D vecV,
const G4Vector3D vecW 
)

Definition at line 39 of file G4ErrorSurfaceTrajParam.cc.

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

◆ G4ErrorSurfaceTrajParam() [3/5]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( const G4Point3D pos,
const G4Vector3D mom,
const G4Plane3D plane 
)

Definition at line 48 of file G4ErrorSurfaceTrajParam.cc.

51{
52 SetParameters(pos, mom, plane);
53}

◆ G4ErrorSurfaceTrajParam() [4/5]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( const G4ErrorSurfaceTrajParam )
default

◆ G4ErrorSurfaceTrajParam() [5/5]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( G4ErrorSurfaceTrajParam &&  )
default

◆ ~G4ErrorSurfaceTrajParam()

virtual G4ErrorSurfaceTrajParam::~G4ErrorSurfaceTrajParam ( )
inlinevirtual

Definition at line 70 of file G4ErrorSurfaceTrajParam.hh.

70{}

Member Function Documentation

◆ GetDirection()

G4Vector3D G4ErrorSurfaceTrajParam::GetDirection ( ) const
inline

Definition at line 86 of file G4ErrorSurfaceTrajParam.hh.

86{ return fDir; }

◆ GetInvP()

G4double G4ErrorSurfaceTrajParam::GetInvP ( ) const
inline

Definition at line 94 of file G4ErrorSurfaceTrajParam.hh.

94{ return fInvP; }

◆ GetPlaneNormal()

G4Vector3D G4ErrorSurfaceTrajParam::GetPlaneNormal ( ) const
inline

Definition at line 87 of file G4ErrorSurfaceTrajParam.hh.

87{ return fVectorV.cross(fVectorW); }
BasicVector3D< T > cross(const BasicVector3D< T > &v) const

◆ GetPV()

G4double G4ErrorSurfaceTrajParam::GetPV ( ) const
inline

Definition at line 90 of file G4ErrorSurfaceTrajParam.hh.

90{ return fPV; }

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetPW()

G4double G4ErrorSurfaceTrajParam::GetPW ( ) const
inline

Definition at line 91 of file G4ErrorSurfaceTrajParam.hh.

91{ return fPW; }

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetV()

G4double G4ErrorSurfaceTrajParam::GetV ( ) const
inline

Definition at line 92 of file G4ErrorSurfaceTrajParam.hh.

92{ return fV; }

◆ GetVectorV()

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorV ( ) const
inline

Definition at line 88 of file G4ErrorSurfaceTrajParam.hh.

88{ return fVectorV; }

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

◆ GetVectorW()

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorW ( ) const
inline

Definition at line 89 of file G4ErrorSurfaceTrajParam.hh.

89{ return fVectorW; }

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

◆ GetW()

G4double G4ErrorSurfaceTrajParam::GetW ( ) const
inline

Definition at line 93 of file G4ErrorSurfaceTrajParam.hh.

93{ return fW; }

◆ operator=() [1/2]

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

◆ operator=() [2/2]

G4ErrorSurfaceTrajParam & G4ErrorSurfaceTrajParam::operator= ( G4ErrorSurfaceTrajParam &&  )
default

◆ SetParameters() [1/2]

void G4ErrorSurfaceTrajParam::SetParameters ( const G4Point3D pos,
const G4Vector3D mom,
const G4Plane3D plane 
)

Definition at line 56 of file G4ErrorSurfaceTrajParam.cc.

59{
60 //--- Get two perpendicular vectors: first parallel X
61 // (unless normal is parallel to X, then take Y)
62
65
66 G4ThreeVector Xvec(1., 0., 0.);
67 G4Vector3D vecV = -Xvec.cross(plane.normal());
68 if(vecV.mag() < kCarTolerance)
69 {
70 G4ThreeVector Zvec(0., 0., 1.);
71 vecV = Zvec.cross(plane.normal());
72 }
73
74 G4Vector3D vecW = plane.normal().cross(vecV);
75
76 SetParameters(pos, mom, vecV, vecW);
77}
const G4double kCarTolerance
double G4double
Definition: G4Types.hh:83
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Normal3D< T > normal() const
Definition: Plane3D.h:97

◆ SetParameters() [2/2]

void G4ErrorSurfaceTrajParam::SetParameters ( const G4Point3D pos,
const G4Vector3D mom,
const G4Vector3D vecV,
const G4Vector3D vecW 
)

Definition at line 80 of file G4ErrorSurfaceTrajParam.cc.

84{
85 if(mom.mag() > 0.)
86 {
87 fDir = mom;
88 fDir /= mom.mag();
89 }
90 else
91 {
92 fDir = G4Vector3D(0., 0., 0.);
93 }
94 fVectorV = vecV / vecV.mag();
95 fVectorW = vecW / vecW.mag();
96 fInvP = 1. / mom.mag();
97 G4ThreeVector momv(mom);
98 // check 3 vectors are ortogonal and right handed
99
100 // now all 4 scalar memeber variables retain the signs
101 // fPV = momv.project( vecV ).mag();
102 // fPW = momv.project( vecW ).mag();
103 fPV = momv.dot(vecV);
104 fPW = momv.dot(vecW);
105
106 G4ThreeVector posv(pos);
107 // fV = posv.project( vecV ).mag();
108 // fW = posv.project( vecW ).mag();
109 fV = posv.dot(vecV);
110 fW = posv.dot(vecW);
111}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:34

Referenced by G4ErrorSurfaceTrajParam(), SetParameters(), and G4ErrorSurfaceTrajState::SetParameters().

Friends And Related Function Documentation

◆ operator<<

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

Definition at line 114 of file G4ErrorSurfaceTrajParam.cc.

115{
116 // long mode = out.setf(std::ios::fixed,std::ios::floatfield);
117
118 // out << tp.theType;
119 // out << std::setprecision(5) << std::setw(10);
120 out << " InvP= " << tp.fInvP << " PV= " << tp.fPV << " PW= " << tp.fPW
121 << " V= " << tp.fV << " W= " << tp.fW << G4endl;
122 out << " vectorV direction= " << tp.fVectorV
123 << " vectorW direction= " << tp.fVectorW << G4endl;
124
125 return out;
126}
#define G4endl
Definition: G4ios.hh:57

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