|
||||
File indexing completed on 2025-01-18 09:58:14
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 // 0027 // 0028 // Class Description: 0029 // 0030 // Represents a free G4ErrorTrajState 0031 // It can be represented by the 5 variables 0032 // 1/p, lambda, phi, y_perp, z_perp 0033 // where lambda and phi are the dip and azimuthal angles related 0034 // to the momentum components in the following way: 0035 // p_x = p cos(lambda) cos(phi) ! lambda = 90 - theta 0036 // p_y = p cos(lambda) sin(phi) 0037 // p_z = p sin(lambda) 0038 // y_perp and z_perp are the coordinates of the trajectory in a 0039 // local orthonormal reference frame with the x_perp axis along the 0040 // particle direction, the y_perp being parallel to the x-y plane. 0041 // 0042 // This class also takes care of propagating the error associated to 0043 // the trajectory 0044 0045 // History: 0046 // - Created: P. Arce 0047 // -------------------------------------------------------------------- 0048 0049 #ifndef G4ErrorFreeTrajState_hh 0050 #define G4ErrorFreeTrajState_hh 0051 0052 #include "globals.hh" 0053 0054 #include "G4ErrorMatrix.hh" 0055 0056 #include "G4ErrorTrajState.hh" 0057 #include "G4ErrorFreeTrajParam.hh" 0058 0059 #include "G4Point3D.hh" 0060 #include "G4Vector3D.hh" 0061 0062 class G4ErrorSurfaceTrajState; 0063 0064 class G4ErrorFreeTrajState : public G4ErrorTrajState 0065 { 0066 public: // with description 0067 G4ErrorFreeTrajState() 0068 : theFirstStep(true) 0069 {} 0070 G4ErrorFreeTrajState(const G4String& partName, const G4Point3D& pos, 0071 const G4Vector3D& mom, 0072 const G4ErrorTrajErr& errmat = G4ErrorTrajErr(5, 0)); 0073 // Constructor by providing particle, position and momentum 0074 0075 G4ErrorFreeTrajState(const G4ErrorSurfaceTrajState& tpOS); 0076 // Constructor by providing G4ErrorSurfaceTrajState 0077 0078 ~G4ErrorFreeTrajState() {} 0079 0080 virtual G4int Update(const G4Track* aTrack); 0081 // update parameters from G4Track 0082 0083 virtual G4int PropagateError(const G4Track* aTrack); 0084 // propagate the error along the step 0085 0086 virtual void Dump(std::ostream& out = G4cout) const; 0087 // dump TrajState parameters 0088 0089 friend std::ostream& operator<<(std::ostream&, 0090 const G4ErrorFreeTrajState& ts); 0091 0092 // Set and Get methods 0093 0094 virtual void SetPosition(const G4Point3D pos) 0095 { 0096 SetParameters(pos, fMomentum); 0097 } 0098 0099 virtual void SetMomentum(const G4Vector3D& mom) 0100 { 0101 SetParameters(fPosition, mom); 0102 } 0103 0104 void SetParameters(const G4Point3D& pos, const G4Vector3D& mom) 0105 { 0106 fPosition = pos; 0107 fMomentum = mom; 0108 fTrajParam.SetParameters(pos, mom); 0109 } 0110 0111 G4ErrorFreeTrajParam GetParameters() const { return fTrajParam; } 0112 0113 G4ErrorMatrix GetTransfMat() const { return theTransfMat; } 0114 0115 private: 0116 void Init(); 0117 // define TrajState type and build charge 0118 0119 G4int PropagateErrorMSC(const G4Track* aTrack); 0120 // add the error associated to multiple scattering 0121 0122 void CalculateEffectiveZandA(const G4Material* mate, double& effZ, 0123 double& effA); 0124 // calculate effective Z and A (needed by PropagateErrorMSC) 0125 0126 G4int PropagateErrorIoni(const G4Track* aTrack); 0127 // add the error associated to ionization energy loss 0128 0129 private: 0130 G4ErrorFreeTrajParam fTrajParam; 0131 0132 G4ErrorMatrix theTransfMat; 0133 0134 G4bool theFirstStep; // to count if transf mat is updated or initialized 0135 }; 0136 0137 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |