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

#include <G4PenelopeComptonModel.hh>

+ Inheritance diagram for G4PenelopeComptonModel:

Public Member Functions

 G4PenelopeComptonModel (const G4ParticleDefinition *p=0, const G4String &processName="PenCompton")
 
virtual ~G4PenelopeComptonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double, G4double, G4double, G4double, G4double)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetVerbosityLevel (G4int lev)
 
G4int GetVerbosityLevel ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
const G4ParticleDefinitionfParticle
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 62 of file G4PenelopeComptonModel.hh.

Constructor & Destructor Documentation

◆ G4PenelopeComptonModel()

G4PenelopeComptonModel::G4PenelopeComptonModel ( const G4ParticleDefinition p = 0,
const G4String processName = "PenCompton" 
)

Definition at line 63 of file G4PenelopeComptonModel.cc.

66 isInitialised(false),fAtomDeexcitation(0),
67 oscManager(0)
68{
69 fIntrinsicLowEnergyLimit = 100.0*eV;
70 fIntrinsicHighEnergyLimit = 100.0*GeV;
71 SetHighEnergyLimit(fIntrinsicHighEnergyLimit);
72 //
74
75 if (part)
76 SetParticle(part);
77
78 verboseLevel= 0;
79 // Verbosity scale:
80 // 0 = nothing
81 // 1 = warning for energy non-conservation
82 // 2 = details of energy budget
83 // 3 = calculation of cross sections, file openings, sampling of atoms
84 // 4 = entering in methods
85
86 //Mark this model as "applicable" for atomic deexcitation
88
89 fTransitionManager = G4AtomicTransitionManager::Instance();
90}
static G4AtomicTransitionManager * Instance()
G4ParticleChangeForGamma * fParticleChange
const G4ParticleDefinition * fParticle
static G4PenelopeOscillatorManager * GetOscillatorManager()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4PenelopeComptonModel()

G4PenelopeComptonModel::~G4PenelopeComptonModel ( )
virtual

Definition at line 94 of file G4PenelopeComptonModel.cc.

95{;}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4PenelopeComptonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  ,
G4double  ,
G4double  ,
G4double  ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 252 of file G4PenelopeComptonModel.cc.

258{
259 G4cout << "*** G4PenelopeComptonModel -- WARNING ***" << G4endl;
260 G4cout << "Penelope Compton model v2008 does not calculate cross section _per atom_ " << G4endl;
261 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
262 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
263 return 0;
264}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ CrossSectionPerVolume()

G4double G4PenelopeComptonModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  kineticEnergy,
G4double  cutEnergy = 0.0,
G4double  maxEnergy = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 178 of file G4PenelopeComptonModel.cc.

183{
184 // Penelope model v2008 to calculate the Compton scattering cross section:
185 // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
186 //
187 // The cross section for Compton scattering is calculated according to the Klein-Nishina
188 // formula for energy > 5 MeV.
189 // For E < 5 MeV it is used a parametrization for the differential cross-section dSigma/dOmega,
190 // which is integrated numerically in cos(theta), G4PenelopeComptonModel::DifferentialCrossSection().
191 // The parametrization includes the J(p)
192 // distribution profiles for the atomic shells, that are tabulated from Hartree-Fock calculations
193 // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201
194 //
195 if (verboseLevel > 3)
196 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeComptonModel" << G4endl;
197 SetupForMaterial(p, material, energy);
198
199
200 G4double cs = 0;
201 //Force null cross-section if below the low-energy edge of the table
202 if (energy < LowEnergyLimit())
203 return cs;
204
205 //Retrieve the oscillator table for this material
206 G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
207
208 if (energy < 5*MeV) //explicit calculation for E < 5 MeV
209 {
210 size_t numberOfOscillators = theTable->size();
211 for (size_t i=0;i<numberOfOscillators;i++)
212 {
213 G4PenelopeOscillator* theOsc = (*theTable)[i];
214 //sum contributions coming from each oscillator
215 cs += OscillatorTotalCrossSection(energy,theOsc);
216 }
217 }
218 else //use Klein-Nishina for E>5 MeV
219 cs = KleinNishinaCrossSection(energy,material);
220
221 //cross sections are in units of pi*classic_electr_radius^2
222 cs *= pi*classic_electr_radius*classic_electr_radius;
223
224 //Now, cs is the cross section *per molecule*, let's calculate the
225 //cross section per volume
226
227 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
228 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
229
230 if (verboseLevel > 3)
231 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
232 "atoms per molecule" << G4endl;
233
234 G4double moleculeDensity = 0.;
235
236 if (atPerMol)
237 moleculeDensity = atomDensity/atPerMol;
238
239 G4double csvolume = cs*moleculeDensity;
240
241 if (verboseLevel > 2)
242 G4cout << "Compton mean free path at " << energy/keV << " keV for material " <<
243 material->GetName() << " = " << (1./csvolume)/mm << " mm" << G4endl;
244 return csvolume;
245}
std::vector< G4PenelopeOscillator * > G4PenelopeOscillatorTable
double G4double
Definition: G4Types.hh:83
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:207
const G4String & GetName() const
Definition: G4Material.hh:175
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.
G4PenelopeOscillatorTable * GetOscillatorTableCompton(const G4Material *)
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:449
G4double energy(const ThreeVector &p, const G4double m)
const G4double pi

