G4MuElecElasticModel Class Reference

#include <G4MuElecElasticModel.hh>

Inheritance diagram for G4MuElecElasticModel:

G4VEmModel

Public Member Functions

 G4MuElecElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="MuElecElasticModel")
virtual ~G4MuElecElasticModel ()
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetKillBelowThreshold (G4double threshold)
G4double GetKillBelowThreshold ()

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma

Detailed Description

Definition at line 51 of file G4MuElecElasticModel.hh.


Constructor & Destructor Documentation

G4MuElecElasticModel::G4MuElecElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MuElecElasticModel" 
)

Definition at line 48 of file G4MuElecElasticModel.cc.

References G4NistManager::FindOrBuildMaterial(), fParticleChangeForGamma, G4cout, G4endl, G4NistManager::Instance(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

00050 :G4VEmModel(nam),isInitialised(false)
00051 {
00052   nistSi = G4NistManager::Instance()->FindOrBuildMaterial("G4_Si");
00053 
00054   killBelowEnergy = 16.7 * eV; // Minimum e- energy for energy loss by excitation
00055   lowEnergyLimit = 0 * eV; 
00056   lowEnergyLimitOfModel = 5 * eV; // The model lower energy is 5 eV
00057   highEnergyLimit = 100. * MeV;
00058   SetLowEnergyLimit(lowEnergyLimit);
00059   SetHighEnergyLimit(highEnergyLimit);
00060 
00061   verboseLevel= 0;
00062   // Verbosity scale:
00063   // 0 = nothing 
00064   // 1 = warning for energy non-conservation 
00065   // 2 = details of energy budget
00066   // 3 = calculation of cross sections, file openings, sampling of atoms
00067   // 4 = entering in methods
00068   
00069   if( verboseLevel>0 ) 
00070   { 
00071     G4cout << "MuElec Elastic model is constructed " << G4endl
00072            << "Energy range: "
00073            << lowEnergyLimit / eV << " eV - "
00074            << highEnergyLimit / keV << " keV"
00075            << G4endl;
00076   }
00077   fParticleChangeForGamma = 0;
00078 }

G4MuElecElasticModel::~G4MuElecElasticModel (  )  [virtual]

Definition at line 82 of file G4MuElecElasticModel.cc.

00083 {  
00084   // For total cross section
00085   
00086   std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
00087   for (pos = tableData.begin(); pos != tableData.end(); ++pos)
00088   {
00089     G4MuElecCrossSectionDataSet* table = pos->second;
00090     delete table;
00091   }
00092 
00093    // For final state
00094    
00095    eVecm.clear();
00096 
00097 }


Member Function Documentation

G4double G4MuElecElasticModel::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 207 of file G4MuElecElasticModel.cc.

References FatalException, G4cout, G4endl, G4Exception(), G4Material::GetBaseMaterial(), G4ParticleDefinition::GetParticleName(), and G4Material::GetTotNbOfAtomsPerVolume().

00212 {
00213   if (verboseLevel > 3)
00214     G4cout << "Calling CrossSectionPerVolume() of G4MuElecElasticModel" << G4endl;
00215 
00216  // Calculate total cross section for model
00217 
00218  G4double sigma=0;
00219  
00220  G4double density = material->GetTotNbOfAtomsPerVolume();
00221 
00222  if (material == nistSi || material->GetBaseMaterial() == nistSi)
00223  {
00224   const G4String& particleName = p->GetParticleName();
00225 
00226   if (ekin < highEnergyLimit)
00227   {
00228       //SI : XS must not be zero otherwise sampling of secondaries method ignored
00229       if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel;
00230       //      
00231       
00232         std::map< G4String,G4MuElecCrossSectionDataSet*,std::less<G4String> >::iterator pos;
00233         pos = tableData.find(particleName);
00234         
00235         if (pos != tableData.end())
00236         {
00237           G4MuElecCrossSectionDataSet* table = pos->second;
00238           if (table != 0)
00239           {
00240             sigma = table->FindValue(ekin);
00241           }
00242         }
00243         else
00244         {
00245             G4Exception("G4MuElecElasticModel::ComputeCrossSectionPerVolume","em0002",FatalException,"Model not applicable to particle type.");
00246         }
00247   }
00248 
00249   if (verboseLevel > 3)
00250   {
00251     G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
00252     G4cout << " - Cross section per Si atom (cm^2)=" << sigma/cm/cm << G4endl;
00253     G4cout << " - Cross section per Si atom (cm^-1)=" << sigma*density/(1./cm) << G4endl;
00254   } 
00255 
00256  } 
00257          
00258  return sigma*density;             
00259 }

G4double G4MuElecElasticModel::GetKillBelowThreshold (  )  [inline]

Definition at line 76 of file G4MuElecElasticModel.hh.

00076 { return killBelowEnergy; }     

void G4MuElecElasticModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
) [virtual]

Implements G4VEmModel.

Definition at line 101 of file G4MuElecElasticModel.cc.

References G4Electron::ElectronDefinition(), FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4MuElecCrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

