G4teoCrossSection Class Reference

#include <G4teoCrossSection.hh>

Inheritance diagram for G4teoCrossSection:

G4VhShellCrossSection

Public Member Functions

 G4teoCrossSection (const G4String &name)
virtual ~G4teoCrossSection ()
std::vector< G4doubleGetCrossSection (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
G4double CrossSection (G4int Z, G4AtomicShellEnumerator shell, G4double incidentEnergy, G4double mass, const G4Material *mat)
std::vector< G4doubleProbabilities (G4int Z, G4double incidentEnergy, G4double mass, G4double deltaEnergy=0, const G4Material *mat=0)
void SetTotalCS (G4double)

Detailed Description

Definition at line 56 of file G4teoCrossSection.hh.


Constructor & Destructor Documentation

G4teoCrossSection::G4teoCrossSection ( const G4String name  ) 

Definition at line 50 of file G4teoCrossSection.cc.

References G4cout, and G4endl.

00051   :G4VhShellCrossSection(nam),totalCS(0.0),ecpssrShellK(0),ecpssrShellLi(0), ecpssrShellMi(0)
00052 { 
00053 
00054   if (nam == "Analytical") 
00055     {
00056       ecpssrShellK  = new G4ecpssrBaseKxsModel();  
00057       ecpssrShellLi = new G4ecpssrBaseLixsModel();      
00058     }
00059   else if (nam == "ECPSSR_FormFactor")
00060     {
00061       ecpssrShellK  = new G4ecpssrFormFactorKxsModel();  
00062       ecpssrShellLi = new G4ecpssrFormFactorLixsModel();
00063       ecpssrShellMi = new G4ecpssrFormFactorMixsModel();
00064     }
00065   else { G4cout << "ERROR" << G4endl;}
00066 
00067 
00068 }

G4teoCrossSection::~G4teoCrossSection (  )  [virtual]

Definition at line 70 of file G4teoCrossSection.cc.

00071 { 
00072   delete ecpssrShellK;
00073   delete ecpssrShellLi;
00074   if (ecpssrShellMi) {delete ecpssrShellMi;}
00075 }


Member Function Documentation

G4double G4teoCrossSection::CrossSection ( G4int  Z,
G4AtomicShellEnumerator  shell,
G4double  incidentEnergy,
G4double  mass,
const G4Material mat 
) [virtual]

Implements G4VhShellCrossSection.

Definition at line 105 of file G4teoCrossSection.cc.

References G4VecpssrKModel::CalculateCrossSection(), G4VecpssrLiModel::CalculateL1CrossSection(), G4VecpssrLiModel::CalculateL2CrossSection(), G4VecpssrLiModel::CalculateL3CrossSection(), fKShell, fL1Shell, fL2Shell, fL3Shell, fM1Shell, fM2Shell, fM3Shell, fM4Shell, and fM5Shell.

00109 {
00110   G4double res = 0.0;
00111   if(shell > 3 && !ecpssrShellMi) {
00112     return res; 
00113   } 
00114 
00115   else if(shell > 8) { 
00116     return res;
00117   } 
00118   
00119   else if(fKShell  == shell) 
00120     { 
00121       res = ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy);
00122     } 
00123   
00124   else if(fL1Shell == shell) 
00125     { 
00126       res = ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy);
00127     } 
00128   
00129   else if(fL2Shell == shell) 
00130     { 
00131       res = ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy);
00132     } 
00133   
00134   else if(fL3Shell == shell) 
00135     { 
00136       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00137     }
00138 
00139   else if(fM1Shell == shell) 
00140     { 
00141       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00142     }
00143 
00144   else if(fM2Shell == shell) 
00145     { 
00146       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00147     }
00148 
00149   else if(fM3Shell == shell) 
00150     { 
00151       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00152     }
00153 
00154   else if(fM4Shell == shell) 
00155     { 
00156       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00157     }
00158 
00159   else if(fM5Shell == shell) 
00160     { 
00161       res = ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy);
00162     }
00163   return res;
00164 }

std::vector< G4double > G4teoCrossSection::GetCrossSection ( G4int  Z,
G4double  incidentEnergy,
G4double  mass,
G4double  deltaEnergy = 0,
const G4Material mat = 0 
) [virtual]

Implements G4VhShellCrossSection.

Definition at line 77 of file G4teoCrossSection.cc.

References G4VecpssrKModel::CalculateCrossSection(), G4VecpssrLiModel::CalculateL1CrossSection(), G4VecpssrLiModel::CalculateL2CrossSection(), and G4VecpssrLiModel::CalculateL3CrossSection().

Referenced by Probabilities().

00082 {
00083   std::vector<G4double> crossSections;
00084 
00085   crossSections.push_back( ecpssrShellK->CalculateCrossSection(Z, mass, incidentEnergy) );
00086   
00087   crossSections.push_back( ecpssrShellLi->CalculateL1CrossSection(Z, mass, incidentEnergy) );
00088   crossSections.push_back( ecpssrShellLi->CalculateL2CrossSection(Z, mass, incidentEnergy) );
00089   crossSections.push_back( ecpssrShellLi->CalculateL3CrossSection(Z, mass, incidentEnergy) );
00090 
00091   if (ecpssrShellMi) {
00092 
00093     crossSections.push_back( ecpssrShellMi->CalculateM1CrossSection(Z, mass, incidentEnergy) );
00094     crossSections.push_back( ecpssrShellMi->CalculateM2CrossSection(Z, mass, incidentEnergy) );
00095     crossSections.push_back( ecpssrShellMi->CalculateM3CrossSection(Z, mass, incidentEnergy) );
00096     crossSections.push_back( ecpssrShellMi->CalculateM4CrossSection(Z, mass, incidentEnergy) );
00097     crossSections.push_back( ecpssrShellMi->CalculateM5CrossSection(Z, mass, incidentEnergy) );
00098 
00099   }
00100 
00101 
00102   return crossSections;
00103 }

std::vector< G4double > G4teoCrossSection::Probabilities ( G4int  Z,
G4double  incidentEnergy,
G4double  mass,
G4double  deltaEnergy = 0,
const G4Material mat = 0 
) [virtual]

Implements G4VhShellCrossSection.

Definition at line 166 of file G4teoCrossSection.cc.

References GetCrossSection().

00171 {
00172   std::vector<G4double> crossSections = 
00173     GetCrossSection(Z, incidentEnergy, mass, deltaEnergy);
00174 
00175   for (size_t i=0; i<crossSections.size(); i++ ) {
00176     
00177     if (totalCS) {
00178       crossSections[i] = crossSections[i]/totalCS;
00179     }
00180     
00181   }
00182   return crossSections;
00183 }

void G4teoCrossSection::SetTotalCS ( G4double   )  [virtual]

Reimplemented from G4VhShellCrossSection.

Definition at line 185 of file G4teoCrossSection.cc.

00185                                               {
00186 
00187   totalCS = val;
00188   //  G4cout << "totalXS set to: " << val / barn << " barns" << G4endl;
00189 }


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