Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ErrorFreeTrajParam.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: G4ErrorFreeTrajParam.cc 69014 2013-04-15 09:42:51Z gcosmo $
27 //
28 // ------------------------------------------------------------
29 // GEANT 4 class implementation file
30 // ------------------------------------------------------------
31 //
32 
33 #include "G4ErrorFreeTrajParam.hh"
34 #include "G4ThreeVector.hh"
35 #include "G4SystemOfUnits.hh"
36 
37 #include <iomanip>
38 
39 //------------------------------------------------------------------------
41  const G4Vector3D& mom )
42 {
43  SetParameters( pos, mom );
44 }
45 
46 
47 //------------------------------------------------------------------------
49  const G4Vector3D& mom )
50 {
51  fInvP = 1./mom.mag();
52  fDir = mom*fInvP;
53  fLambda = 90.*deg - mom.theta();
54  fPhi = mom.phi();
55  G4Vector3D vxPerp(0.,0.,0.);
56  if( mom.mag() > 0.) {
57  vxPerp = mom/mom.mag();
58  }
59  G4Vector3D vyPerp = G4Vector3D( -vxPerp.y(), vxPerp.x(), 0.);
60  vyPerp /= vyPerp.mag();
61  G4Vector3D vzPerp = vxPerp.cross( vyPerp );
62  vzPerp /= vzPerp.mag();
63  // check if right handed
64  // fXPerp = pos.proj( mom );
65  G4ThreeVector posv(pos);
66  if( vyPerp.mag() != 0. ) {
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  } else {
73  fYPerp = 0.;
74  fZPerp = 0.;
75  }
76 }
77 
78 //------------------------------------------------------------------------
80 {
81  SetParameters( aTrack->GetPosition(), aTrack->GetMomentum() );
82 
83 }
84 
85 
86 //------------------------------------------------------------------------
87 std::ostream& operator<<(std::ostream& out, const G4ErrorFreeTrajParam& tp)
88 {
89  G4int oldprc = out.precision(8);
90  out << " InvP= " << tp.fInvP << " Theta= "
91  << tp.fLambda << " 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 }
double dot(const Hep3Vector &) const
const G4ThreeVector & GetPosition() const
HepGeom::Vector3D< G4double > G4Vector3D
Definition: G4Vector3D.hh:35
void SetParameters(const G4Point3D &pos, const G4Vector3D &mom)
int G4int
Definition: G4Types.hh:78
void Update(const G4Track *aTrack)
G4ThreeVector GetMomentum() const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
#define G4endl
Definition: G4ios.hh:61
static constexpr double deg
Definition: G4SIunits.hh:152
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
static const G4double pos