Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ErrorSurfaceTrajState.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: G4ErrorSurfaceTrajState.cc 69014 2013-04-15 09:42:51Z gcosmo $
27 //
28 // ------------------------------------------------------------
29 // GEANT 4 class implementation file
30 // ------------------------------------------------------------
31 //
32 
34 #include "G4ErrorPropagatorData.hh"
35 
36 #include "G4PhysicalConstants.hh"
37 #include "G4SystemOfUnits.hh"
38 #include "G4Field.hh"
39 #include "G4FieldManager.hh"
41 
42 #include "G4ErrorMatrix.hh"
43 
44 #include <iomanip>
45 
46 //------------------------------------------------------------------------
48 G4ErrorSurfaceTrajState( const G4String& partType, const G4Point3D& pos,
49  const G4Vector3D& mom, const G4Vector3D& vecU,
50  const G4Vector3D& vecV, const G4ErrorTrajErr& errmat)
51  : G4ErrorTrajState( partType, pos, mom, errmat )
52 {
53  Init();
54  fTrajParam = G4ErrorSurfaceTrajParam( pos, mom, vecU, vecV );
55 }
56 
57 
58 //------------------------------------------------------------------------
60 G4ErrorSurfaceTrajState( const G4String& partType, const G4Point3D& pos,
61  const G4Vector3D& mom, const G4Plane3D& plane,
62  const G4ErrorTrajErr& errmat )
63  : G4ErrorTrajState( partType, pos, mom, errmat )
64 {
65  Init();
66  fTrajParam = G4ErrorSurfaceTrajParam( pos, mom, plane );
67 
68 }
69 
70 
71 //------------------------------------------------------------------------
74  : G4ErrorTrajState( tpSC.GetParticleType(), tpSC.GetPosition(),
75  tpSC.GetMomentum() )
76 {
77  // fParticleType = tpSC.GetParticleType();
78  // fPosition = tpSC.GetPosition();
79  // fMomentum = tpSC.GetMomentum();
80  fTrajParam = G4ErrorSurfaceTrajParam( fPosition, fMomentum, plane );
81  Init();
82 
83  //----- Get the error matrix in SC coordinates
85 }
86 
87 
88 //------------------------------------------------------------------------
91  const G4Vector3D& vecV , G4ErrorMatrix &transfM)
92  : G4ErrorTrajState( tpSC.GetParticleType(), tpSC.GetPosition(),
93  tpSC.GetMomentum() )
94 {
95  Init(); // needed to define charge sign
96  fTrajParam = G4ErrorSurfaceTrajParam( fPosition, fMomentum, vecU, vecV );
97  //----- Get the error matrix in SC coordinates
98  transfM= BuildErrorMatrix( tpSC, vecU, vecV );
99 }
100 
101 
102 //------------------------------------------------------------------------
105  const G4Vector3D& )
106 {
107  G4double sclambda = tpSC.GetParameters().GetLambda();
108  G4double scphi = tpSC.GetParameters().GetPhi();
110  sclambda *= -1;
111  scphi += CLHEP::pi;
112  }
113  G4double cosLambda = std::cos( sclambda );
114  G4double sinLambda = std::sin( sclambda );
115  G4double sinPhi = std::sin( scphi );
116  G4double cosPhi = std::cos( scphi );
117 
118 #ifdef G4EVERBOSE
119  if( iverbose >= 4) G4cout << " PM " << fMomentum.mag() << " pLambda " << sclambda << " pPhi " << scphi << G4endl;
120 #endif
121 
122  G4ThreeVector vTN( cosLambda*cosPhi, cosLambda*sinPhi,sinLambda );
123  G4ThreeVector vUN( -sinPhi, cosPhi, 0. );
124  G4ThreeVector vVN( -vTN.z()*vUN.y(),vTN.z()*vUN.x(), cosLambda );
125 
126 #ifdef G4EVERBOSE
127  if( iverbose >= 4) G4cout << " SC2SD: vTN " << vTN << " vUN " << vUN << " vVN " << vVN << G4endl;
128 #endif
129  G4double UJ = vUN*GetVectorV();
130  G4double UK = vUN*GetVectorW();
131  G4double VJ = vVN*GetVectorV();
132  G4double VK = vVN*GetVectorW();
133 
134 
135  //--- Get transformation first
136  G4ErrorMatrix transfM(5, 5, 0 );
137  //--- Get magnetic field
139 
140  G4Vector3D vectorU = GetVectorV().cross( GetVectorW() );
141  G4double T1R = 1. / ( vTN * vectorU );
142 
143 #ifdef G4EVERBOSE
144  if( iverbose >= 4) G4cout << "surf vectors:U,V,W " << vectorU << " " << GetVectorV() << " " << GetVectorW() << " T1R:"<<T1R<<G4endl;
145 #endif
146 
147 
148  if( fCharge != 0 && field ) {
149  G4double pos[3]; pos[0] = fPosition.x()*cm; pos[1] = fPosition.y()*cm; pos[2] = fPosition.z()*cm;
150  G4double Hd[3];
151  field->GetFieldValue( pos, Hd );
152  G4ThreeVector H = G4ThreeVector( Hd[0], Hd[1], Hd[2] ) / tesla *10.; //in kilogauus
153  G4double magH = H.mag();
154  G4double invP = 1./(fMomentum.mag()/GeV);
155  G4double magHM = magH * invP;
156  if( magH != 0. ) {
157  G4double magHM2 = fCharge / magH;
158  G4double Q = -magHM * c_light/ (km/ns);
159 #ifdef G4EVERBOSE
160  if( iverbose >= 4) G4cout << GeV << " Q " << Q << " magHM " << magHM << " c_light/(km/ns) " << c_light/(km/ns) << G4endl;
161 #endif
162 
163  G4double sinz = -H*vUN * magHM2;
164  G4double cosz = H*vVN * magHM2;
165  G4double T3R = Q * std::pow(T1R,3);
166  G4double UI = vUN * vectorU;
167  G4double VI = vVN * vectorU;
168 #ifdef G4EVERBOSE
169  if( iverbose >= 4) {
170  G4cout << " T1R " << T1R << " T3R " << T3R << G4endl;
171  G4cout << " UI " << UI << " VI " << VI << " vectorU " << vectorU << G4endl;
172  G4cout << " UJ " << UJ << " VJ " << VJ << G4endl;
173  G4cout << " UK " << UK << " VK " << VK << G4endl;
174  }
175 #endif
176 
177  transfM[1][3] = -UI*( VK*cosz-UK*sinz)*T3R;
178  transfM[1][4] = -VI*( VK*cosz-UK*sinz)*T3R;
179  transfM[2][3] = UI*( VJ*cosz-UJ*sinz)*T3R;
180  transfM[2][4] = VI*( VJ*cosz-UJ*sinz)*T3R;
181  }
182  }
183 
184  G4double T2R = T1R * T1R;
185  transfM[0][0] = 1.;
186  transfM[1][1] = -UK*T2R;
187  transfM[1][2] = VK*cosLambda*T2R;
188  transfM[2][1] = UJ*T2R;
189  transfM[2][2] = -VJ*cosLambda*T2R;
190  transfM[3][3] = VK*T1R;
191  transfM[3][4] = -UK*T1R;
192  transfM[4][3] = -VJ*T1R;
193  transfM[4][4] = UJ*T1R;
194 
195 #ifdef G4EVERBOSE
196  if( iverbose >= 4) G4cout << " SC2SD transf matrix A= " << transfM << G4endl;
197 #endif
198  fError = G4ErrorTrajErr( tpSC.GetError().similarity( transfM ) );
199 
200 #ifdef G4EVERBOSE
201  if( iverbose >= 1) G4cout << "G4EP: error matrix SC2SD " << fError << G4endl;
202  if( iverbose >= 4) G4cout << "G4ErrorSurfaceTrajState from SC " << *this << G4endl;
203 #endif
204 
205  return transfM; // want to use trasnfM for the reverse transform
206 }
207 
208 
209 //------------------------------------------------------------------------
210 void G4ErrorSurfaceTrajState::Init()
211 {
213  BuildCharge();
214 }
215 
216 
217 //------------------------------------------------------------------------
218 void G4ErrorSurfaceTrajState::Dump( std::ostream& out ) const
219 {
220  out << *this;
221 }
222 
223 
224 //------------------------------------------------------------------------
225 std::ostream& operator<<(std::ostream& out, const G4ErrorSurfaceTrajState& ts)
226 {
227  std::ios::fmtflags oldFlags = out.flags();
228  out.setf(std::ios::fixed,std::ios::floatfield);
229 
230  ts.DumpPosMomError( out );
231 
232  out << " G4ErrorSurfaceTrajState: Params: " << ts.fTrajParam << G4endl;
233  out.flags(oldFlags);
234  return out;
235 }
G4ErrorSymMatrix G4ErrorTrajErr
G4ErrorSymMatrix similarity(const G4ErrorMatrix &m1) const
CLHEP::Hep3Vector G4ThreeVector
virtual void GetFieldValue(const double Point[4], double *fieldArr) const =0
double x() const
G4double GetPhi() const
G4ErrorMatrix BuildErrorMatrix(G4ErrorFreeTrajState &tpSC, const G4Vector3D &vecV, const G4Vector3D &vecW)
G4double GetLambda() const
G4ErrorTrajErr fError
virtual void Dump(std::ostream &out=G4cout) const
double z() const
G4GLOB_DLL std::ostream G4cout
G4ErrorSurfaceTrajState(const G4String &partType, const G4Point3D &pos, const G4Vector3D &mom, const G4Plane3D &plane, const G4ErrorTrajErr &errmat=G4ErrorTrajErr(5, 0))
void DumpPosMomError(std::ostream &out=G4cout) const
static G4TransportationManager * GetTransportationManager()
G4FieldManager * GetFieldManager() const
double y() const
std::ostream & operator<<(std::ostream &, const BasicVector3D< float > &)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static G4ErrorPropagatorData * GetErrorPropagatorData()
const G4Field * GetDetectorField() const
double mag() const
#define ns
Definition: xmlparse.cc:597
BasicVector3D< T > cross(const BasicVector3D< T > &v) const
float c_light
Definition: hepunit.py:257
G4ErrorTrajErr GetError() const
G4ErrorFreeTrajParam GetParameters() const