G4UAtomicDeexcitation Class Reference

#include <G4UAtomicDeexcitation.hh>

Inheritance diagram for G4UAtomicDeexcitation:

G4VAtomDeexcitation

Public Member Functions

 G4UAtomicDeexcitation ()
virtual ~G4UAtomicDeexcitation ()
virtual void InitialiseForNewRun ()
virtual void InitialiseForExtraAtom (G4int Z)
void SetCutForSecondaryPhotons (G4double cut)
void SetCutForAugerElectrons (G4double cut)
virtual const G4AtomicShellGetAtomicShell (G4int Z, G4AtomicShellEnumerator shell)
virtual void GenerateParticles (std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4double gammaCut, G4double eCut)
virtual G4double GetShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)
virtual G4double ComputeShellIonisationCrossSectionPerAtom (const G4ParticleDefinition *, G4int Z, G4AtomicShellEnumerator shell, G4double kinE, const G4Material *mat=0)

Detailed Description

Definition at line 62 of file G4UAtomicDeexcitation.hh.


Constructor & Destructor Documentation

G4UAtomicDeexcitation::G4UAtomicDeexcitation (  ) 

Definition at line 75 of file G4UAtomicDeexcitation.cc.

References G4Electron::Electron(), G4LossTableManager::EmCorrections(), G4LossTableManager::Instance(), and G4Positron::Positron().

00075                                             :
00076   G4VAtomDeexcitation("UAtomDeexcitation"),
00077   minGammaEnergy(DBL_MAX), 
00078   minElectronEnergy(DBL_MAX),
00079   emcorr(0)
00080 {
00081   PIXEshellCS    = 0;
00082   ePIXEshellCS   = 0;
00083   emcorr = G4LossTableManager::Instance()->EmCorrections();
00084   theElectron = G4Electron::Electron();
00085   thePositron = G4Positron::Positron();
00086   transitionManager = 0;
00087   anaPIXEshellCS = 0;
00088 }

G4UAtomicDeexcitation::~G4UAtomicDeexcitation (  )  [virtual]

Definition at line 90 of file G4UAtomicDeexcitation.cc.

00091 {
00092   delete PIXEshellCS;
00093   delete anaPIXEshellCS;
00094   delete ePIXEshellCS;
00095 }


Member Function Documentation

G4double G4UAtomicDeexcitation::ComputeShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
) [virtual]

Implements G4VAtomDeexcitation.

Definition at line 390 of file G4UAtomicDeexcitation.cc.

References GetShellIonisationCrossSectionPerAtom().

00396 {
00397   return GetShellIonisationCrossSectionPerAtom(p,Z,shell,kinE,mat);
00398 }

void G4UAtomicDeexcitation::GenerateParticles ( std::vector< G4DynamicParticle * > *  secVect,
const G4AtomicShell ,
G4int  Z,
G4double  gammaCut,
G4double  eCut 
) [virtual]

Implements G4VAtomDeexcitation.

Definition at line 233 of file G4UAtomicDeexcitation.cc.

References G4Exception(), JustWarning, and G4AtomicShell::ShellId().

