Geant4 9.6.0
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)
 
virtual ~G4ErrorSurfaceTrajParam ()
 
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
 

Friends

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

Detailed Description

Definition at line 52 of file G4ErrorSurfaceTrajParam.hh.

Constructor & Destructor Documentation

◆ G4ErrorSurfaceTrajParam() [1/3]

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam ( )
inline

Definition at line 56 of file G4ErrorSurfaceTrajParam.hh.

57 : fInvP(0.), fPV(0.), fPW(0.), fV(0.), fW(0.) {}

◆ G4ErrorSurfaceTrajParam() [2/3]

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

Definition at line 40 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/3]

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

Definition at line 49 of file G4ErrorSurfaceTrajParam.cc.

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

◆ ~G4ErrorSurfaceTrajParam()

virtual G4ErrorSurfaceTrajParam::~G4ErrorSurfaceTrajParam ( )
inlinevirtual

Definition at line 62 of file G4ErrorSurfaceTrajParam.hh.

62{}

Member Function Documentation

◆ GetDirection()

G4Vector3D G4ErrorSurfaceTrajParam::GetDirection ( ) const
inline

Definition at line 74 of file G4ErrorSurfaceTrajParam.hh.

74{ return fDir; }

◆ GetPlaneNormal()

G4Vector3D G4ErrorSurfaceTrajParam::GetPlaneNormal ( ) const
inline

Definition at line 75 of file G4ErrorSurfaceTrajParam.hh.

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

◆ GetPV()

G4double G4ErrorSurfaceTrajParam::GetPV ( ) const
inline

Definition at line 78 of file G4ErrorSurfaceTrajParam.hh.

78{ return fPV; }

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetPW()

G4double G4ErrorSurfaceTrajParam::GetPW ( ) const
inline

Definition at line 79 of file G4ErrorSurfaceTrajParam.hh.

79{ return fPW; }

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

◆ GetV()

G4double G4ErrorSurfaceTrajParam::GetV ( ) const
inline

Definition at line 80 of file G4ErrorSurfaceTrajParam.hh.

80{ return fV; }

◆ GetVectorV()

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorV ( ) const
inline

Definition at line 76 of file G4ErrorSurfaceTrajParam.hh.

76{ return fVectorV; }

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

◆ GetVectorW()

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorW ( ) const
inline

Definition at line 77 of file G4ErrorSurfaceTrajParam.hh.

77{ return fVectorW; }

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

◆ GetW()

G4double G4ErrorSurfaceTrajParam::GetW ( ) const
inline

Definition at line 81 of file G4ErrorSurfaceTrajParam.hh.

81{ return fW; }

◆ SetParameters() [1/2]

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

Definition at line 57 of file G4ErrorSurfaceTrajParam.cc.

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

◆ SetParameters() [2/2]

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

Definition at line 82 of file G4ErrorSurfaceTrajParam.cc.

85{
86 if( mom.mag() > 0. ) {
87 fDir = mom;
88 fDir /= mom.mag();
89 } else {
90 fDir = G4Vector3D(0.,0.,0.);
91 }
92 fVectorV = vecV / vecV.mag();
93 fVectorW = vecW / vecW.mag();
94 fInvP = 1./mom.mag();
95 G4ThreeVector momv(mom);
96 //check 3 vectors are ortogonal and right handed
97
98 fPV = momv.project( vecV ).mag();
99 fPW = momv.project( vecW ).mag();
100
101 G4ThreeVector posv(pos);
102 fV = posv.project( vecV ).mag();
103 fW = posv.project( vecW ).mag();
104}
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35

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 108 of file G4ErrorSurfaceTrajParam.cc.

109{
110 // long mode = out.setf(std::ios::fixed,std::ios::floatfield);
111
112 // out << tp.theType;
113 // out << std::setprecision(5) << std::setw(10);
114 out << " InvP= " << tp.fInvP << " PV= " << tp.fPV
115 << " PW= " << tp.fPW << " V= " << tp.fV << " W= " << tp.fW << G4endl;
116 out << " vectorV direction= " << tp.fVectorV
117 << " vectorW direction= " << tp.fVectorW << G4endl;
118
119 return out;
120}
#define G4endl
Definition: G4ios.hh:52

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