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

#include <G4ParticleChangeForGamma.hh>

+ Inheritance diagram for G4ParticleChangeForGamma:

Public Member Functions

 G4ParticleChangeForGamma ()
 
 ~G4ParticleChangeForGamma () override=default
 
 G4ParticleChangeForGamma (const G4ParticleChangeForGamma &right)=delete
 
G4ParticleChangeForGammaoperator= (const G4ParticleChangeForGamma &right)=delete
 
G4StepUpdateStepForAtRest (G4Step *pStep) final
 
G4StepUpdateStepForPostStep (G4Step *Step) final
 
void InitializeForPostStep (const G4Track &)
 
void AddSecondary (G4DynamicParticle *aParticle)
 
G4double GetProposedKineticEnergy () const
 
void SetProposedKineticEnergy (G4double proposedKinEnergy)
 
const G4ThreeVectorGetProposedMomentumDirection () const
 
void ProposeMomentumDirection (const G4ThreeVector &Pfinal)
 
const G4ThreeVectorGetProposedPolarization () const
 
void ProposePolarization (const G4ThreeVector &dir)
 
void ProposePolarization (G4double Px, G4double Py, G4double Pz)
 
void DumpInfo () const override
 
- Public Member Functions inherited from G4VParticleChange
 G4VParticleChange ()
 
virtual ~G4VParticleChange ()=default
 
 G4VParticleChange (const G4VParticleChange &right)=delete
 
G4VParticleChangeoperator= (const G4VParticleChange &right)=delete
 
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)
 
const G4TrackGetCurrentTrack () const
 
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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VParticleChange
G4StepUpdateStepInfo (G4Step *Step)
 
void InitializeLocalEnergyDeposit ()
 
void InitializeSteppingControl ()
 
void InitializeParentWeight (const G4Track &)
 
void InitializeStatusChange (const G4Track &)
 
void InitializeSecondaries ()
 
void InitializeFromStep (const G4Step *)
 
G4double ComputeBeta (G4double kinEnergy)
 
G4bool CheckSecondary (G4Track &)
 
G4double GetAccuracyForWarning () const
 
G4double GetAccuracyForException () const
 
- Protected Attributes inherited from G4VParticleChange
const G4TracktheCurrentTrack = nullptr
 
std::vector< G4Track * > theListOfSecondaries
 
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
 
G4int nError = 0
 
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
 
static const G4int maxError = 10
 

Detailed Description