00239 {
00240 
00241   // Defined initial conditions
00242   G4int givenShellId = atomicShell->ShellId();
00243   //G4cout << "generating particles for vacancy in shellId: " << givenShellId << G4endl; // debug
00244   minGammaEnergy = gammaCut;
00245   minElectronEnergy = eCut;
00246 
00247   // V.I. short-cut
00248   //  if(!IsAugerActive()) {  minElectronEnergy = DBL_MAX; }
00249 
00250   // generation secondaries
00251   G4DynamicParticle* aParticle=0;
00252   G4int provShellId = 0;
00253   G4int counter = 0;
00254   
00255   // let's check that 5<Z<100
00256 
00257   if (Z>5 && Z<100) {
00258 
00259   // The aim of this loop is to generate more than one fluorecence photon 
00260   // from the same ionizing event 
00261     do
00262       {
00263         if (counter == 0) 
00264           // First call to GenerateParticles(...):
00265           // givenShellId is given by the process
00266           {
00267             provShellId = SelectTypeOfTransition(Z, givenShellId);
00268 
00269             if  ( provShellId >0) 
00270               {
00271                 aParticle = GenerateFluorescence(Z,givenShellId,provShellId);
00272                 //if (aParticle != 0) { G4cout << "****FLUO!_1**** " << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl ;} //debug  
00273               }
00274             else if ( provShellId == -1)
00275               {
00276                 //              G4cout << "Try to generate Auger 1" << G4endl; //debug
00277                 aParticle = GenerateAuger(Z, givenShellId);
00278                 //              if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
00279               }
00280             else
00281               {
00282                 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
00283               }
00284           }
00285         else 
00286           // Following calls to GenerateParticles(...):
00287           // newShellId is given by GenerateFluorescence(...)
00288           {
00289             provShellId = SelectTypeOfTransition(Z,newShellId);
00290             if  (provShellId >0)
00291               {
00292                 aParticle = GenerateFluorescence(Z,newShellId,provShellId);
00293                 //if (aParticle != 0) { G4cout << "****FLUO!_2****" << aParticle->GetParticleDefinition()->GetParticleType() << " " << aParticle->GetKineticEnergy()/keV << G4endl;} //debug
00294               }
00295             else if ( provShellId == -1)
00296               {
00297                 //              G4cout << "Try to generate Auger 2" << G4endl; //debug
00298                 aParticle = GenerateAuger(Z, newShellId);
00299                 //              if (aParticle != 0) { G4cout << "****AUGER!****" << G4endl;} //debug
00300               }
00301             else
00302               {
00303                 G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0002",JustWarning, "Energy deposited locally");
00304               }
00305           }
00306         counter++;
00307         if (aParticle != 0) 
00308           {
00309             vectorOfParticles->push_back(aParticle);
00310             //G4cout << "Deexcitation Occurred!" << G4endl; //debug
00311           }
00312         else {provShellId = -2;}
00313       }  
00314     while (provShellId > -2); 
00315   }
00316   else
00317     {
00318       G4Exception("G4UAtomicDeexcitation::GenerateParticles()","de0001",JustWarning, "Energy deposited locally");
00319     }
00320   
00321   //G4cout << "---------FATTO!---------" << G4endl; //debug 
00322 
00323 }

const G4AtomicShell * G4UAtomicDeexcitation::GetAtomicShell ( G4int  Z,
G4AtomicShellEnumerator  shell 
) [virtual]

Implements G4VAtomDeexcitation.

Definition at line 228 of file G4UAtomicDeexcitation.cc.

References G4AtomicTransitionManager::Shell().

00229 {
00230   return transitionManager->Shell(Z, size_t(shell));
00231 }

G4double G4UAtomicDeexcitation::GetShellIonisationCrossSectionPerAtom ( const G4ParticleDefinition ,
G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  kinE,
const G4Material mat = 0 
) [virtual]

Implements G4VAtomDeexcitation.

Definition at line 326 of file G4UAtomicDeexcitation.cc.

References G4VhShellCrossSection::CrossSection(), G4EmCorrections::EffectiveChargeSquareRatio(), G4AtomicShells::GetNumberOfShells(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), and G4ParticleDefinition::GetPDGMass().

Referenced by ComputeShellIonisationCrossSectionPerAtom().

