G4DNAChampionElasticModel Class Reference

#include <G4DNAChampionElasticModel.hh>

Inheritance diagram for G4DNAChampionElasticModel:

G4VEmModel

Public Member Functions

 G4DNAChampionElasticModel (const G4ParticleDefinition *p=0, const G4String &nam="DNAChampionElasticModel")
virtual ~G4DNAChampionElasticModel ()
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 41 of file G4DNAChampionElasticModel.hh.


Constructor & Destructor Documentation

G4DNAChampionElasticModel::G4DNAChampionElasticModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNAChampionElasticModel" 
)

Definition at line 40 of file G4DNAChampionElasticModel.cc.

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

00042 :G4VEmModel(nam),isInitialised(false)
00043 {
00044 //  nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
00045 
00046   killBelowEnergy = 7.4*eV; 
00047   lowEnergyLimit = 0 * eV; 
00048   highEnergyLimit = 1. * MeV;
00049   SetLowEnergyLimit(lowEnergyLimit);
00050   SetHighEnergyLimit(highEnergyLimit);
00051 
00052   verboseLevel= 0;
00053   // Verbosity scale:
00054   // 0 = nothing 
00055   // 1 = warning for energy non-conservation 
00056   // 2 = details of energy budget
00057   // 3 = calculation of cross sections, file openings, sampling of atoms
00058   // 4 = entering in methods
00059   
00060   if( verboseLevel>0 ) 
00061   { 
00062     G4cout << "Champion Elastic model is constructed " << G4endl
00063            << "Energy range: "
00064            << lowEnergyLimit / eV << " eV - "
00065            << highEnergyLimit / MeV << " MeV"
00066            << G4endl;
00067   }
00068   fParticleChangeForGamma = 0;
00069   fpMolWaterDensity = 0;
00070 }

G4DNAChampionElasticModel::~G4DNAChampionElasticModel (  )  [virtual]

Definition at line 74 of file G4DNAChampionElasticModel.cc.

00075 {  
00076   // For total cross section
00077   
00078   std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
00079   for (pos = tableData.begin(); pos != tableData.end(); ++pos)
00080   {
00081     G4DNACrossSectionDataSet* table = pos->second;
00082     delete table;
00083   }
00084 
00085    // For final state
00086    
00087    eVecm.clear();
00088 
00089 }


Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 203 of file G4DNAChampionElasticModel.cc.

References DBL_MAX, FatalException, G4cout, G4endl, G4Exception(), G4Material::GetIndex(), and G4ParticleDefinition::GetParticleName().

00208 {
00209   if (verboseLevel > 3)
00210     G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel" << G4endl;
00211 
00212  // Calculate total cross section for model
00213 
00214  G4double sigma=0;
00215 
00216  G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
00217 
00218  if(waterDensity!= 0.0)
00219 //  if (material == nistwater || material->GetBaseMaterial() == nistwater)
00220   {
00221   const G4String& particleName = p->GetParticleName();
00222 
00223   if (ekin < highEnergyLimit)
00224   {
00225       //SI : XS must not be zero otherwise sampling of secondaries method ignored
00226       if (ekin < killBelowEnergy) return DBL_MAX;
00227       //      
00228       
00229         std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
00230         pos = tableData.find(particleName);
00231         
00232         if (pos != tableData.end())
00233         {
00234           G4DNACrossSectionDataSet* table = pos->second;
00235           if (table != 0)
00236           {
00237             sigma = table->FindValue(ekin);
00238           }
00239         }
00240         else
00241         {
00242            G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume","em0002",
00243                       FatalException,"Model not applicable to particle type.");
00244         }
00245   }
00246 
00247   if (verboseLevel > 2)
00248   {
00249     G4cout << "__________________________________" << G4endl;
00250     G4cout << "°°° G4DNAChampionElasticModel - XS INFO START" << G4endl;
00251     G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleName << G4endl;
00252     G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
00253     G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
00254     //      G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
00255     G4cout << "°°° G4DNAChampionElasticModel - XS INFO END" << G4endl;
00256   } 
00257 
00258  } 
00259          
00260  return sigma*waterDensity;
00261 // return sigma*material->GetAtomicNumDensityVector()[1];
00262 }

G4double G4DNAChampionElasticModel::GetKillBelowThreshold (  )  [inline]

Definition at line 66 of file G4DNAChampionElasticModel.hh.

00066 { return killBelowEnergy; }              

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

Implements G4VEmModel.

Definition at line 93 of file G4DNAChampionElasticModel.cc.

References G4Electron::ElectronDefinition(), FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4DNAMolecularMaterial::Instance(), G4DNACrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

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

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

Implements G4VEmModel.

Definition at line 266 of file G4DNAChampionElasticModel.cc.

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

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

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

Definition at line 133 of file G4DNAChampionElasticModel.hh.

References G4Exception(), and JustWarning.

Referenced by G4EmDNAPhysicsChemistry::ConstructProcess().

00134 { 
00135 
00136 // SI - commented on 19/06/2013
00137 /*
00138   killBelowEnergy = threshold; 
00139     
00140   if (threshold < 1*eV)
00141      G4Exception ("*** WARNING : the G4DNAChampionElasticModel class is not validated below 1 eV !","",JustWarning,"") ;
00142 
00143   if (threshold < 0.025*eV) threshold = 0.025*eV;
00144 */
00145 
00146   G4Exception ("*** WARNING : G4DNAChampionElasticModel::SetKillBelowThreshold INACTIVE for now","",JustWarning,"") ;
00147             
00148 }                


Field Documentation

G4ParticleChangeForGamma* G4DNAChampionElasticModel::fParticleChangeForGamma [protected]

Definition at line 70 of file G4DNAChampionElasticModel.hh.

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


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