File indexing completed on 2025-01-18 09:59:08
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035 #ifndef G4StepPoint_hh
0036 #define G4StepPoint_hh 1
0037
0038 #include <cmath> // Include from 'system'
0039 #include <CLHEP/Units/PhysicalConstants.h>
0040
0041 #include "globals.hh" // Include from 'global'
0042 #include "G4Allocator.hh" // Include from 'global'
0043 #include "G4ThreeVector.hh" // Include from 'geometry'
0044 #include "G4VPhysicalVolume.hh" // Include from 'geometry'
0045 #include "G4SteppingControl.hh"
0046 #include "G4StepStatus.hh" // Include from 'track'
0047 #include "G4TouchableHandle.hh" // Include from 'geometry'
0048 #include "G4Material.hh"
0049 #include "G4LogicalVolume.hh"
0050
0051 class G4VProcess;
0052 class G4MaterialCutsCouple;
0053 class G4VSensitiveDetector;
0054
0055 class G4StepPoint
0056 {
0057 public:
0058
0059 G4StepPoint() = default;
0060 ~G4StepPoint()= default;
0061
0062
0063 G4StepPoint(const G4StepPoint&) = default;
0064 G4StepPoint& operator=(const G4StepPoint&);
0065
0066
0067 const G4ThreeVector& GetPosition() const;
0068 void SetPosition(const G4ThreeVector& aValue);
0069 void AddPosition(const G4ThreeVector& aValue);
0070
0071
0072 G4double GetLocalTime() const;
0073 void SetLocalTime(const G4double aValue);
0074 void AddLocalTime(const G4double aValue);
0075
0076
0077 G4double GetGlobalTime() const;
0078 void SetGlobalTime(const G4double aValue);
0079 void AddGlobalTime(const G4double aValue);
0080
0081
0082 G4double GetProperTime() const;
0083 void SetProperTime(const G4double aValue);
0084 void AddProperTime(const G4double aValue);
0085
0086
0087 const G4ThreeVector& GetMomentumDirection() const;
0088 void SetMomentumDirection(const G4ThreeVector& aValue);
0089 void AddMomentumDirection(const G4ThreeVector& aValue);
0090
0091
0092 G4ThreeVector GetMomentum() const;
0093
0094
0095 G4double GetTotalEnergy() const;
0096
0097
0098 G4double GetKineticEnergy() const;
0099 void SetKineticEnergy(const G4double aValue);
0100 void AddKineticEnergy(const G4double aValue);
0101
0102
0103 G4double GetVelocity() const;
0104 void SetVelocity(G4double v);
0105
0106
0107 G4double GetBeta() const;
0108
0109
0110 G4double GetGamma() const;
0111
0112
0113 G4VPhysicalVolume* GetPhysicalVolume() const;
0114
0115 const G4VTouchable* GetTouchable() const;
0116 const G4TouchableHandle& GetTouchableHandle() const;
0117 void SetTouchableHandle(const G4TouchableHandle& apValue);
0118
0119 G4Material* GetMaterial() const;
0120 void SetMaterial(G4Material*);
0121
0122 const G4MaterialCutsCouple* GetMaterialCutsCouple() const;
0123 void SetMaterialCutsCouple(const G4MaterialCutsCouple*);
0124
0125 G4VSensitiveDetector* GetSensitiveDetector() const;
0126 void SetSensitiveDetector(G4VSensitiveDetector*);
0127
0128 G4double GetSafety() const;
0129 void SetSafety(const G4double aValue);
0130
0131 const G4ThreeVector& GetPolarization() const;
0132 void SetPolarization(const G4ThreeVector& aValue);
0133 void AddPolarization(const G4ThreeVector& aValue);
0134
0135 G4StepStatus GetStepStatus() const;
0136 void SetStepStatus(const G4StepStatus aValue);
0137
0138 const G4VProcess* GetProcessDefinedStep() const;
0139
0140
0141 void SetProcessDefinedStep(const G4VProcess* aValue);
0142
0143 G4double GetMass() const;
0144 void SetMass(G4double value);
0145
0146 G4double GetCharge() const;
0147 void SetCharge(G4double value);
0148
0149 G4double GetMagneticMoment() const;
0150 void SetMagneticMoment(G4double value);
0151
0152 void SetWeight(G4double aValue);
0153 G4double GetWeight() const;
0154
0155 private:
0156
0157 G4ThreeVector fPosition;
0158 G4double fGlobalTime = 0.0;
0159
0160 G4double fLocalTime = 0.0;
0161
0162 G4double fProperTime = 0.0;
0163
0164 G4ThreeVector fMomentumDirection;
0165 G4double fKineticEnergy = 0.0;
0166 G4double fVelocity = 0.0;
0167
0168 G4TouchableHandle fpTouchable;
0169
0170 G4Material* fpMaterial = nullptr;
0171
0172 const G4MaterialCutsCouple* fpMaterialCutsCouple = nullptr;
0173
0174 G4VSensitiveDetector* fpSensitiveDetector = nullptr;
0175 G4double fSafety = 0.0;
0176 G4ThreeVector fPolarization;
0177 G4StepStatus fStepStatus = fUndefined;
0178
0179 const G4VProcess* fpProcessDefinedStep = nullptr;
0180
0181 G4double fMass = 0.0;
0182
0183 G4double fCharge = 0.0;
0184
0185 G4double fMagneticMoment = 0.0;
0186
0187 G4double fWeight = 0.0;
0188
0189 };
0190
0191 #include "G4StepPoint.icc"
0192
0193 #endif