Geant4 9.6.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 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 61 of file G4PenelopeComptonModel.cc.

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

◆ ~G4PenelopeComptonModel()

G4PenelopeComptonModel::~G4PenelopeComptonModel ( )
virtual

Definition at line 89 of file G4PenelopeComptonModel.cc.

90{;}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 198 of file G4PenelopeComptonModel.cc.

204{
205 G4cout << "*** G4PenelopeComptonModel -- WARNING ***" << G4endl;
206 G4cout << "Penelope Compton model v2008 does not calculate cross section _per atom_ " << G4endl;
207 G4cout << "so the result is always zero. For physics values, please invoke " << G4endl;
208 G4cout << "GetCrossSectionPerVolume() or GetMeanFreePath() via the G4EmCalculator" << G4endl;
209 return 0;
210}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT 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 128 of file G4PenelopeComptonModel.cc.

133{
134 // Penelope model v2008 to calculate the Compton scattering cross section:
135 // D. Brusa et al., Nucl. Instrum. Meth. A 379 (1996) 167
136 //
137 // The cross section for Compton scattering is calculated according to the Klein-Nishina
138 // formula for energy > 5 MeV.
139 // For E < 5 MeV it is used a parametrization for the differential cross-section dSigma/dOmega,
140 // which is integrated numerically in cos(theta), G4PenelopeComptonModel::DifferentialCrossSection().
141 // The parametrization includes the J(p)
142 // distribution profiles for the atomic shells, that are tabulated from Hartree-Fock calculations
143 // from F. Biggs et al., At. Data Nucl. Data Tables 16 (1975) 201
144 //
145 if (verboseLevel > 3)
146 G4cout << "Calling CrossSectionPerVolume() of G4PenelopeComptonModel" << G4endl;
147 SetupForMaterial(p, material, energy);
148
149 //Retrieve the oscillator table for this material
150 G4PenelopeOscillatorTable* theTable = oscManager->GetOscillatorTableCompton(material);
151
152 G4double cs = 0;
153
154 if (energy < 5*MeV) //explicit calculation for E < 5 MeV
155 {
156 size_t numberOfOscillators = theTable->size();
157 for (size_t i=0;i<numberOfOscillators;i++)
158 {
159 G4PenelopeOscillator* theOsc = (*theTable)[i];
160 //sum contributions coming from each oscillator
161 cs += OscillatorTotalCrossSection(energy,theOsc);
162 }
163 }
164 else //use Klein-Nishina for E>5 MeV
165 cs = KleinNishinaCrossSection(energy,material);
166
167 //cross sections are in units of pi*classic_electr_radius^2
168 cs *= pi*classic_electr_radius*classic_electr_radius;
169
170 //Now, cs is the cross section *per molecule*, let's calculate the
171 //cross section per volume
172
173 G4double atomDensity = material->GetTotNbOfAtomsPerVolume();
174 G4double atPerMol = oscManager->GetAtomsPerMolecule(material);
175
176 if (verboseLevel > 3)
177 G4cout << "Material " << material->GetName() << " has " << atPerMol <<
178 "atoms per molecule" << G4endl;
179
180 G4double moleculeDensity = 0.;
181
182 if (atPerMol)
183 moleculeDensity = atomDensity/atPerMol;
184
185 G4double csvolume = cs*moleculeDensity;
186
187 if (verboseLevel > 2)
188 G4cout << "Compton mean free path at " << energy/keV << " keV for material " <<
189 material->GetName() << " = " << (1./csvolume)/mm << " mm" << G4endl;
190 return csvolume;
191}
std::vector< G4PenelopeOscillator * > G4PenelopeOscillatorTable
double G4double
Definition: G4Types.hh:64
G4double GetTotNbOfAtomsPerVolume() const
Definition: G4Material.hh:208
const G4String & GetName() const
Definition: G4Material.hh:177
G4double GetAtomsPerMolecule(const G4Material *)
Returns the total number of atoms per molecule.
G4PenelopeOscillatorTable * GetOscillatorTableCompton(const G4Material *)
virtual void SetupForMaterial(const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
Definition: G4VEmModel.cc:311
const G4double pi

◆ GetVerbosityLevel()

G4int G4PenelopeComptonModel::GetVerbosityLevel ( )
inline

Definition at line 97 of file G4PenelopeComptonModel.hh.

97{return verboseLevel;};

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 94 of file G4PenelopeComptonModel.cc.

96{
97 if (verboseLevel > 3)
98 G4cout << "Calling G4PenelopeComptonModel::Initialise()" << G4endl;
99
100 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
101 //Issue warning if the AtomicDeexcitation has not been declared
102 if (!fAtomDeexcitation)
103 {
104 G4cout << G4endl;
105 G4cout << "WARNING from G4PenelopeComptonModel " << G4endl;
106 G4cout << "Atomic de-excitation module is not instantiated, so there will not be ";
107 G4cout << "any fluorescence/Auger emission." << G4endl;
108 G4cout << "Please make sure this is intended" << G4endl;
109 }
110
111
112 if (verboseLevel > 0)
113 {
114 G4cout << "Penelope Compton model v2008 is initialized " << G4endl
115 << "Energy range: "
116 << LowEnergyLimit() / keV << " keV - "
117 << HighEnergyLimit() / GeV << " GeV";
118 }
119
120 if(isInitialised) return;
122 isInitialised = true;
123
124}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 214 of file G4PenelopeComptonModel.cc.

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

96{verboseLevel = lev;};

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4PenelopeComptonModel::fParticleChange
protected

Definition at line 101 of file G4PenelopeComptonModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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