◆ GetVerbosityLevel()

G4int G4PenelopeComptonModel::GetVerbosityLevel ( )
inline

Definition at line 99 of file G4PenelopeComptonModel.hh.

99{return verboseLevel;};

◆ Initialise()

void G4PenelopeComptonModel::Initialise ( const G4ParticleDefinition part,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 99 of file G4PenelopeComptonModel.cc.

101{
102 if (verboseLevel > 3)
103 G4cout << "Calling G4PenelopeComptonModel::Initialise()" << G4endl;
104
105 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
106 //Issue warning if the AtomicDeexcitation has not been declared
107 if (!fAtomDeexcitation)
108 {
109 G4cout << G4endl;
110 G4cout << "WARNING from G4PenelopeComptonModel " << G4endl;
111 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
112 G4cout << "any fluorescence/Auger emission." << G4endl;
113 G4cout << "Please make sure this is intended" << G4endl;
114 }
115
116 SetParticle(part);
117
118 if (IsMaster() && part == fParticle)
119 {
120
121 if (verboseLevel > 0)
122 {
123 G4cout << "Penelope Compton model v2008 is initialized " << G4endl
124 << "Energy range: "
125 << LowEnergyLimit() / keV << " keV - "
126 << HighEnergyLimit() / GeV << " GeV";
127 }
128 //Issue a warning, if the model is going to be used down to a
129 //energy which is outside the validity of the model itself
130 if (LowEnergyLimit() < fIntrinsicLowEnergyLimit)
131 {
133 ed << "Using the Penelope Compton model outside its intrinsic validity range. "
134 << G4endl;
135 ed << "-> LowEnergyLimit() in process = " << LowEnergyLimit()/keV << "keV " << G4endl;
136 ed << "-> Instrinsic low-energy limit = " << fIntrinsicLowEnergyLimit/keV << "keV "
137 << G4endl;
138 ed << "Result of the simulation have to be taken with care" << G4endl;
139 G4Exception("G4PenelopeComptonModel::Initialise()",
140 "em2100",JustWarning,ed);
141 }
142 }
143
144 if(isInitialised) return;
146 isInitialised = true;
147
148}
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4bool IsMaster() const
Definition: G4VEmModel.hh:736
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645

◆ InitialiseLocal()

void G4PenelopeComptonModel::InitialiseLocal ( const G4ParticleDefinition part,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 152 of file G4PenelopeComptonModel.cc.

154{
155 if (verboseLevel > 3)
156 G4cout << "Calling G4PenelopeComptonModel::InitialiseLocal()" << G4endl;
157
158 //
159 //Check that particle matches: one might have multiple master models (e.g.
160 //for e+ and e-).
161 //
162 if (part == fParticle)
163 {
164 //Get the const table pointers from the master to the workers
165 const G4PenelopeComptonModel* theModel =
166 static_cast<G4PenelopeComptonModel*> (masterModel);
167
168 //Same verbosity for all workers, as the master
169 verboseLevel = theModel->verboseLevel;
170 }
171
172 return;
173}

◆ SampleSecondaries()

void G4PenelopeComptonModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 268 of file G4PenelopeComptonModel.cc.

273{
274
275 // Penelope model v2008 to sample the Compton scattering final state.
276 // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
277 // The model determines also the original shell from which the electron is expelled,
278 // in order to produce fluorescence de-excitation (from G4DeexcitationManager)
279 //
280 // The final state for Compton scattering is calculated according to the Klein-Nishina
281 // formula for energy > 5 MeV. In this case, the Doppler broadening is negligible and
282 // one can assume that the target electron is at rest.
283 // For E < 5 MeV it is used the parametrization for the differential cross-section dSigma/dOmega,
284 // to sample the scattering angle and the energy of the emerging electron, which is
285 // G4PenelopeComptonModel::DifferentialCrossSection(). The rejection method is
286 // used to sample cos(theta). The efficiency increases monotonically with photon energy and is
287 // nearly independent on the Z; typical values are 35%, 80% and 95% for 1 keV, 1 MeV and 10 MeV,
288 // respectively.
289 // The parametrization includes the J(p) distribution profiles for the atomic shells, that are
290 // tabulated
291 // from Hartree-Fock calculations from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201.
292 // Doppler broadening is included.
293 //
294
295 if (verboseLevel > 3)
296 G4cout << "Calling SampleSecondaries() of G4PenelopeComptonModel" << G4endl;
297
298 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
299
300 // do nothing below the threshold
301 // should never get here because the XS is zero below the limit
302 if(photonEnergy0 < LowEnergyLimit())
303 return;
304
305 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
306 const G4Material* material = couple->GetMaterial();
307
308 G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
309
310 const G4int nmax = 64;
311 G4double rn[nmax]={0.0};
312 G4double pac[nmax]={0.0};
313
314 G4double S=0.0;
315 G4double epsilon = 0.0;
316 G4double cosTheta = 1.0;
317 G4double hartreeFunc = 0.0;
318 G4double oscStren = 0.0;
319 size_t numberOfOscillators = theTable->size();
320 size_t targetOscillator = 0;
321 G4double ionEnergy = 0.0*eV;
322
323 G4double ek = photonEnergy0/electron_mass_c2;
324 G4double ek2 = 2.*ek+1.0;
325 G4double eks = ek*ek;
326 G4double ek1 = eks-ek2-1.0;
327
328 G4double taumin = 1.0/ek2;
329 G4double a1 = G4Log(ek2);
330 G4double a2 = a1+2.0*ek*(1.0+ek)/(ek2*ek2);
331
332 G4double TST = 0;
333 G4double tau = 0.;
334
335 //If the incoming photon is above 5 MeV, the quicker approach based on the
336 //pure Klein-Nishina formula is used
337 if (photonEnergy0 > 5*MeV)
338 {
339 do{
340 do{
341 if ((a2*G4UniformRand()) < a1)
342 tau = std::pow(taumin,G4UniformRand());
343 else
344 tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
345 //rejection function
346 TST = (1.0+tau*(ek1+tau*(ek2+tau*eks)))/(eks*tau*(1.0+tau*tau));
347 }while (G4UniformRand()> TST);
348 epsilon=tau;
349 cosTheta = 1.0 - (1.0-tau)/(ek*tau);
350
351 //Target shell electrons
352 TST = oscManager->GetTotalZ(material)*G4UniformRand();
353 targetOscillator = numberOfOscillators-1; //last level
354 S=0.0;
355 G4bool levelFound = false;
356 for (size_t j=0;j<numberOfOscillators && !levelFound; j++)
357 {
358 S += (*theTable)[j]->GetOscillatorStrength();
359 if (S > TST)
360 {
361 targetOscillator = j;
362 levelFound = true;
363 }
364 }
365 //check whether the level is valid
366 ionEnergy = (*theTable)[targetOscillator]->GetIonisationEnergy();
367 }while((epsilon*photonEnergy0-photonEnergy0+ionEnergy) >0);
368 }
369 else //photonEnergy0 < 5 MeV
370 {
371 //Incoherent scattering function for theta=PI
372 G4double s0=0.0;
373 G4double pzomc=0.0;
374 G4double rni=0.0;
375 G4double aux=0.0;
376 for (size_t i=0;i<numberOfOscillators;i++)
377 {
378 ionEnergy = (*theTable)[i]->GetIonisationEnergy();
379 if (photonEnergy0 > ionEnergy)
380 {
381 G4double aux2 = photonEnergy0*(photonEnergy0-ionEnergy)*2.0;
382 hartreeFunc = (*theTable)[i]->GetHartreeFactor();
383 oscStren = (*theTable)[i]->GetOscillatorStrength();
384 pzomc = hartreeFunc*(aux2-electron_mass_c2*ionEnergy)/
385 (electron_mass_c2*std::sqrt(2.0*aux2+ionEnergy*ionEnergy));
386 if (pzomc > 0)
387 rni = 1.0-0.5*G4Exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*
388 (std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
389 else
390 rni = 0.5*G4Exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*
391 (std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
392 s0 += oscStren*rni;
393 }
394 }
395 //Sampling tau
396 G4double cdt1 = 0.;
397 do
398 {
399 if ((G4UniformRand()*a2) < a1)
400 tau = std::pow(taumin,G4UniformRand());
401 else
402 tau = std::sqrt(1.0+G4UniformRand()*(taumin*taumin-1.0));
403 cdt1 = (1.0-tau)/(ek*tau);
404 //Incoherent scattering function
405 S = 0.;
406 for (size_t i=0;i<numberOfOscillators;i++)
407 {
408 ionEnergy = (*theTable)[i]->GetIonisationEnergy();
409 if (photonEnergy0 > ionEnergy) //sum only on excitable levels
410 {
411 aux = photonEnergy0*(photonEnergy0-ionEnergy)*cdt1;
412 hartreeFunc = (*theTable)[i]->GetHartreeFactor();
413 oscStren = (*theTable)[i]->GetOscillatorStrength();
414 pzomc = hartreeFunc*(aux-electron_mass_c2*ionEnergy)/
415 (electron_mass_c2*std::sqrt(2.0*aux+ionEnergy*ionEnergy));
416 if (pzomc > 0)
417 rn[i] = 1.0-0.5*G4Exp(0.5-(std::sqrt(0.5)+std::sqrt(2.0)*pzomc)*
418 (std::sqrt(0.5)+std::sqrt(2.0)*pzomc));
419 else
420 rn[i] = 0.5*G4Exp(0.5-(std::sqrt(0.5)-std::sqrt(2.0)*pzomc)*
421 (std::sqrt(0.5)-std::sqrt(2.0)*pzomc));
422 S += oscStren*rn[i];
423 pac[i] = S;
424 }
425 else
426 pac[i] = S-1e-6;
427 }
428 //Rejection function
429 TST = S*(1.0+tau*(ek1+tau*(ek2+tau*eks)))/(eks*tau*(1.0+tau*tau));
430 }while ((G4UniformRand()*s0) > TST);
431
432 cosTheta = 1.0 - cdt1;
433 G4double fpzmax=0.0,fpz=0.0;
434 G4double A=0.0;
435 //Target electron shell
436 do
437 {
438 do
439 {
440 TST = S*G4UniformRand();
441 targetOscillator = numberOfOscillators-1; //last level
442 G4bool levelFound = false;
443 for (size_t i=0;i<numberOfOscillators && !levelFound;i++)
444 {
445 if (pac[i]>TST)
446 {
447 targetOscillator = i;
448 levelFound = true;
449 }
450 }
451 A = G4UniformRand()*rn[targetOscillator];
452 hartreeFunc = (*theTable)[targetOscillator]->GetHartreeFactor();
453 oscStren = (*theTable)[targetOscillator]->GetOscillatorStrength();
454 if (A < 0.5)
455 pzomc = (std::sqrt(0.5)-std::sqrt(0.5-G4Log(2.0*A)))/
456 (std::sqrt(2.0)*hartreeFunc);
457 else
458 pzomc = (std::sqrt(0.5-G4Log(2.0-2.0*A))-std::sqrt(0.5))/
459 (std::sqrt(2.0)*hartreeFunc);
460 } while (pzomc < -1);
461
462 // F(EP) rejection
463 G4double XQC = 1.0+tau*(tau-2.0*cosTheta);
464 G4double AF = std::sqrt(XQC)*(1.0+tau*(tau-cosTheta)/XQC);
465 if (AF > 0)
466 fpzmax = 1.0+AF*0.2;
467 else
468 fpzmax = 1.0-AF*0.2;
469 fpz = 1.0+AF*std::max(std::min(pzomc,0.2),-0.2);
470 }while ((fpzmax*G4UniformRand())>fpz);
471
472 //Energy of the scattered photon
473 G4double T = pzomc*pzomc;
474 G4double b1 = 1.0-T*tau*tau;
475 G4double b2 = 1.0-T*tau*cosTheta;
476 if (pzomc > 0.0)
477 epsilon = (tau/b1)*(b2+std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
478 else
479 epsilon = (tau/b1)*(b2-std::sqrt(std::abs(b2*b2-b1*(1.0-T))));
480 } //energy < 5 MeV
481
482 //Ok, the kinematics has been calculated.
483 G4double sinTheta = std::sqrt(1-cosTheta*cosTheta);
484 G4double phi = twopi * G4UniformRand() ;
485 G4double dirx = sinTheta * std::cos(phi);
486 G4double diry = sinTheta * std::sin(phi);
487 G4double dirz = cosTheta ;
488
489 // Update G4VParticleChange for the scattered photon
490 G4ThreeVector photonDirection1(dirx,diry,dirz);
491 photonDirection1.rotateUz(photonDirection0);
492 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
493
494 G4double photonEnergy1 = epsilon * photonEnergy0;
495
496 if (photonEnergy1 > 0.)
498 else
499 {
502 }
503
504 //Create scattered electron
505 G4double diffEnergy = photonEnergy0*(1-epsilon);
506 ionEnergy = (*theTable)[targetOscillator]->GetIonisationEnergy();
507
508 G4double Q2 =
509 photonEnergy0*photonEnergy0+photonEnergy1*(photonEnergy1-2.0*photonEnergy0*cosTheta);
510 G4double cosThetaE = 0.; //scattering angle for the electron
511
512 if (Q2 > 1.0e-12)
513 cosThetaE = (photonEnergy0-photonEnergy1*cosTheta)/std::sqrt(Q2);
514 else
515 cosThetaE = 1.0;
516 G4double sinThetaE = std::sqrt(1-cosThetaE*cosThetaE);
517
518 //Now, try to handle fluorescence
519 //Notice: merged levels are indicated with Z=0 and flag=30
520 G4int shFlag = (*theTable)[targetOscillator]->GetShellFlag();
521 G4int Z = (G4int) (*theTable)[targetOscillator]->GetParentZ();
522
523 //initialize here, then check photons created by Atomic-Deexcitation, and the final state e-
524 G4double bindingEnergy = 0.*eV;
525 const G4AtomicShell* shell = 0;
526
527 //Real level
528 if (Z > 0 && shFlag<30)
529 {
530 shell = fTransitionManager->Shell(Z,shFlag-1);
531 bindingEnergy = shell->BindingEnergy();
532 }
533
534 G4double ionEnergyInPenelopeDatabase = ionEnergy;
535 //protection against energy non-conservation
536 ionEnergy = std::max(bindingEnergy,ionEnergyInPenelopeDatabase);
537
538 //subtract the excitation energy. If not emitted by fluorescence
539 //the ionization energy is deposited as local energy deposition
540 G4double eKineticEnergy = diffEnergy - ionEnergy;
541 G4double localEnergyDeposit = ionEnergy;
542 G4double energyInFluorescence = 0.; //testing purposes only
543 G4double energyInAuger = 0; //testing purposes
544
545 if (eKineticEnergy < 0)
546 {
547 //It means that there was some problem/mismatch between the two databases.
548 //Try to make it work
549 //In this case available Energy (diffEnergy) < ionEnergy
550 //Full residual energy is deposited locally
551 localEnergyDeposit = diffEnergy;
552 eKineticEnergy = 0.0;
553 }
554
555 //the local energy deposit is what remains: part of this may be spent for fluorescence.
556 //Notice: shell might be nullptr (invalid!) if shFlag=30. Must be protected
557 //Now, take care of fluorescence, if required
558 if (fAtomDeexcitation && shell)
559 {
560 G4int index = couple->GetIndex();
561 if (fAtomDeexcitation->CheckDeexcitationActiveRegion(index))
562 {
563 size_t nBefore = fvect->size();
564 fAtomDeexcitation->GenerateParticles(fvect,shell,Z,index);
565 size_t nAfter = fvect->size();
566
567 if (nAfter > nBefore) //actual production of fluorescence
568 {
569 for (size_t j=nBefore;j<nAfter;j++) //loop on products
570 {
571 G4double itsEnergy = ((*fvect)[j])->GetKineticEnergy();
572 if (itsEnergy < localEnergyDeposit) // valid secondary, generate it
573 {
574 localEnergyDeposit -= itsEnergy;
575 if (((*fvect)[j])->GetParticleDefinition() == G4Gamma::Definition())
576 energyInFluorescence += itsEnergy;
577 else if (((*fvect)[j])->GetParticleDefinition() ==
579 energyInAuger += itsEnergy;
580 }
581 else //invalid secondary: takes more than the available energy: delete it
582 {
583 delete (*fvect)[j];
584 (*fvect)[j] = nullptr;
585 }
586 }
587 }
588
589 }
590 }
591
592
593 /*
594 if(DeexcitationFlag() && Z > 5 && shellId>0) {
595
596 const G4ProductionCutsTable* theCoupleTable=
597 G4ProductionCutsTable::GetProductionCutsTable();
598
599 size_t index = couple->GetIndex();
600 G4double cutg = (*(theCoupleTable->GetEnergyCutsVector(0)))[index];
601 G4double cute = (*(theCoupleTable->GetEnergyCutsVector(1)))[index];
602
603 // Generation of fluorescence
604 // Data in EADL are available only for Z > 5
605 // Protection to avoid generating photons in the unphysical case of
606 // shell binding energy > photon energy
607 if (localEnergyDeposit > cutg || localEnergyDeposit > cute)
608 {
609 G4DynamicParticle* aPhoton;
610 deexcitationManager.SetCutForSecondaryPhotons(cutg);
611 deexcitationManager.SetCutForAugerElectrons(cute);
612
613 photonVector = deexcitationManager.GenerateParticles(Z,shellId);
614 if(photonVector)
615 {
616 size_t nPhotons = photonVector->size();
617 for (size_t k=0; k<nPhotons; k++)
618 {
619 aPhoton = (*photonVector)[k];
620 if (aPhoton)
621 {
622 G4double itsEnergy = aPhoton->GetKineticEnergy();
623 G4bool keepIt = false;
624 if (itsEnergy <= localEnergyDeposit)
625 {
626 //check if good!
627 if(aPhoton->GetDefinition() == G4Gamma::Gamma()
628 && itsEnergy >= cutg)
629 {
630 keepIt = true;
631 energyInFluorescence += itsEnergy;
632 }
633 if (aPhoton->GetDefinition() == G4Electron::Electron() &&
634 itsEnergy >= cute)
635 {
636 energyInAuger += itsEnergy;
637 keepIt = true;
638 }
639 }
640 //good secondary, register it
641 if (keepIt)
642 {
643 localEnergyDeposit -= itsEnergy;
644 fvect->push_back(aPhoton);
645 }
646 else
647 {
648 delete aPhoton;
649 (*photonVector)[k] = 0;
650 }
651 }
652 }
653 delete photonVector;
654 }
655 }
656 }
657 */
658
659
660 //Always produce explicitly the electron
662
663 G4double xEl = sinThetaE * std::cos(phi+pi);
664 G4double yEl = sinThetaE * std::sin(phi+pi);
665 G4double zEl = cosThetaE;
666 G4ThreeVector eDirection(xEl,yEl,zEl); //electron direction
667 eDirection.rotateUz(photonDirection0);
669 eDirection,eKineticEnergy) ;
670 fvect->push_back(electron);
671
672
673 if (localEnergyDeposit < 0) //Should not be: issue a G4Exception (warning)
674 {
675 G4Exception("G4PenelopeComptonModel::SampleSecondaries()",
676 "em2099",JustWarning,"WARNING: Negative local energy deposit");
677 localEnergyDeposit=0.;
678 }
679 fParticleChange->ProposeLocalEnergyDeposit(localEnergyDeposit);
680
681 G4double electronEnergy = 0.;
682 if (electron)
683 electronEnergy = eKineticEnergy;
684 if (verboseLevel > 1)
685 {
686 G4cout << "-----------------------------------------------------------" << G4endl;
687 G4cout << "Energy balance from G4PenelopeCompton" << G4endl;
688 G4cout << "Incoming photon energy: " << photonEnergy0/keV << " keV" << G4endl;
689 G4cout << "-----------------------------------------------------------" << G4endl;
690 G4cout << "Scattered photon: " << photonEnergy1/keV << " keV" << G4endl;
691 G4cout << "Scattered electron " << electronEnergy/keV << " keV" << G4endl;
692 if (energyInFluorescence)
693 G4cout << "Fluorescence x-rays: " << energyInFluorescence/keV << " keV" << G4endl;
694 if (energyInAuger)
695 G4cout << "Auger electrons: " << energyInAuger/keV << " keV" << G4endl;
696 G4cout << "Local energy deposit " << localEnergyDeposit/keV << " keV" << G4endl;
697 G4cout << "Total final state: " << (photonEnergy1+electronEnergy+energyInFluorescence+
698 localEnergyDeposit+energyInAuger)/keV <<
699 " keV" << G4endl;
700 G4cout << "-----------------------------------------------------------" << G4endl;
701 }
702 if (verboseLevel > 0)
703 {
704 G4double energyDiff = std::fabs(photonEnergy1+
705 electronEnergy+energyInFluorescence+
706 localEnergyDeposit+energyInAuger-photonEnergy0);
707 if (energyDiff > 0.05*keV)
708 G4cout << "Warning from G4PenelopeCompton: problem with energy conservation: " <<
709 (photonEnergy1+electronEnergy+energyInFluorescence+energyInAuger+localEnergyDeposit)/keV <<
710 " keV (final) vs. " <<
711 photonEnergy0/keV << " keV (initial)" << G4endl;
712 }
713}
double epsilon(double density, double temperature)
double S(double temp)
double A(double temperature)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
G4double BindingEnergy() const
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
static G4Electron * Definition()
Definition: G4Electron.cc:48
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Definition()
Definition: G4Gamma.cc:48
const G4Material * GetMaterial() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetTotalZ(const G4Material *)
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SetVerbosityLevel()

void G4PenelopeComptonModel::SetVerbosityLevel ( G4int  lev)
inline

Definition at line 98 of file G4PenelopeComptonModel.hh.

98{verboseLevel = lev;};

Member Data Documentation

◆ fParticle

const G4ParticleDefinition* G4PenelopeComptonModel::fParticle
protected

Definition at line 104 of file G4PenelopeComptonModel.hh.

Referenced by Initialise(), and InitialiseLocal().

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeComptonModel::fParticleChange
protected

Definition at line 103 of file G4PenelopeComptonModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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