Geant4 9.6.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
 
G4int theNumberOfSecondaries
 
G4int theSizeOftheListOfSecondaries
 
G4TrackStatus theStatusChange
 
G4SteppingControl theSteppingControlFlag
 
G4double theLocalEnergyDeposit
 
G4double theNonIonizingEnergyDeposit
 
G4double theTrueStepLength
 
G4bool theFirstStepInVolume
 
G4bool theLastStepInVolume
 
G4double theParentWeight
 
G4bool isParentWeightProposed
 
G4bool fSetSecondaryWeightByProcess
 
G4double theParentGlobalTime
 
G4int verboseLevel
 
G4bool debugFlag
 
- Static Protected Attributes inherited from G4VParticleChange
static const G4double accuracyForWarning = 1.0e-9
 
static const G4double accuracyForException = 0.001
 

Detailed Description

Definition at line 58 of file G4ParticleChangeForGamma.hh.

Constructor & Destructor Documentation

◆ G4ParticleChangeForGamma() [1/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( )

Definition at line 49 of file G4ParticleChangeForGamma.cc.

50 : G4VParticleChange(), currentTrack(0), proposedKinEnergy(0.)
51{
53 debugFlag = false;
54#ifdef G4VERBOSE
55 if (verboseLevel>2) {
56 G4cout << "G4ParticleChangeForGamma::G4ParticleChangeForGamma() " << G4endl;
57 }
58#endif
59}
@ NormalCondition
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4SteppingControl theSteppingControlFlag

◆ ~G4ParticleChangeForGamma()

G4ParticleChangeForGamma::~G4ParticleChangeForGamma ( )
virtual

Definition at line 61 of file G4ParticleChangeForGamma.cc.

62{
63#ifdef G4VERBOSE
64 if (verboseLevel>2) {
65 G4cout << "G4ParticleChangeForGamma::~G4ParticleChangeForGamma() " << G4endl;
66 }
67#endif
68}

◆ G4ParticleChangeForGamma() [2/2]

G4ParticleChangeForGamma::G4ParticleChangeForGamma ( const G4ParticleChangeForGamma right)
protected

Definition at line 70 of file G4ParticleChangeForGamma.cc.

71 : G4VParticleChange(right)
72{
73 if (verboseLevel>1)
74 G4cout << "G4ParticleChangeForGamma:: copy constructor is called " << G4endl;
75
76 currentTrack = right.currentTrack;
77 proposedKinEnergy = right.proposedKinEnergy;
78 proposedMomentumDirection = right.proposedMomentumDirection;
79 proposedPolarization = right.proposedPolarization;
80}

Member Function Documentation

◆ AddSecondary()

void G4ParticleChangeForGamma::AddSecondary ( G4DynamicParticle aParticle)

Definition at line 155 of file G4ParticleChangeForGamma.cc.

156{
157 // create track
158 G4Track* aTrack = new G4Track(aParticle, currentTrack->GetGlobalTime(),
159 currentTrack->GetPosition());
160
161 // Touchable handle is copied to keep the pointer
162 aTrack->SetTouchableHandle(currentTrack->GetTouchableHandle());
163
164 // add a secondary
166}
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4TouchableHandle & GetTouchableHandle() const
void AddSecondary(G4Track *aSecondary)

Referenced by G4eplusPolarizedAnnihilation::AtRestDoIt(), and G4eplusAnnihilation::AtRestDoIt().

◆ CheckIt()

G4bool G4ParticleChangeForGamma::CheckIt ( const G4Track aTrack)
virtual

Reimplemented from G4VParticleChange.

Definition at line 190 of file G4ParticleChangeForGamma.cc.

191{
192 G4bool itsOK = true;
193 G4bool exitWithError = false;
194
195 G4double accuracy;
196
197 // Energy should not be lager than initial value
198 accuracy = ( proposedKinEnergy - aTrack.GetKineticEnergy())/MeV;
199 if (accuracy > accuracyForWarning) {
200 itsOK = false;
201 exitWithError = (accuracy > accuracyForException);
202#ifdef G4VERBOSE
203 G4cout << "G4ParticleChangeForGamma::CheckIt: ";
204 G4cout << "KinEnergy become larger than the initial value!"
205 << " Difference: " << accuracy << "[MeV] " <<G4endl;
207 << " E=" << aTrack.GetKineticEnergy()/MeV
208 << " pos=" << aTrack.GetPosition().x()/m
209 << ", " << aTrack.GetPosition().y()/m
210 << ", " << aTrack.GetPosition().z()/m
211 << G4endl;
212#endif
213 }
214
215 // dump out information of this particle change
216#ifdef G4VERBOSE
217 if (!itsOK) DumpInfo();
218#endif
219
220 // Exit with error
221 if (exitWithError) {
222 G4Exception("G4ParticleChangeForGamma::CheckIt",
223 "TRACK004",EventMustBeAborted,
224 "energy was illegal");
225 }
226
227 //correction
228 if (!itsOK) {
229 proposedKinEnergy = aTrack.GetKineticEnergy();
230 }
231
232 itsOK = (itsOK) && G4VParticleChange::CheckIt(aTrack);
233 return itsOK;
234}
@ EventMustBeAborted
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
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
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ DumpInfo()

void G4ParticleChangeForGamma::DumpInfo ( ) const
virtual

Reimplemented from G4VParticleChange.

Definition at line 172 of file G4ParticleChangeForGamma.cc.

173{
174// use base-class DumpInfo
176
177 G4int oldprc = G4cout.precision(3);
178 G4cout << " Kinetic Energy (MeV): "
179 << std::setw(20) << proposedKinEnergy/MeV
180 << G4endl;
181 G4cout << " Momentum Direction: "
182 << std::setw(20) << proposedMomentumDirection
183 << G4endl;
184 G4cout << " Polarization: "
185 << std::setw(20) << proposedPolarization
186 << G4endl;
187 G4cout.precision(oldprc);
188}
int G4int
Definition: G4Types.hh:66
virtual void DumpInfo() const

Referenced by CheckIt().

◆ GetCurrentTrack()

◆ GetProposedKineticEnergy()

G4double G4ParticleChangeForGamma::GetProposedKineticEnergy ( ) const
inline

Definition at line 124 of file G4ParticleChangeForGamma.hh.

125{
126 return proposedKinEnergy;
127}

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

◆ GetProposedMomentumDirection()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedMomentumDirection ( ) const
inline

Definition at line 135 of file G4ParticleChangeForGamma.hh.

136{
137 return proposedMomentumDirection;
138}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing().

◆ GetProposedPolarization()

const G4ThreeVector & G4ParticleChangeForGamma::GetProposedPolarization ( ) const
inline

Definition at line 160 of file G4ParticleChangeForGamma.hh.

161{
162 return proposedPolarization;
163}

◆ InitializeForPostStep()

void G4ParticleChangeForGamma::InitializeForPostStep ( const G4Track track)
inline

Definition at line 179 of file G4ParticleChangeForGamma.hh.

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

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

◆ operator=()

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

Definition at line 83 of file G4ParticleChangeForGamma.cc.

85{
86#ifdef G4VERBOSE
87 if (verboseLevel>1)
88 G4cout << "G4ParticleChangeForGamma:: assignment operator is called " << G4endl;
89#endif
90
91 if (this != &right) {
93#ifdef G4VERBOSE
94 if (verboseLevel>0) {
95 G4cout << "G4ParticleChangeForGamma: assignment operator Warning ";
96 G4cout << "theListOfSecondaries is not empty ";
97 }
98#endif
99 for (G4int index= 0; index<theNumberOfSecondaries; index++){
100 if ( (*theListOfSecondaries)[index] ) delete (*theListOfSecondaries)[index] ;
101 }
102 }
103 delete theListOfSecondaries;
106 for (G4int index = 0; index<theNumberOfSecondaries; index++){
107 G4Track* newTrack = new G4Track(*((*right.theListOfSecondaries)[index] ));
108 theListOfSecondaries->SetElement(index, newTrack); }
109
114
115 currentTrack = right.currentTrack;
116 proposedKinEnergy = right.proposedKinEnergy;
117 proposedMomentumDirection = right.proposedMomentumDirection;
118 proposedPolarization = right.proposedPolarization;
119 }
120 return *this;
121}
G4FastVector< G4Track, G4TrackFastVectorSize > G4TrackFastVector
void SetElement(G4int anIndex, Type *anElement)
Definition: G4FastVector.hh:76
G4TrackFastVector * theListOfSecondaries

