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

#include <G4LivermorePolarizedComptonModel.hh>

+ Inheritance diagram for G4LivermorePolarizedComptonModel:

Public Member Functions

 G4LivermorePolarizedComptonModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedCompton")
 
virtual ~G4LivermorePolarizedComptonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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
 

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 *)
 
- 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
 

Detailed Description

Definition at line 52 of file G4LivermorePolarizedComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedComptonModel()

G4LivermorePolarizedComptonModel::G4LivermorePolarizedComptonModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermorePolarizedCompton" 
)

Definition at line 69 of file G4LivermorePolarizedComptonModel.cc.

70 :G4VEmModel(nam),isInitialised(false)
71{
72 verboseLevel= 1;
73 // Verbosity scale:
74 // 0 = nothing
75 // 1 = warning for energy non-conservation
76 // 2 = details of energy budget
77 // 3 = calculation of cross sections, file openings, sampling of atoms
78 // 4 = entering in methods
79
80 if( verboseLevel>1 )
81 G4cout << "Livermore Polarized Compton is constructed " << G4endl;
82
83 //Mark this model as "applicable" for atomic deexcitation
85
86 fParticleChange = 0;
87 fAtomDeexcitation = 0;
88}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4LivermorePolarizedComptonModel()

G4LivermorePolarizedComptonModel::~G4LivermorePolarizedComptonModel ( )
virtual

Definition at line 92 of file G4LivermorePolarizedComptonModel.cc.

93{
94 if(IsMaster()) {
95 delete shellData;
96 shellData = 0;
97 delete profileData;
98 profileData = 0;
99 delete scatterFunctionData;
100 scatterFunctionData = 0;
101 for(G4int i=0; i<maxZ; ++i) {
102 if(data[i]) {
103 delete data[i];
104 data[i] = 0;
105 }
106 }
107 }
108}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LivermorePolarizedComptonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 252 of file G4LivermorePolarizedComptonModel.cc.

257{
258 if (verboseLevel > 3)
259 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermorePolarizedComptonModel" << G4endl;
260
261 G4double cs = 0.0;
262
263 if (GammaEnergy < LowEnergyLimit())
264 return 0.0;
265
266 G4int intZ = G4lrint(Z);
267 if(intZ < 1 || intZ > maxZ) { return cs; }
268
269 G4LPhysicsFreeVector* pv = data[intZ];
270
271 // if element was not initialised
272 // do initialisation safely for MT mode
273 if(!pv)
274 {
275 InitialiseForElement(0, intZ);
276 pv = data[intZ];
277 if(!pv) { return cs; }
278 }
279
280 G4int n = pv->GetVectorLength() - 1;
281 G4double e1 = pv->Energy(0);
282 G4double e2 = pv->Energy(n);
283
284 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
285 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
286 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
287
288 return cs;
289
290}
double G4double
Definition: G4Types.hh:83
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double Energy(std::size_t index) const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
int G4lrint(double ad)
Definition: templates.hh:134

◆ Initialise()

void G4LivermorePolarizedComptonModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 112 of file G4LivermorePolarizedComptonModel.cc.

114{
115 if (verboseLevel > 1)
116 G4cout << "Calling G4LivermorePolarizedComptonModel::Initialise()" << G4endl;
117
118 // Initialise element selector
119
120 if(IsMaster()) {
121
122 // Access to elements
123
124 char* path = std::getenv("G4LEDATA");
125
126 G4ProductionCutsTable* theCoupleTable =
128
129 G4int numOfCouples = theCoupleTable->GetTableSize();
130
131 for(G4int i=0; i<numOfCouples; ++i) {
132 const G4Material* material =
133 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
134 const G4ElementVector* theElementVector = material->GetElementVector();
135 G4int nelm = material->GetNumberOfElements();
136
137 for (G4int j=0; j<nelm; ++j) {
138 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
139 if(Z < 1) { Z = 1; }
140 else if(Z > maxZ){ Z = maxZ; }
141
142 if( (!data[Z]) ) { ReadData(Z, path); }
143 }
144 }
145
146 // For Doppler broadening
147 if(!shellData) {
148 shellData = new G4ShellData();
149 shellData->SetOccupancyData();
150 G4String file = "/doppler/shell-doppler";
151 shellData->LoadData(file);
152 }
153 if(!profileData) { profileData = new G4DopplerProfile(); }
154
155 // Scattering Function
156
157 if(!scatterFunctionData)
158 {
159
160 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
161 G4String scatterFile = "comp/ce-sf-";
162 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
163 scatterFunctionData->LoadData(scatterFile);
164 }
165
166 InitialiseElementSelectors(particle, cuts);
167 }
168
169 if (verboseLevel > 2) {
170 G4cout << "Loaded cross section files" << G4endl;
171 }
172
173 if( verboseLevel>1 ) {
174 G4cout << "G4LivermoreComptonModel is initialized " << G4endl
175 << "Energy range: "
176 << LowEnergyLimit() / eV << " eV - "
177 << HighEnergyLimit() / GeV << " GeV"
178 << G4endl;
179 }
180 //
181 if(isInitialised) { return; }
182
183 fParticleChange = GetParticleChangeForGamma();
184 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
185 isInitialised = true;
186}
std::vector< G4Element * > G4ElementVector
virtual G4bool LoadData(const G4String &fileName)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetOccupancyData()
Definition: G4ShellData.hh:69
void LoadData(const G4String &fileName)
Definition: G4ShellData.cc:233
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

void G4LivermorePolarizedComptonModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 874 of file G4LivermorePolarizedComptonModel.cc.

876{
877 G4AutoLock l(&LivermorePolarizedComptonModelMutex);
878 // G4cout << "G4LivermoreComptonModel::InitialiseForElement Z= "
879 // << Z << G4endl;
880 if(!data[Z]) { ReadData(Z); }
881 l.unlock();
882}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4LivermorePolarizedComptonModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 189 of file G4LivermorePolarizedComptonModel.cc.

191{
193}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 294 of file G4LivermorePolarizedComptonModel.cc.

299{
300 // The scattered gamma energy is sampled according to Klein - Nishina formula.
301 // The random number techniques of Butcher & Messel are used (Nuc Phys 20(1960),15).
302 // GEANT4 internal units
303 //
304 // Note : Effects due to binding of atomic electrons are negliged.
305
306 if (verboseLevel > 3)
307 G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedComptonModel" << G4endl;
308
309 G4double gammaEnergy0 = aDynamicGamma->GetKineticEnergy();
310
311 // do nothing below the threshold
312 // should never get here because the XS is zero below the limit
313 if (gammaEnergy0 < LowEnergyLimit())
314 return ;
315
316
317 G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
318
319 // Protection: a polarisation parallel to the
320 // direction causes problems;
321 // in that case find a random polarization
322
323 G4ThreeVector gammaDirection0 = aDynamicGamma->GetMomentumDirection();
324
325 // Make sure that the polarization vector is perpendicular to the
326 // gamma direction. If not
327
328 if(!(gammaPolarization0.isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.mag()==0))
329 { // only for testing now
330 gammaPolarization0 = GetRandomPolarization(gammaDirection0);
331 }
332 else
333 {
334 if ( gammaPolarization0.howOrthogonal(gammaDirection0) != 0)
335 {
336 gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
337 }
338 }
339
340 // End of Protection
341
342 G4double E0_m = gammaEnergy0 / electron_mass_c2 ;
343
344 // Select randomly one element in the current material
345 //G4int Z = crossSectionHandler->SelectRandomAtom(couple,gammaEnergy0);
346 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
347 const G4Element* elm = SelectRandomAtom(couple,particle,gammaEnergy0);
348 G4int Z = (G4int)elm->GetZ();
349
350 // Sample the energy and the polarization of the scattered photon
351
352 G4double epsilon, epsilonSq, onecost, sinThetaSqr, greject ;
353
354 G4double epsilon0Local = 1./(1. + 2*E0_m);
355 G4double epsilon0Sq = epsilon0Local*epsilon0Local;
356 G4double alpha1 = - G4Log(epsilon0Local);
357 G4double alpha2 = 0.5*(1.- epsilon0Sq);
358
359 G4double wlGamma = h_Planck*c_light/gammaEnergy0;
360 G4double gammaEnergy1;
361 G4ThreeVector gammaDirection1;
362
363 do {
364 if ( alpha1/(alpha1+alpha2) > G4UniformRand() )
365 {
366 epsilon = G4Exp(-alpha1*G4UniformRand());
367 epsilonSq = epsilon*epsilon;
368 }
369 else
370 {
371 epsilonSq = epsilon0Sq + (1.- epsilon0Sq)*G4UniformRand();
372 epsilon = std::sqrt(epsilonSq);
373 }
374
375 onecost = (1.- epsilon)/(epsilon*E0_m);
376 sinThetaSqr = onecost*(2.-onecost);
377
378 // Protection
379 if (sinThetaSqr > 1.)
380 {
381 G4cout
382 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
383 << "sin(theta)**2 = "
384 << sinThetaSqr
385 << "; set to 1"
386 << G4endl;
387 sinThetaSqr = 1.;
388 }
389 if (sinThetaSqr < 0.)
390 {
391 G4cout
392 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
393 << "sin(theta)**2 = "
394 << sinThetaSqr
395 << "; set to 0"
396 << G4endl;
397 sinThetaSqr = 0.;
398 }
399 // End protection
400
401 G4double x = std::sqrt(onecost/2.) / (wlGamma/cm);;
402 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
403 greject = (1. - epsilon*sinThetaSqr/(1.+ epsilonSq))*scatteringFunction;
404
405 } while(greject < G4UniformRand()*Z);
406
407
408 // ****************************************************
409 // Phi determination
410 // ****************************************************
411
412 G4double phi = SetPhi(epsilon,sinThetaSqr);
413
414 //
415 // scattered gamma angles. ( Z - axis along the parent gamma)
416 //
417
418 G4double cosTheta = 1. - onecost;
419
420 // Protection
421
422 if (cosTheta > 1.)
423 {
424 G4cout
425 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
426 << "cosTheta = "
427 << cosTheta
428 << "; set to 1"
429 << G4endl;
430 cosTheta = 1.;
431 }
432 if (cosTheta < -1.)
433 {
434 G4cout
435 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
436 << "cosTheta = "
437 << cosTheta
438 << "; set to -1"
439 << G4endl;
440 cosTheta = -1.;
441 }
442 // End protection
443
444
445 G4double sinTheta = std::sqrt (sinThetaSqr);
446
447 // Protection
448 if (sinTheta > 1.)
449 {
450 G4cout
451 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
452 << "sinTheta = "
453 << sinTheta
454 << "; set to 1"
455 << G4endl;
456 sinTheta = 1.;
457 }
458 if (sinTheta < -1.)
459 {
460 G4cout
461 << " -- Warning -- G4LivermorePolarizedComptonModel::SampleSecondaries "
462 << "sinTheta = "
463 << sinTheta
464 << "; set to -1"
465 << G4endl;
466 sinTheta = -1.;
467 }
468 // End protection
469
470
471 G4double dirx = sinTheta*std::cos(phi);
472 G4double diry = sinTheta*std::sin(phi);
473 G4double dirz = cosTheta ;
474
475
476 // oneCosT , eom
477
478 // Doppler broadening - Method based on:
479 // Y. Namito, S. Ban and H. Hirayama,
480 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon Into the EGS4 Code"
481 // NIM A 349, pp. 489-494, 1994
482
483 // Maximum number of sampling iterations
484
485 static G4int maxDopplerIterations = 1000;
486 G4double bindingE = 0.;
487 G4double photonEoriginal = epsilon * gammaEnergy0;
488 G4double photonE = -1.;
489 G4int iteration = 0;
490 G4double eMax = gammaEnergy0;
491
492 G4int shellIdx = 0;
493
494 if (verboseLevel > 3) {
495 G4cout << "Started loop to sample broading" << G4endl;
496 }
497
498 do
499 {
500 iteration++;
501 // Select shell based on shell occupancy
502 shellIdx = shellData->SelectRandomShell(Z);
503 bindingE = shellData->BindingEnergy(Z,shellIdx);
504
505 if (verboseLevel > 3) {
506 G4cout << "Shell ID= " << shellIdx
507 << " Ebind(keV)= " << bindingE/keV << G4endl;
508 }
509 eMax = gammaEnergy0 - bindingE;
510
511 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units)
512 G4double pSample = profileData->RandomSelectMomentum(Z,shellIdx);
513
514 if (verboseLevel > 3) {
515 G4cout << "pSample= " << pSample << G4endl;
516 }
517 // Rescale from atomic units
518 G4double pDoppler = pSample * fine_structure_const;
519 G4double pDoppler2 = pDoppler * pDoppler;
520 G4double var2 = 1. + onecost * E0_m;
521 G4double var3 = var2*var2 - pDoppler2;
522 G4double var4 = var2 - pDoppler2 * cosTheta;
523 G4double var = var4*var4 - var3 + pDoppler2 * var3;
524 if (var > 0.)
525 {
526 G4double varSqrt = std::sqrt(var);
527 G4double scale = gammaEnergy0 / var3;
528 // Random select either root
529 if (G4UniformRand() < 0.5) photonE = (var4 - varSqrt) * scale;
530 else photonE = (var4 + varSqrt) * scale;
531 }
532 else
533 {
534 photonE = -1.;
535 }
536 } while ( iteration <= maxDopplerIterations &&
537 (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
538 //while (iteration <= maxDopplerIterations && photonE > eMax); ???
539
540
541 // End of recalculation of photon energy with Doppler broadening
542 // Revert to original if maximum number of iterations threshold has been reached
543 if (iteration >= maxDopplerIterations)
544 {
545 photonE = photonEoriginal;
546 bindingE = 0.;
547 }
548
549 gammaEnergy1 = photonE;
550
551 //
552 // update G4VParticleChange for the scattered photon
553 //
554
555
556
557 // gammaEnergy1 = epsilon*gammaEnergy0;
558
559
560 // New polarization
561
562 G4ThreeVector gammaPolarization1 = SetNewPolarization(epsilon,
563 sinThetaSqr,
564 phi,
565 cosTheta);
566
567 // Set new direction
568 G4ThreeVector tmpDirection1( dirx,diry,dirz );
569 gammaDirection1 = tmpDirection1;
570
571 // Change reference frame.
572
573 SystemOfRefChange(gammaDirection0,gammaDirection1,
574 gammaPolarization0,gammaPolarization1);
575
576 if (gammaEnergy1 > 0.)
577 {
578 fParticleChange->SetProposedKineticEnergy( gammaEnergy1 ) ;
579 fParticleChange->ProposeMomentumDirection( gammaDirection1 );
580 fParticleChange->ProposePolarization( gammaPolarization1 );
581 }
582 else
583 {
584 gammaEnergy1 = 0.;
585 fParticleChange->SetProposedKineticEnergy(0.) ;
586 fParticleChange->ProposeTrackStatus(fStopAndKill);
587 }
588
589 //
590 // kinematic of the scattered electron
591 //
592
593 G4double ElecKineEnergy = gammaEnergy0 - gammaEnergy1 -bindingE;
594
595 // SI -protection against negative final energy: no e- is created
596 // like in G4LivermoreComptonModel.cc
597 if(ElecKineEnergy < 0.0) {
598 fParticleChange->ProposeLocalEnergyDeposit(gammaEnergy0 - gammaEnergy1);
599 return;
600 }
601
602 // SI - Removed range test
603
604 G4double ElecMomentum = std::sqrt(ElecKineEnergy*(ElecKineEnergy+2.*electron_mass_c2));
605
606 G4ThreeVector ElecDirection((gammaEnergy0 * gammaDirection0 -
607 gammaEnergy1 * gammaDirection1) * (1./ElecMomentum));
608
609 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),ElecDirection.unit(),ElecKineEnergy) ;
610 fvect->push_back(dp);
611
612 // sample deexcitation
613 //
614
615 if (verboseLevel > 3) {
616 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
617 }
618
619 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
620 G4int index = couple->GetIndex();
621 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
622 size_t nbefore = fvect->size();
624 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
625 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
626 size_t nafter = fvect->size();
627 if(nafter > nbefore) {
628 for (size_t i=nbefore; i<nafter; ++i) {
629 //Check if there is enough residual energy
630 if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
631 {
632 //Ok, this is a valid secondary: keep it
633 bindingE -= ((*fvect)[i])->GetKineticEnergy();
634 }
635 else
636 {
637 //Invalid secondary: not enough energy to create it!
638 //Keep its energy in the local deposit
639 delete (*fvect)[i];
640 (*fvect)[i]=0;
641 }
642 }
643 }
644 }
645 }
646 //This should never happen
647 if(bindingE < 0.0)
648 G4Exception("G4LivermoreComptonModel::SampleSecondaries()",
649 "em2050",FatalException,"Negative local energy deposit");
650
651 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
652
653}
G4AtomicShellEnumerator
double epsilon(double density, double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:233
double mag() const
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:215
virtual G4double FindValue(G4double x, G4int componentId=0) const
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4ThreeVector & GetPolarization() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4double GetZ() const
Definition: G4Element.hh:130
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposePolarization(const G4ThreeVector &dir)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
Definition: G4ShellData.cc:165
G4int SelectRandomShell(G4int Z) const
Definition: G4ShellData.cc:362
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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