00103 {
00104 
00105   if (verboseLevel > 3)
00106     G4cout << "Calling G4MuElecElasticModel::Initialise()" << G4endl;
00107 
00108   // Energy limits
00109   
00110   if (LowEnergyLimit() < lowEnergyLimit)
00111   {
00112     G4cout << "G4MuElecElasticModel: low energy limit increased from " << 
00113         LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
00114     SetLowEnergyLimit(lowEnergyLimit);
00115     }
00116 
00117   if (HighEnergyLimit() > highEnergyLimit)
00118   {
00119     G4cout << "G4MuElecElasticModel: high energy limit decreased from " << 
00120         HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
00121     SetHighEnergyLimit(highEnergyLimit);
00122   }
00123 
00124   // Reading of data files 
00125   
00126   G4double scaleFactor = 1e-18 * cm * cm;
00127 
00128   G4String fileElectron("muelec/sigma_elastic_e_Si");
00129 
00130   G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
00131   G4String electron;
00132 
00133     // For total cross section
00134     
00135     electron = electronDef->GetParticleName();
00136 
00137     tableFile[electron] = fileElectron;
00138 
00139     G4MuElecCrossSectionDataSet* tableE = new G4MuElecCrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
00140     tableE->LoadData(fileElectron);
00141     tableData[electron] = tableE;
00142     
00143     // For final state
00144     
00145     char *path = getenv("G4LEDATA");
00146  
00147     if (!path)
00148     {
00149       G4Exception("G4MuElecElasticModel::Initialise","em0006",FatalException,"G4LEDATA environment variable not set.");
00150       return;
00151     }
00152 
00153     std::ostringstream eFullFileName;
00154     eFullFileName << path << "/muelec/sigmadiff_elastic_e_Si.dat";
00155     std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
00156      
00157     if (!eDiffCrossSection) 
00158         G4Exception("G4MuElecElasticModel::Initialise","em0003",FatalException,"Missing data file: /muelec/sigmadiff_elastic_e_Si.dat");
00159       
00160     eTdummyVec.push_back(0.);
00161 
00162     while(!eDiffCrossSection.eof())
00163     {
00164         double tDummy;
00165         double eDummy;
00166         eDiffCrossSection>>tDummy>>eDummy;
00167 
00168         // SI : mandatory eVecm initialization
00169         if (tDummy != eTdummyVec.back()) 
00170         { 
00171           eTdummyVec.push_back(tDummy); 
00172           eVecm[tDummy].push_back(0.);
00173         }
00174           
00175         eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
00176 
00177         // SI : only if not end of file reached !
00178         if (!eDiffCrossSection.eof()) eDiffCrossSectionData[tDummy][eDummy]*=scaleFactor;
00179           
00180         if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
00181           
00182     }
00183 
00184     // End final state
00185   
00186   if (verboseLevel > 2) 
00187     G4cout << "Loaded cross section files for MuElec Elastic model" << G4endl;
00188 
00189   if( verboseLevel>0 ) 
00190   { 
00191     G4cout << "MuElec Elastic model is initialized " << G4endl
00192            << "Energy range: "
00193            << LowEnergyLimit() / eV << " eV - "
00194            << HighEnergyLimit() / keV << " keV"
00195            << G4endl;
00196   }
00197 
00198   if (isInitialised) { return; }
00199   fParticleChangeForGamma = GetParticleChangeForGamma();
00200   isInitialised = true;
00201 
00202   // InitialiseElementSelectors(particle,cuts);
00203 }

void G4MuElecElasticModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  maxEnergy 
) [virtual]

Implements G4VEmModel.

Definition at line 263 of file G4MuElecElasticModel.cc.

References fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4INCL::Math::pi, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

00268 {
00269 
00270   if (verboseLevel > 3)
00271     G4cout << "Calling SampleSecondaries() of G4MuElecElasticModel" << G4endl;
00272 
00273   G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
00274   
00275   if (electronEnergy0 < killBelowEnergy)
00276   {
00277     fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
00278     fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
00279     return ;
00280   }
00281 
00282   if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
00283   {  
00284     G4double cosTheta = RandomizeCosTheta(electronEnergy0);
00285   
00286     G4double phi = 2. * pi * G4UniformRand();
00287 
00288     G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
00289     G4ThreeVector xVers = zVers.orthogonal();
00290     G4ThreeVector yVers = zVers.cross(xVers);
00291 
00292     G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
00293     G4double yDir = xDir;
00294     xDir *= std::cos(phi);
00295     yDir *= std::sin(phi);
00296 
00297     G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
00298 
00299     fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit()) ;
00300 
00301     fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
00302   }
00303 
00304 }

void G4MuElecElasticModel::SetKillBelowThreshold ( G4double  threshold  )  [inline]

Definition at line 138 of file G4MuElecElasticModel.hh.

References G4Exception(), and JustWarning.

00139 { 
00140     killBelowEnergy = threshold; 
00141     
00142     if (threshold < 5*CLHEP::eV)
00143     {
00144        G4Exception ("*** WARNING : the G4MuElecElasticModel class is not validated below 5 eV !","",JustWarning,"") ;
00145        threshold = 0.025*CLHEP::eV;
00146     }
00147              
00148 }                


Field Documentation

G4ParticleChangeForGamma* G4MuElecElasticModel::fParticleChangeForGamma [protected]

Definition at line 80 of file G4MuElecElasticModel.hh.

Referenced by G4MuElecElasticModel(), Initialise(), and SampleSecondaries().


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