◆ ProposeMomentumDirection() [1/2]

void G4ParticleChangeForGamma::ProposeMomentumDirection ( const G4ThreeVector Pfinal)
inline

Definition at line 141 of file G4ParticleChangeForGamma.hh.

142{
143 proposedMomentumDirection = dir;
144}

◆ ProposeMomentumDirection() [2/2]

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

Definition at line 147 of file G4ParticleChangeForGamma.hh.

148{
149 proposedMomentumDirection.setX(Px);
150 proposedMomentumDirection.setY(Py);
151 proposedMomentumDirection.setZ(Pz);
152}
void setY(double)
void setZ(double)
void setX(double)

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4DNABornExcitationModel::SampleSecondaries(), G4DNABornIonisationModel::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4MuElecElasticModel::SampleSecondaries(), G4MuElecInelasticModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), and G4XrayRayleighModel::SampleSecondaries().

◆ ProposePolarization() [1/2]

void G4ParticleChangeForGamma::ProposePolarization ( const G4ThreeVector dir)
inline

◆ ProposePolarization() [2/2]

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

Definition at line 172 of file G4ParticleChangeForGamma.hh.

173{
174 proposedPolarization.setX(Px);
175 proposedPolarization.setY(Py);
176 proposedPolarization.setZ(Pz);
177}

