Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4ErrorSurfaceTrajParam.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// $Id$
27//
28// ------------------------------------------------------------
29// GEANT 4 class implementation file
30// ------------------------------------------------------------
31//
32
34#include <iomanip>
35
36#include "G4ThreeVector.hh"
38
39//------------------------------------------------------------------------
41G4ErrorSurfaceTrajParam( const G4Point3D& pos, const G4Vector3D& mom,
42 const G4Vector3D& vecV, const G4Vector3D& vecW )
43{
44 SetParameters( pos, mom, vecV, vecW );
45}
46
47
48//------------------------------------------------------------------------
50G4ErrorSurfaceTrajParam( const G4Point3D& pos, const G4Vector3D& mom,
51 const G4Plane3D& plane )
52{
53 SetParameters( pos, mom, plane );
54}
55
56//------------------------------------------------------------------------
58SetParameters( const G4Point3D& pos, const G4Vector3D& mom,
59 const G4Plane3D& plane )
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}
79
80
81//------------------------------------------------------------------------
83SetParameters( const G4Point3D& pos, const G4Vector3D& mom,
84 const G4Vector3D& vecV, const G4Vector3D& vecW )
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}
105
106
107//------------------------------------------------------------------------
108std::ostream& operator<<(std::ostream& out, const G4ErrorSurfaceTrajParam& tp)
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}
std::ostream & operator<<(std::ostream &out, const G4ErrorSurfaceTrajParam &tp)
double G4double
Definition: G4Types.hh:64
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
#define G4endl
Definition: G4ios.hh:52
Hep3Vector cross(const Hep3Vector &) const
double mag() const
Hep3Vector project() const
void SetParameters(const G4Point3D &pos, const G4Vector3D &mom, const G4Vector3D &vecV, const G4Vector3D &vecW)
G4double GetSurfaceTolerance() const
static G4GeometryTolerance * GetInstance()
Normal3D< T > normal() const
Definition: Plane3D.h:90