source: trunk/source/error_propagation/include/G4ErrorFreeTrajState.hh@ 941

Last change on this file since 941 was 850, checked in by garnier, 17 years ago

geant4.8.2 beta

File size: 4.7 KB
RevLine 
[815]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//
27// $Id: G4ErrorFreeTrajState.hh,v 1.4 2007/05/31 20:27:06 arce Exp $
[850]28// GEANT4 tag $Name: HEAD $
[815]29//
30// Class Description:
31//
32// Represents a free G4ErrorTrajState
33// It can be represented by the 5 variables
34// 1/p, lambda, phi, y_perp, z_perp
35// where lambda and phi are the dip and azimuthal angles related
36// to the momentum components in the following way:
37// p_x = p cos(lambda) cos(phi) ! lambda = 90 - theta
38// p_y = p cos(lambda) sin(phi)
39// p_z = p sin(lambda)
40// y_perp and z_perp are the coordinates of the trajectory in a
41// local orthonormal reference frame with the x_perp axis along the
42// particle direction, the y_perp being parallel to the x-y plane.
43//
44// This class also takes care of propagating the error associated to
45// the trajectory
46
47// History:
48// - Created: P. Arce
49// --------------------------------------------------------------------
50
51#ifndef G4ErrorFreeTrajState_hh
52#define G4ErrorFreeTrajState_hh
53
54#include "globals.hh"
55
56#include "G4ErrorMatrix.hh"
57
58#include "G4ErrorTrajState.hh"
59#include "G4ErrorFreeTrajParam.hh"
60
61#include "G4Point3D.hh"
62#include "G4Vector3D.hh"
63
64class G4ErrorSurfaceTrajState;
65
66class G4ErrorFreeTrajState : public G4ErrorTrajState
67{
68 public: // with description
69
70 G4ErrorFreeTrajState(){}
71 G4ErrorFreeTrajState( const G4String& partType,
72 const G4Point3D& pos,
73 const G4Vector3D& mom,
74 const G4ErrorTrajErr& errmat = G4ErrorTrajErr(5,0) );
75 // Constructor by providing particle, position and momentum
76
77 G4ErrorFreeTrajState( const G4ErrorSurfaceTrajState& tpOS );
78 // Constructor by providing G4ErrorSurfaceTrajState
79
80 ~G4ErrorFreeTrajState(){}
81
82 virtual G4int Update( const G4Track* aTrack );
83 // update parameters from G4Track
84
85 virtual G4int PropagateError( const G4Track* aTrack );
86 // propagate the error along the step
87
88 virtual void Dump( std::ostream& out = G4cout ) const;
89 // dump TrajState parameters
90
91 friend
92 std::ostream& operator<<(std::ostream&, const G4ErrorFreeTrajState& ts);
93
94 // Set and Get methods
95
96 virtual void SetPosition( const G4Point3D pos )
97 { SetParameters( pos, fMomentum ); }
98
99 virtual void SetMomentum( const G4Vector3D& mom )
100 { SetParameters( fPosition, mom ); }
101
102 void SetParameters( const G4Point3D& pos, const G4Vector3D& mom )
103 {
104 fPosition = pos;
105 fMomentum = mom;
106 fTrajParam.SetParameters( pos, mom );
107 }
108
109 G4ErrorFreeTrajParam GetParameters() const
110 { return fTrajParam; }
111
112 private:
113
114 void Init();
115 // define TrajState type and build charge
116
117 G4int PropagateErrorMSC( const G4Track* aTrack );
118 // add the error associated to multiple scattering
119
120 void CalculateEffectiveZandA( const G4Material* mate, double& effZ, double& effA );
121 // calculate effective Z and A (needed by PropagateErrorMSC)
122
123 G4int PropagateErrorIoni( const G4Track* aTrack );
124 // add the error associated to ionization energy loss
125
126
127 private:
128
129 G4ErrorFreeTrajParam fTrajParam;
130
131 G4ErrorMatrix theTransfMat;
132
133 G4bool theFirstStep; // to count if transf mat is updated or initialized
134};
135
136#endif
Note: See TracBrowser for help on using the repository browser.