Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ParticleChangeForGamma Class Reference

#include <G4ParticleChangeForGamma.hh>

+ Inheritance diagram for G4ParticleChangeForGamma:

Public Member Functions

 G4ParticleChangeForGamma ()
 
virtual ~G4ParticleChangeForGamma ()
 
G4StepUpdateStepForAtRest (G4Step *pStep)
 
G4StepUpdateStepForPostStep (G4Step *Step)
 
void InitializeForPostStep (const G4Track &)
 
void AddSecondary (G4DynamicParticle *aParticle)
 
G4double GetProposedKineticEnergy () const
 
void SetProposedKineticEnergy (G4double proposedKinEnergy)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void ProposeMomentumDirection (G4double Px, G4double Py, G4double Pz)
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedPolarization () const
 
void ProposePolarization (const G4ThreeVector &dir)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
const G4TrackGetCurrentTrack () const
 
virtual void DumpInfo () const
 
virtual G4bool CheckIt (const G4Track &)
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()
 
G4bool operator== (const G4VParticleChange &right) const
 
G4bool operator!= (const G4VParticleChange &right) const
 
virtual G4StepUpdateStepForAtRest (G4Step *Step)
 
virtual G4StepUpdateStepForAlongStep (G4Step *Step)
 
virtual G4StepUpdateStepForPostStep (G4Step *Step)
 
virtual void Initialize (const G4Track &)
 
G4double GetTrueStepLength () const
 
void ProposeTrueStepLength (G4double truePathLength)
 
G4double GetLocalEnergyDeposit () const
 
void ProposeLocalEnergyDeposit (G4double anEnergyPart)
 
G4double GetNonIonizingEnergyDeposit () const
 
void ProposeNonIonizingEnergyDeposit (G4double anEnergyPart)
 
G4TrackStatus GetTrackStatus () const
 
void ProposeTrackStatus (G4TrackStatus status)
 
G4SteppingControl GetSteppingControl () const
 
void ProposeSteppingControl (G4SteppingControl StepControlFlag)
 
G4bool GetFirstStepInVolume () const
 
G4bool GetLastStepInVolume () const
 
void ProposeFirstStepInVolume (G4bool flag)
 
void ProposeLastStepInVolume (G4bool flag)
 
void Clear ()
 
void SetNumberOfSecondaries (G4int totSecondaries)
 
G4int GetNumberOfSecondaries () const
 
G4TrackGetSecondary (G4int anIndex) const
 
void AddSecondary (G4Track *aSecondary)
 
G4double GetWeight () const
 
G4double GetParentWeight () const
 
void ProposeWeight (G4double finalWeight)
 
void ProposeParentWeight (G4double finalWeight)
 
void SetSecondaryWeightByProcess (G4bool)
 
G4bool IsSecondaryWeightSetByProcess () const
 
void SetParentWeightByProcess (G4bool)
 
G4bool IsParentWeightSetByProcess () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int vLevel)
 
G4int GetVerboseLevel () const
 
virtual G4bool CheckIt (const G4Track &)
 
void ClearDebugFlag ()
 
void SetDebugFlag ()
 
G4bool GetDebugFlag () const
 

Protected Member Functions

 G4ParticleChangeForGamma (const G4ParticleChangeForGamma &right)
 
G4ParticleChangeForGammaoperator= (const G4ParticleChangeForGamma &right)
 
- Protected Member Functions inherited from G4VParticleChange
 G4VParticleChange (const G4VParticleChange &right)
 
G4VParticleChangeoperator= (const G4VParticleChange &right)
 
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeTrueStepLength (const G4Track &)
 
void InitializeLocalEnergyDeposit (const G4Track &)
 
void InitializeSteppingControl (const G4Track &)
 
void InitializeParentWeight (const G4Track &)
 