Definition at line 45 of file G4ParticleChangeForGamma.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForGamma() [1/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( )

Definition at line 40 of file G4ParticleChangeForGamma.cc.

41{
42 // Disable flag to avoid check of each secondary at each step
43 debugFlag = false;
44}

◆ ~G4ParticleChangeForGamma()

G4ParticleChangeForGamma::~G4ParticleChangeForGamma ( )
overridedefault

◆ G4ParticleChangeForGamma() [2/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( const G4ParticleChangeForGamma right)
delete

Member Function Documentation

◆ AddSecondary()

void G4ParticleChangeForGamma::AddSecondary ( G4DynamicParticle aParticle)

Definition at line 47 of file G4ParticleChangeForGamma.cc.

48{
49 // create track
50 G4Track* aTrack = new G4Track(aParticle, theCurrentTrack->GetGlobalTime(),
52
53 // touchable handle is copied to keep the pointer
55
56 // add a secondary
58}
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void AddSecondary(G4Track *aSecondary)
const G4Track * theCurrentTrack

◆ DumpInfo()

void G4ParticleChangeForGamma::DumpInfo ( ) const
overridevirtual

Reimplemented from G4VParticleChange.

Definition at line 111 of file G4ParticleChangeForGamma.cc.

112{
113 // use base-class DumpInfo
115
116 G4long oldprc = G4cout.precision(8);
117 G4cout << " -----------------------------------------------" << G4endl;
118 G4cout << " G4ParticleChangeForGamma proposes: " << G4endl;
119 G4cout << " Kinetic Energy (MeV): " << std::setw(20)
120 << proposedKinEnergy / MeV << G4endl;
121 G4cout << " Momentum Direction: " << std::setw(20)
122 << proposedMomentumDirection << G4endl;
123 G4cout << " Polarization: " << std::setw(20) << proposedPolarization
124 << G4endl;
125 G4cout.precision(oldprc);
126}
long G4long
Definition: G4Types.hh:87
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void DumpInfo() const

◆ GetProposedKineticEnergy()

G4double G4ParticleChangeForGamma::GetProposedKineticEnergy ( ) const
inline

Definition at line 100 of file G4ParticleChangeForGamma.hh.

101{
102 return proposedKinEnergy;
103}

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

◆ GetProposedMomentumDirection()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedMomentumDirection ( ) const
inline

Definition at line 112 of file G4ParticleChangeForGamma.hh.

114{
115 return proposedMomentumDirection;
116}

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

◆ GetProposedPolarization()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedPolarization ( ) const
inline

Definition at line 126 of file G4ParticleChangeForGamma.hh.

127{
128 return proposedPolarization;
129}

◆ InitializeForPostStep()

void G4ParticleChangeForGamma::InitializeForPostStep ( const G4Track track)
inline

Definition at line 148 of file G4ParticleChangeForGamma.hh.

149{
154 proposedKinEnergy = track.GetKineticEnergy();
155 proposedMomentumDirection = track.GetMomentumDirection();
156 proposedPolarization = track.GetPolarization();
157}
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
void InitializeLocalEnergyDeposit()
void InitializeStatusChange(const G4Track &)
void InitializeSecondaries()
void InitializeParentWeight(const G4Track &)

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

◆ operator=()

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

◆ ProposeMomentumDirection()

void G4ParticleChangeForGamma::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

Definition at line 119 of file G4ParticleChangeForGamma.hh.

121{
122 proposedMomentumDirection = dir;
123}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4DNABornExcitationModel1::SampleSecondaries(), G4DNABornExcitationModel2::SampleSecondaries(), G4DNABornIonisationModel1::SampleSecondaries(), G4DNABornIonisationModel2::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNACPA100ElasticModel::SampleSecondaries(), G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNADiracRMatrixExcitationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4DNARPWBAExcitationModel::SampleSecondaries(), G4DNARPWBAIonisationModel::SampleSecondaries(), G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MicroElecLOPhononModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eDPWACoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4XrayRayleighModel::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 138 of file G4ParticleChangeForGamma.hh.

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

◆ SetProposedKineticEnergy()

void G4ParticleChangeForGamma::SetProposedKineticEnergy ( G4double  proposedKinEnergy)
inline

Definition at line 106 of file G4ParticleChangeForGamma.hh.

107{
108 proposedKinEnergy = energy;
109}
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(), G4DNADiracRMatrixExcitationModel::SampleSecondaries(), G4DNAELSEPAElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAEmfietzoglouIonisationModel::SampleSecondaries(), G4DNAIonElasticModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNAQuinnPlasmonExcitationModel::SampleSecondaries(), G4DNARelativisticIonisationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4DNARPWBAExcitationModel::SampleSecondaries(), G4DNARPWBAIonisationModel::SampleSecondaries(), G4DNAUeharaScreenedRutherfordElasticModel::SampleSecondaries(), G4JAEAElasticScatteringModel::SampleSecondaries(), G4JAEAPolarizedElasticScatteringModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4LowEPPolarizedComptonModel::SampleSecondaries(), G4MicroElecElasticModel::SampleSecondaries(), G4MicroElecElasticModel_new::SampleSecondaries(), G4MicroElecInelasticModel::SampleSecondaries(), G4MicroElecInelasticModel_new::SampleSecondaries(), G4MicroElecLOPhononModel::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PenelopeRayleighModelMI::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4XrayRayleighModel::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)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 61 of file G4ParticleChangeForGamma.cc.

62{
64 pStep->SetStepLength(0.0);
65
67 {
69 }
70
71 return pStep;
72}
void SetWeight(G4double aValue)
void SetStepLength(G4double value)
void AddTotalEnergyDeposit(G4double value)
G4StepPoint * GetPostStepPoint() const

◆ UpdateStepForPostStep()

G4Step * G4ParticleChangeForGamma::UpdateStepForPostStep ( G4Step Step)
finalvirtual

Reimplemented from G4VParticleChange.

Definition at line 75 of file G4ParticleChangeForGamma.cc.

76{
77 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
78
79 pPostStepPoint->SetMomentumDirection(proposedMomentumDirection);
80 pPostStepPoint->SetPolarization(proposedPolarization);
81
82 // update velocity for scattering process and particles with mass
83 if(proposedKinEnergy > 0.0)
84 {
85 pPostStepPoint->SetKineticEnergy(proposedKinEnergy);
87 G4double v = CLHEP::c_light;
88 if(mass > 0.0) {
89 v *= std::sqrt(proposedKinEnergy*(proposedKinEnergy + 2*mass))/
90 (proposedKinEnergy + mass);
91 }
92 pPostStepPoint->SetVelocity(v);
93 }
94 else
95 {
96 pPostStepPoint->SetKineticEnergy(0.0);
97 pPostStepPoint->SetVelocity(0.0);
98 }
99
101 {
102 pPostStepPoint->SetWeight(theParentWeight);
103 }
104
105 pStep->AddTotalEnergyDeposit(theLocalEnergyDeposit);
106 pStep->AddNonIonizingEnergyDeposit(theNonIonizingEnergyDeposit);
107 return pStep;
108}
double G4double
Definition: G4Types.hh:83
void SetKineticEnergy(const G4double aValue)
void SetVelocity(G4double v)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)
G4ParticleDefinition * GetDefinition() const
G4double theNonIonizingEnergyDeposit

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