00332 {
00333 
00334   // we must put a control on the shell that are passed: 
00335   // some shells should not pass (line "0" or "2")
00336 
00337   // check atomic number
00338   G4double xsec = 0.0;
00339   if(Z > 93 || Z < 6 ) { return xsec; } //corrected by alf - Z<6 missing
00340   G4int idx = G4int(shellEnum);
00341   if(idx >= G4AtomicShells::GetNumberOfShells(Z)) { return xsec; }
00342 
00343   // 
00344   if(pdef == theElectron || pdef == thePositron) {
00345     xsec = ePIXEshellCS->CrossSection(Z,shellEnum,kineticEnergy,0.0,mat);
00346     return xsec;
00347   }
00348 
00349   G4double mass = pdef->GetPDGMass();
00350   G4double escaled = kineticEnergy;
00351   G4double q2 = 0.0;
00352 
00353   // scaling to protons
00354   if ((pdef->GetParticleName() != "proton" && pdef->GetParticleName() != "alpha" ) )
00355   {
00356     mass = proton_mass_c2;
00357     escaled = kineticEnergy*mass/(pdef->GetPDGMass());
00358 
00359     if(mat) {
00360       q2 = emcorr->EffectiveChargeSquareRatio(pdef,mat,kineticEnergy);
00361     } else {
00362       G4double q = pdef->GetPDGCharge()/eplus;
00363       q2 = q*q;
00364     }
00365   }
00366   
00367   if(PIXEshellCS) { xsec = PIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); }
00368   if(xsec < 1e-100) { 
00369     
00370     xsec = anaPIXEshellCS->CrossSection(Z,shellEnum,escaled,mass,mat); 
00371     
00372   }
00373 
00374   if (q2)  {xsec *= q2;}
00375 
00376   return xsec;
00377 }

void G4UAtomicDeexcitation::InitialiseForExtraAtom ( G4int  Z  )  [virtual]

Implements G4VAtomDeexcitation.

Definition at line 225 of file G4UAtomicDeexcitation.cc.

00226 {}

void G4UAtomicDeexcitation::InitialiseForNewRun (  )  [virtual]

Implements G4VAtomDeexcitation.

Definition at line 97 of file G4UAtomicDeexcitation.cc.

References G4cout, G4endl, G4VhShellCrossSection::GetName(), G4AtomicTransitionManager::Instance(), G4VAtomDeexcitation::IsFluoActive(), G4VAtomDeexcitation::IsPIXEActive(), G4VAtomDeexcitation::PIXECrossSectionModel(), G4VAtomDeexcitation::PIXEElectronCrossSectionModel(), G4VAtomDeexcitation::SetPIXECrossSectionModel(), and G4VAtomDeexcitation::SetPIXEElectronCrossSectionModel().

