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
G4double GetInvP () const

Friends

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


Detailed Description

Definition at line 52 of file G4ErrorSurfaceTrajParam.hh.


Constructor & Destructor Documentation

G4ErrorSurfaceTrajParam::G4ErrorSurfaceTrajParam (  )  [inline]

Definition at line 56 of file G4ErrorSurfaceTrajParam.hh.

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

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

Definition at line 41 of file G4ErrorSurfaceTrajParam.cc.

References SetParameters().

00043 {
00044   SetParameters( pos, mom, vecV, vecW );
00045 }

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

Definition at line 50 of file G4ErrorSurfaceTrajParam.cc.

References SetParameters().

00052 {
00053   SetParameters( pos, mom, plane );
00054 }

virtual G4ErrorSurfaceTrajParam::~G4ErrorSurfaceTrajParam (  )  [inline, virtual]

Definition at line 62 of file G4ErrorSurfaceTrajParam.hh.

00062 {}


Member Function Documentation

G4Vector3D G4ErrorSurfaceTrajParam::GetDirection (  )  const [inline]

Definition at line 74 of file G4ErrorSurfaceTrajParam.hh.

00074 { return fDir; }

G4double G4ErrorSurfaceTrajParam::GetInvP (  )  const [inline]

Definition at line 82 of file G4ErrorSurfaceTrajParam.hh.

00082 { return fInvP; }

G4Vector3D G4ErrorSurfaceTrajParam::GetPlaneNormal (  )  const [inline]

Definition at line 75 of file G4ErrorSurfaceTrajParam.hh.

00075 { return fVectorV.cross(fVectorW); }

G4double G4ErrorSurfaceTrajParam::GetPV (  )  const [inline]

Definition at line 78 of file G4ErrorSurfaceTrajParam.hh.

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

00078 { return fPV; }

G4double G4ErrorSurfaceTrajParam::GetPW (  )  const [inline]

Definition at line 79 of file G4ErrorSurfaceTrajParam.hh.

Referenced by G4ErrorFreeTrajState::G4ErrorFreeTrajState().

00079 { return fPW; }

G4double G4ErrorSurfaceTrajParam::GetV (  )  const [inline]

Definition at line 80 of file G4ErrorSurfaceTrajParam.hh.

00080 { return fV; }

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorV (  )  const [inline]

Definition at line 76 of file G4ErrorSurfaceTrajParam.hh.

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

00076 { return fVectorV; }

G4Vector3D G4ErrorSurfaceTrajParam::GetVectorW (  )  const [inline]

Definition at line 77 of file G4ErrorSurfaceTrajParam.hh.

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

00077 { return fVectorW; }

G4double G4ErrorSurfaceTrajParam::GetW (  )  const [inline]

Definition at line 81 of file G4ErrorSurfaceTrajParam.hh.

00081 { return fW; }

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

Definition at line 58 of file G4ErrorSurfaceTrajParam.cc.

References G4GeometryTolerance::GetInstance(), G4GeometryTolerance::GetSurfaceTolerance(), and SetParameters().

00060 {
00061   //--- Get two perpendicular vectors: first parallel X
00062   //    (unless normal is parallel to X, then take Y)
00063 
00064   G4double kCarTolerance =
00065     G4GeometryTolerance::GetInstance()->GetSurfaceTolerance();
00066 
00067   G4ThreeVector Xvec(1.,0.,0.);
00068   G4Vector3D vecV = -Xvec.cross(plane.normal());
00069   if( vecV.mag() < kCarTolerance )
00070   {
00071     G4ThreeVector Zvec(0.,0.,1.);
00072     vecV = Zvec.cross(plane.normal());
00073   }
00074 
00075   G4Vector3D vecW = plane.normal().cross( vecV );
00076 
00077   SetParameters( pos, mom, vecV, vecW );
00078 }

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

Definition at line 83 of file G4ErrorSurfaceTrajParam.cc.

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

00085 {
00086   if( mom.mag() > 0. ) {
00087     fDir = mom;
00088     fDir /= mom.mag();
00089   } else {
00090     fDir = G4Vector3D(0.,0.,0.);
00091   }
00092   fVectorV = vecV / vecV.mag();
00093   fVectorW = vecW / vecW.mag();
00094   fInvP = 1./mom.mag();
00095   G4ThreeVector momv(mom);
00096   //check 3 vectors are ortogonal and right handed
00097 
00098   // now all 4 scalar memeber variables retain the signs
00099   //  fPV = momv.project( vecV ).mag();
00100   //  fPW = momv.project( vecW ).mag();
00101   fPV = momv.dot( vecV );
00102   fPW = momv.dot( vecW );
00103 
00104 
00105   G4ThreeVector posv(pos);
00106   // fV = posv.project( vecV ).mag();
00107   // fW = posv.project( vecW ).mag();
00108   fV = posv.dot( vecV );
00109   fW = posv.dot( vecW );
00110 }


Friends And Related Function Documentation

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

Definition at line 114 of file G4ErrorSurfaceTrajParam.cc.

00115 {
00116   //  long mode = out.setf(std::ios::fixed,std::ios::floatfield);
00117   
00118   //  out << tp.theType;
00119   //  out << std::setprecision(5) << std::setw(10);
00120   out << " InvP= " << tp.fInvP << " PV= " << tp.fPV
00121       << " PW= " << tp.fPW << " V= " << tp.fV << " W= " << tp.fW << G4endl;
00122   out << " vectorV direction= " << tp.fVectorV
00123       << " vectorW direction= " << tp.fVectorW << G4endl;
00124     
00125   return out;
00126 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:55 2013 for Geant4 by  doxygen 1.4.7