◆ SetProposedKineticEnergy()

void G4ParticleChangeForGamma::SetProposedKineticEnergy ( G4double  proposedKinEnergy)
inline

Definition at line 129 of file G4ParticleChangeForGamma.hh.

130{
131 proposedKinEnergy = energy;
132}

Referenced by G4EmBiasingManager::ApplySecondaryBiasing(), G4DNABornExcitationModel::SampleSecondaries(), G4DNABornIonisationModel::SampleSecondaries(), G4DNAChampionElasticModel::SampleSecondaries(), G4DNAEmfietzoglouExcitationModel::SampleSecondaries(), G4DNAMeltonAttachmentModel::SampleSecondaries(), G4DNAMillerGreenExcitationModel::SampleSecondaries(), G4DNARuddIonisationExtendedModel::SampleSecondaries(), G4DNARuddIonisationModel::SampleSecondaries(), G4DNASancheExcitationModel::SampleSecondaries(), G4DNASancheSolvatationModel::SampleSecondaries(), G4DNAScreenedRutherfordElasticModel::SampleSecondaries(), G4BoldyshevTripletModel::SampleSecondaries(), G4LivermoreComptonModel::SampleSecondaries(), G4LivermoreComptonModifiedModel::SampleSecondaries(), G4LivermoreGammaConversionModel::SampleSecondaries(), G4LivermoreGammaConversionModelRC::SampleSecondaries(), G4LivermoreNuclearGammaConversionModel::SampleSecondaries(), G4LivermorePhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedComptonModel::SampleSecondaries(), G4LivermorePolarizedGammaConversionModel::SampleSecondaries(), G4LivermorePolarizedPhotoElectricModel::SampleSecondaries(), G4LivermorePolarizedRayleighModel::SampleSecondaries(), G4LivermoreRayleighModel::SampleSecondaries(), G4LowEPComptonModel::SampleSecondaries(), G4MuElecElasticModel::SampleSecondaries(), G4MuElecInelasticModel::SampleSecondaries(), G4PenelopeAnnihilationModel::SampleSecondaries(), G4PenelopeComptonModel::SampleSecondaries(), G4PenelopeGammaConversionModel::SampleSecondaries(), G4PenelopePhotoElectricModel::SampleSecondaries(), G4PenelopeRayleighModel::SampleSecondaries(), G4PolarizedAnnihilationModel::SampleSecondaries(), G4PolarizedComptonModel::SampleSecondaries(), G4BetheHeitlerModel::SampleSecondaries(), G4eCoulombScatteringModel::SampleSecondaries(), G4eeToTwoGammaModel::SampleSecondaries(), G4eSingleCoulombScatteringModel::SampleSecondaries(), G4hCoulombScatteringModel::SampleSecondaries(), G4HeatedKleinNishinaCompton::SampleSecondaries(), G4IonCoulombScatteringModel::SampleSecondaries(), G4KleinNishinaCompton::SampleSecondaries(), G4KleinNishinaModel::SampleSecondaries(), G4PairProductionRelModel::SampleSecondaries(), G4PEEffectFluoModel::SampleSecondaries(), G4PEEffectModel::SampleSecondaries(), and G4XrayRayleighModel::SampleSecondaries().

◆ UpdateStepForAtRest()

G4Step * G4ParticleChangeForGamma::UpdateStepForAtRest ( G4Step pStep)
virtual

Reimplemented from G4VParticleChange.

Definition at line 127 of file G4ParticleChangeForGamma.cc.

128{
130 pStep->SetStepLength( 0.0 );
131
134 }
135
136 return pStep;
137}
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 139 of file G4ParticleChangeForGamma.cc.

140{
141 G4StepPoint* pPostStepPoint = pStep->GetPostStepPoint();
142 pPostStepPoint->SetKineticEnergy( proposedKinEnergy );
143 pPostStepPoint->SetMomentumDirection( proposedMomentumDirection );
144 pPostStepPoint->SetPolarization( proposedPolarization );
145
147 pPostStepPoint->SetWeight( theParentWeight );
148 }
149
150 pStep->AddTotalEnergyDeposit( theLocalEnergyDeposit );
151 pStep->AddNonIonizingEnergyDeposit( theNonIonizingEnergyDeposit );
152 return pStep;
153}
void SetKineticEnergy(const G4double aValue)
void SetMomentumDirection(const G4ThreeVector &aValue)
void SetPolarization(const G4ThreeVector &aValue)

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