00098 {
00099   if(!IsFluoActive()) { return; }
00100   transitionManager = G4AtomicTransitionManager::Instance();
00101   if(IsPIXEActive()) {
00102     G4cout << G4endl;
00103     G4cout << "### === G4UAtomicDeexcitation::InitialiseForNewRun()" << G4endl;
00104     anaPIXEshellCS = new G4teoCrossSection("Analytical");
00105  
00106   }
00107   else  {return;}
00108   // initializing PIXE x-section name
00109   // 
00110   if (PIXECrossSectionModel() == "" ||
00111       PIXECrossSectionModel() == "Empirical" ||
00112       PIXECrossSectionModel() == "empirical") 
00113     {
00114       SetPIXECrossSectionModel("Empirical");
00115     }
00116   else if (PIXECrossSectionModel() == "ECPSSR_Analytical" ||
00117            PIXECrossSectionModel() == "Analytical" || 
00118            PIXECrossSectionModel() == "analytical") 
00119     {
00120       SetPIXECrossSectionModel("Analytical");
00121     }
00122   else if (PIXECrossSectionModel() == "ECPSSR_FormFactor" ||
00123            PIXECrossSectionModel() == "ECPSSR_Tabulated" ||
00124            PIXECrossSectionModel() == "Analytical_Tabulated") 
00125     {
00126       SetPIXECrossSectionModel("ECPSSR_FormFactor");
00127     }
00128   else 
00129     {
00130       G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
00131              << G4endl;
00132       G4cout << "    PIXE cross section name " << PIXECrossSectionModel()
00133              << " is unknown, Analytical cross section will be used" << G4endl; 
00134       SetPIXECrossSectionModel("Analytical");
00135     }
00136     
00137   // Check if old model should be deleted 
00138   if(PIXEshellCS) 
00139     {
00140       if(PIXECrossSectionModel() != PIXEshellCS->GetName()) 
00141         {
00142           delete PIXEshellCS;
00143           PIXEshellCS = 0;
00144         }
00145     }
00146 
00147   // Instantiate empirical model
00148   if(!PIXEshellCS) {
00149     if (PIXECrossSectionModel() == "Empirical")
00150       {
00151         PIXEshellCS = new G4empCrossSection("Empirical");
00152       }
00153 
00154     if (PIXECrossSectionModel() == "ECPSSR_FormFactor")
00155       {
00156         PIXEshellCS = new G4teoCrossSection("ECPSSR_FormFactor");
00157       }
00158   }
00159 
00160   // Electron cross section
00161   // initializing PIXE x-section name
00162   // 
00163   if (PIXEElectronCrossSectionModel() == "" ||
00164       PIXEElectronCrossSectionModel() == "Livermore")
00165     {
00166       SetPIXEElectronCrossSectionModel("Livermore");
00167     }
00168   else if (PIXEElectronCrossSectionModel() == "ProtonAnalytical" ||
00169            PIXEElectronCrossSectionModel() == "Analytical" || 
00170            PIXEElectronCrossSectionModel() == "analytical") 
00171     {
00172       SetPIXEElectronCrossSectionModel("ProtonAnalytical");
00173     }
00174   else if (PIXEElectronCrossSectionModel() == "ProtonEmpirical" ||
00175            PIXEElectronCrossSectionModel() == "Empirical" || 
00176            PIXEElectronCrossSectionModel() == "empirical") 
00177     {
00178       SetPIXEElectronCrossSectionModel("ProtonEmpirical");
00179     }
00180   else if (PIXEElectronCrossSectionModel() == "Penelope")
00181     SetPIXEElectronCrossSectionModel("Penelope");
00182   else 
00183     {
00184       G4cout << "### G4UAtomicDeexcitation::InitialiseForNewRun WARNING "
00185              << G4endl;
00186       G4cout << "    PIXE e- cross section name " << PIXEElectronCrossSectionModel()
00187              << " is unknown, PIXE is disabled" << G4endl; 
00188       SetPIXEElectronCrossSectionModel("Livermore");
00189     }
00190     
00191   // Check if old model should be deleted 
00192   if(ePIXEshellCS) 
00193     {
00194       if(PIXEElectronCrossSectionModel() != ePIXEshellCS->GetName()) 
00195         {
00196           delete ePIXEshellCS;
00197           ePIXEshellCS = 0;
00198         }
00199     }
00200 
00201   // Instantiate empirical model
00202   if(!ePIXEshellCS) 
00203     {
00204       if(PIXEElectronCrossSectionModel() == "Empirical")
00205         {
00206           ePIXEshellCS = new G4empCrossSection("Empirical");
00207         }
00208 
00209       else if(PIXEElectronCrossSectionModel() == "Analytical") 
00210         {
00211           ePIXEshellCS = new G4teoCrossSection("Analytical");
00212         }
00213 
00214       else if(PIXEElectronCrossSectionModel() == "Livermore")
00215         {
00216           ePIXEshellCS = new G4LivermoreIonisationCrossSection();
00217         }
00218       else if (PIXEElectronCrossSectionModel() == "Penelope")
00219         {
00220           ePIXEshellCS = new G4PenelopeIonisationCrossSection();
00221         }
00222     } 
00223 }

void G4UAtomicDeexcitation::SetCutForAugerElectrons ( G4double  cut  ) 

Definition at line 384 of file G4UAtomicDeexcitation.cc.

00385 {
00386   minElectronEnergy = cut;
00387 }

void G4UAtomicDeexcitation::SetCutForSecondaryPhotons ( G4double  cut  ) 

Definition at line 379 of file G4UAtomicDeexcitation.cc.

00380 {
00381   minGammaEnergy = cut;
00382 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:35 2013 for Geant4 by  doxygen 1.4.7