File indexing completed on 2025-01-18 09:58:51
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
0036 #ifndef G4ParticleChangeForLoss_hh
0037 #define G4ParticleChangeForLoss_hh 1
0038
0039 #include "G4VParticleChange.hh"
0040 #include "G4DynamicParticle.hh"
0041
0042 class G4ParticleChangeForLoss final : public G4VParticleChange
0043 {
0044 public:
0045
0046 G4ParticleChangeForLoss();
0047
0048 ~G4ParticleChangeForLoss() override = default;
0049
0050 G4ParticleChangeForLoss(const G4ParticleChangeForLoss& right) = delete;
0051 G4ParticleChangeForLoss& operator=(const G4ParticleChangeForLoss& right) = delete;
0052
0053
0054
0055 G4Step* UpdateStepForAlongStep(G4Step* step) final;
0056 G4Step* UpdateStepForPostStep(G4Step* step) final;
0057
0058
0059 inline void InitializeForAlongStep(const G4Track&);
0060 inline void InitializeForPostStep(const G4Track&);
0061
0062
0063 inline G4double GetProposedCharge() const;
0064 inline void SetProposedCharge(G4double theCharge);
0065
0066
0067 inline G4double GetProposedKineticEnergy() const;
0068 inline void SetProposedKineticEnergy(G4double proposedKinEnergy);
0069
0070
0071
0072 inline const G4ThreeVector& GetProposedMomentumDirection() const;
0073 inline void SetProposedMomentumDirection(const G4ThreeVector& dir);
0074 inline void ProposeMomentumDirection(const G4ThreeVector& Pfinal);
0075
0076 inline const G4ThreeVector& GetProposedPolarization() const;
0077 inline void ProposePolarization(const G4ThreeVector& dir);
0078 inline void ProposePolarization(G4double Px, G4double Py, G4double Pz);
0079
0080 void DumpInfo() const final;
0081
0082 private:
0083
0084 G4double proposedKinEnergy = 0.0;
0085
0086
0087 G4double currentCharge = 0.0;
0088
0089
0090 G4ThreeVector proposedMomentumDirection;
0091
0092
0093 G4ThreeVector proposedPolarization;
0094
0095 };
0096
0097
0098
0099
0100
0101 inline
0102 G4double G4ParticleChangeForLoss::GetProposedKineticEnergy() const
0103 {
0104 return proposedKinEnergy;
0105 }
0106
0107 inline
0108 void G4ParticleChangeForLoss::SetProposedKineticEnergy(G4double energy)
0109 {
0110 proposedKinEnergy = energy;
0111 }
0112
0113 inline G4double G4ParticleChangeForLoss::GetProposedCharge() const
0114 {
0115 return currentCharge;
0116 }
0117
0118 inline
0119 void G4ParticleChangeForLoss::SetProposedCharge(G4double theCharge)
0120 {
0121 currentCharge = theCharge;
0122 }
0123
0124 inline
0125 const G4ThreeVector&
0126 G4ParticleChangeForLoss::GetProposedMomentumDirection() const
0127 {
0128 return proposedMomentumDirection;
0129 }
0130
0131 inline
0132 void G4ParticleChangeForLoss::ProposeMomentumDirection(const G4ThreeVector& dir)
0133 {
0134 proposedMomentumDirection = dir;
0135 }
0136
0137 inline
0138 void
0139 G4ParticleChangeForLoss::SetProposedMomentumDirection(const G4ThreeVector& dir)
0140 {
0141 proposedMomentumDirection = dir;
0142 }
0143
0144 inline
0145 const G4ThreeVector& G4ParticleChangeForLoss::GetProposedPolarization() const
0146 {
0147 return proposedPolarization;
0148 }
0149
0150 inline
0151 void G4ParticleChangeForLoss::ProposePolarization(const G4ThreeVector& dir)
0152 {
0153 proposedPolarization = dir;
0154 }
0155
0156 inline void G4ParticleChangeForLoss::ProposePolarization(G4double Px,
0157 G4double Py,
0158 G4double Pz)
0159 {
0160 proposedPolarization.set(Px, Py, Pz);
0161 }
0162
0163 inline
0164 void G4ParticleChangeForLoss::InitializeForAlongStep(const G4Track& track)
0165 {
0166 InitializeSecondaries();
0167 InitializeLocalEnergyDeposit();
0168 InitializeParentWeight(track);
0169 InitializeStatusChange(track);
0170 proposedKinEnergy = track.GetKineticEnergy();
0171 currentCharge = track.GetDynamicParticle()->GetCharge();
0172 }
0173
0174 inline
0175 void G4ParticleChangeForLoss::InitializeForPostStep(const G4Track& track)
0176 {
0177 InitializeForAlongStep(track);
0178 proposedMomentumDirection = track.GetMomentumDirection();
0179 proposedPolarization = track.GetPolarization();
0180 }
0181
0182 #endif