void InitializeParentGlobalTime (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries (const G4Track &)
 
void InitializeStepInVolumeFlags (const G4Track &)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VParticleChange
G4TrackFastVectortheListOfSecondaries = nullptr
 
G4TrackStatus theStatusChange = fAlive
 
G4SteppingControl theSteppingControlFlag = NormalCondition
 
G4double theLocalEnergyDeposit = 0.0
 
G4double theNonIonizingEnergyDeposit = 0.0
 
G4double theTrueStepLength = 0.0
 
G4double theParentWeight = 1.0
 
G4double theParentGlobalTime = 0.0
 
G4int theNumberOfSecondaries = 0
 
G4int theSizeOftheListOfSecondaries = 0
 
G4int verboseLevel = 1
 
G4bool theFirstStepInVolume = false
 
G4bool theLastStepInVolume = false
 
G4bool isParentWeightProposed = false
 
G4bool fSetSecondaryWeightByProcess = false
 
G4bool debugFlag = false
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 

Detailed Description

Definition at line 44 of file G4ParticleChangeForGamma.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForGamma() [1/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( )

Definition at line 40 of file G4ParticleChangeForGamma.cc.

◆ ~G4ParticleChangeForGamma()

G4ParticleChangeForGamma::~G4ParticleChangeForGamma ( )
virtual

Definition at line 48 of file G4ParticleChangeForGamma.cc.

49{
50}

◆ G4ParticleChangeForGamma() [2/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( const G4ParticleChangeForGamma right)
protected

Definition at line 53 of file G4ParticleChangeForGamma.cc.

55 : G4VParticleChange(right)
56{
57 currentTrack = right.currentTrack;
58 proposedKinEnergy = right.proposedKinEnergy;
59 proposedMomentumDirection = right.proposedMomentumDirection;
60 proposedPolarization = right.proposedPolarization;
61}

Member Function Documentation

◆ AddSecondary()

void G4ParticleChangeForGamma::AddSecondary ( G4DynamicParticle aParticle)

Definition at line 100 of file G4ParticleChangeForGamma.cc.

101{
102 // create track
103 G4Track* aTrack = new G4Track(aParticle, currentTrack->GetGlobalTime(),
104 currentTrack->GetPosition());
105
106 // touchable handle is copied to keep the pointer
107 aTrack->SetTouchableHandle(currentTrack->GetTouchableHandle());
108
109 // add a secondary
111}
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void AddSecondary(G4Track *aSecondary)

◆ CheckIt()

G4bool G4ParticleChangeForGamma::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 173 of file G4ParticleChangeForGamma.cc.

174{
175 G4bool itsOK = true;
176 G4bool exitWithError = false;
177
178 G4double accuracy;
179
180 // Energy should not be lager than initial value
181 accuracy = (proposedKinEnergy - aTrack.GetKineticEnergy()) / MeV;
182 if(accuracy > accuracyForWarning)
183 {
184 itsOK = false;
185 exitWithError = (accuracy > accuracyForException);
186#ifdef G4VERBOSE
187 G4cout << "G4ParticleChangeForGamma::CheckIt: ";
188 G4cout << "KinEnergy become larger than the initial value!"
189 << " Difference: " << accuracy << "[MeV] " << G4endl;
191 << " E=" << aTrack.GetKineticEnergy() / MeV
192 << " pos=" << aTrack.GetPosition().x() / m << ", "
193 << aTrack.GetPosition().y() / m << ", "
194 << aTrack.GetPosition().z() / m << G4endl;
195#endif
196 }
197
198 // dump out information of this particle change
199#ifdef G4VERBOSE
200 if(!itsOK) { DumpInfo(); }
201#endif
202
203 // Exit with error
204 if(exitWithError)
205 {
206 G4Exception("G4ParticleChangeForGamma::CheckIt()", "TRACK004",
207 EventMustBeAborted, "energy was illegal");
208 }
209
210 // correction
211 if(!itsOK)
212 {
213 proposedKinEnergy = aTrack.GetKineticEnergy();
214 }
215
216 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
217 return itsOK;
218}
@ EventMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
double x() const
double y() const
const G4String & GetParticleName() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
virtual G4bool CheckIt(const G4Track &)
static const G4double accuracyForException
static const G4double accuracyForWarning

◆ DumpInfo()

void G4ParticleChangeForGamma::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 157 of file G4ParticleChangeForGamma.cc.

158{
159 // use base-class DumpInfo
161
162 G4int oldprc = G4cout.precision(3);
163 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
164 << proposedKinEnergy / MeV << G4endl;
165 G4cout << " Momentum Direction: " << std::setw(20)
166 << proposedMomentumDirection << G4endl;
167 G4cout << " Polarization: " << std::setw(20) << proposedPolarization
168 << G4endl;
169 G4cout.precision(oldprc);
170}
int G4int
Definition: G4Types.hh:85
virtual void DumpInfo() const

Referenced by CheckIt().

◆ GetCurrentTrack()

◆ GetProposedKineticEnergy()

G4double G4ParticleChangeForGamma::GetProposedKineticEnergy ( ) const
inline

Definition at line 113 of file G4ParticleChangeForGamma.hh.

114{
115 return proposedKinEnergy;
116}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), and G4VEmProcess::PostStepDoIt().

◆ GetProposedMomentumDirection()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedMomentumDirection ( ) const
inline

Definition at line 125 of file G4ParticleChangeForGamma.hh.

127{
128 return proposedMomentumDirection;
129}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), and G4LEPTSElasticModel::SampleSecondaries().

◆ GetProposedPolarization()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedPolarization ( ) const
inline

Definition at line 155 of file G4ParticleChangeForGamma.hh.

156{
157 return proposedPolarization;
158}

◆ InitializeForPostStep()

void G4ParticleChangeForGamma::InitializeForPostStep ( const G4Track track)
inline

Definition at line 177 of file G4ParticleChangeForGamma.hh.

178{
183 theParentWeight = track.GetWeight();
185 proposedKinEnergy = track.GetKineticEnergy();
186 proposedMomentumDirection = track.GetMomentumDirection();
187 proposedPolarization = track.GetPolarization();
188 currentTrack = &track;
189}
G4TrackStatus GetTrackStatus() const
G4double GetWeight() const
const G4ThreeVector & GetMomentumDirection() const
const G4ThreeVector & GetPolarization() const
G4TrackStatus theStatusChange
G4double theNonIonizingEnergyDeposit
void InitializeSecondaries(const G4Track &)

Referenced by G4eplusAnnihilation::AtRestDoIt(), G4GammaGeneralProcess::PostStepDoIt(), and G4VEmProcess::PostStepDoIt().

◆ operator=()

G4ParticleChangeForGamma & G4ParticleChangeForGamma::operator= ( const G4ParticleChangeForGamma right)
protected

Definition at line 65 of file G4ParticleChangeForGamma.cc.

66{
67 if(this != &right)
68 {
70 {
71 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
72 {
73 if((*theListOfSecondaries)[index])
74 delete(*theListOfSecondaries)[index];
75 }
76 }
80 for(G4int index = 0; index < theNumberOfSecondaries; ++index)
81 {
82 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index]));
83 theListOfSecondaries->SetElement(index, newTrack);
84 }
85
90
91 currentTrack = right.currentTrack;
92 proposedKinEnergy = right.proposedKinEnergy;
93 proposedMomentumDirection = right.proposedMomentumDirection;
94 proposedPolarization = right.proposedPolarization;
95 }
96 return *this;
97}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:72
G4TrackFastVector * theListOfSecondaries

◆ ProposeMomentumDirection() [1/2]

void G4ParticleChangeForGamma::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

Definition at line 132 of file G4ParticleChangeForGamma.hh.

134{
135 proposedMomentumDirection = dir;
136}

◆ ProposeMomentumDirection() [2/2]

void G4ParticleChangeForGamma::ProposeMomentumDirection ( G4double  Px,
G4double  Py,
G4double  Pz 
)
inline

Definition at line 139 of file G4ParticleChangeForGamma.hh.

142{
143 proposedMomentumDirection.setX(Px);
144 proposedMomentumDirection.setY(Py);
145 proposedMomentumDirection.setZ(Pz);
146}
void setY(double)
void setZ(double)
void setX(double)

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4XrayRayleighModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MicroElecLOPhononModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eDPWACoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), G4LEPTSVibExcitationModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ProposePolarization() [1/2]

◆ ProposePolarization() [2/2]

void G4ParticleChangeForGamma::ProposePolarization ( G4double  Px,
G4double  Py,
G4double  Pz 
)
inline

Definition at line 167 of file G4ParticleChangeForGamma.hh.

170{
171 proposedPolarization.setX(Px);
172 proposedPolarization.setY(Py);
173 proposedPolarization.setZ(Pz);
174}

◆ SetProposedKineticEnergy()

void G4ParticleChangeForGamma::SetProposedKineticEnergy ( G4double  proposedKinEnergy)
inline

Definition at line 119 of file G4ParticleChangeForGamma.hh.

120{
121 proposedKinEnergy = energy;
122}
G4double energy(const ThreeVector &p, const G4double m)

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricGDModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4XrayRayleighModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MicroElecLOPhononModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4eplusTo2GammaOKVIModel::SampleSecondaries(), G4eplusTo3GammaOKVIModel::SampleSecondaries(), G4LEPTSAttachmentModel::SampleSecondaries(), G4LEPTSDissociationModel::SampleSecondaries(), G4LEPTSElasticModel::SampleSecondaries(), G4LEPTSExcitationModel::SampleSecondaries(), G4LEPTSIonisationModel::SampleSecondaries(), G4LEPTSPositroniumModel::SampleSecondaries(), G4LEPTSRotExcitationModel::SampleSecondaries(), G4LEPTSVibExcitationModel::SampleSecondaries(), G4DNAPTBElasticModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), G4DNAPTBIonisationModel::SampleSecondaries(), and G4BetheHeitler5DModel::SampleSecondaries().

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForGamma::UpdateStepForAtRest ( G4Step pStep)
virtual

Reimplemented from G4VParticleChange.

Definition at line 114 of file G4ParticleChangeForGamma.cc.

115{
117 pStep->SetStepLength(0.0);
118
120 {
122 }
123
124 return pStep;
125}
void SetWeight(G4double aValue)
void SetStepLength(G4double value)
void AddTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForGamma::UpdateStepForPostStep ( G4Step Step)
virtual

Reimplemented from G4VParticleChange.

Definition at line 128 of file G4ParticleChangeForGamma.cc.

129{
130 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
131 G4Track* pTrack = pStep->GetTrack();
132
133 pPostStepPoint->SetKineticEnergy(proposedKinEnergy);
134 pPostStepPoint->SetMomentumDirection(proposedMomentumDirection);
135 pPostStepPoint->SetPolarization(proposedPolarization);
136
137 // update velocity for scattering process and particles with mass
138 if(proposedKinEnergy > 0.0)
139 {
140 if(pTrack->GetParticleDefinition()->GetPDGMass() > 0.0)
141 {
142 pPostStepPoint->SetVelocity(pTrack->CalculateVelocity());
143 }
144 }
145
147 {
148 pPostStepPoint->SetWeight(theParentWeight);
149 }
150
151 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
152 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
153 return pStep;
154}
void SetKineticEnergy(const G4double aValue)
void SetVelocity(G4double v)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
const G4ParticleDefinition * GetParticleDefinition() const
G4double CalculateVelocity() const

The documentation for this class was generated from the following files: