#include <G4DNASancheExcitationModel.hh>
Inheritance diagram for G4DNASancheExcitationModel:
Public Member Functions | |
G4DNASancheExcitationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNASancheExcitationModel") | |
virtual | ~G4DNASancheExcitationModel () |
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) |
G4double | PartialCrossSection (G4double energy, G4int level) |
void | ExtendLowEnergyLimit (G4double) |
void | SetVerboseLevel (int verbose) |
Protected Attributes | |
G4ParticleChangeForGamma * | fParticleChangeForGamma |
Definition at line 45 of file G4DNASancheExcitationModel.hh.
G4DNASancheExcitationModel::G4DNASancheExcitationModel | ( | const G4ParticleDefinition * | p = 0 , |
|
const G4String & | nam = "DNASancheExcitationModel" | |||
) |
Definition at line 41 of file G4DNASancheExcitationModel.cc.
References fParticleChangeForGamma, G4cout, G4endl, G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().
00043 :G4VEmModel(nam),isInitialised(false) 00044 { 00045 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER"); 00046 fpWaterDensity = 0; 00047 00048 lowEnergyLimit = 2 * eV; 00049 highEnergyLimit = 100 * eV; 00050 SetLowEnergyLimit(lowEnergyLimit); 00051 SetHighEnergyLimit(highEnergyLimit); 00052 nLevels = 9; 00053 00054 verboseLevel= 0; 00055 // Verbosity scale: 00056 // 0 = nothing 00057 // 1 = warning for energy non-conservation 00058 // 2 = details of energy budget 00059 // 3 = calculation of cross sections, file openings, sampling of atoms 00060 // 4 = entering in methods 00061 00062 if (verboseLevel > 0) 00063 { 00064 G4cout << "Sanche Excitation model is constructed " << G4endl 00065 << "Energy range: " 00066 << lowEnergyLimit / eV << " eV - " 00067 << highEnergyLimit / eV << " eV" 00068 << G4endl; 00069 } 00070 fParticleChangeForGamma = 0; 00071 fpWaterDensity = 0; 00072 }
G4DNASancheExcitationModel::~G4DNASancheExcitationModel | ( | ) | [virtual] |
G4double G4DNASancheExcitationModel::CrossSectionPerVolume | ( | const G4Material * | material, | |
const G4ParticleDefinition * | p, | |||
G4double | ekin, | |||
G4double | emin, | |||
G4double | emax | |||
) | [virtual] |
Reimplemented from G4VEmModel.
Definition at line 146 of file G4DNASancheExcitationModel.cc.
References G4Electron::ElectronDefinition(), G4cout, G4endl, G4Material::GetIndex(), and G4ParticleDefinition::GetParticleName().
00151 { 00152 if (verboseLevel > 3) 00153 G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel" << G4endl; 00154 00155 // Calculate total cross section for model 00156 00157 G4double sigma=0; 00158 00159 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()]; 00160 00161 if(waterDensity!= 0.0) 00162 // if (material == nistwater || material->GetBaseMaterial() == nistwater) 00163 { 00164 00165 if (particleDefinition == G4Electron::ElectronDefinition()) 00166 { 00167 if (ekin >= lowEnergyLimit && ekin < highEnergyLimit) 00168 { 00169 sigma = Sum(ekin); 00170 } 00171 } 00172 00173 if (verboseLevel > 2) 00174 { 00175 G4cout << "__________________________________" << G4endl; 00176 G4cout << "°°° G4DNASancheExcitationModel - XS INFO START" << G4endl; 00177 G4cout << "°°° Kinetic energy(eV)=" << ekin/eV << " particle : " << particleDefinition->GetParticleName() << G4endl; 00178 G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; 00179 G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl; 00180 // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl; 00181 G4cout << "°°° G4DNASancheExcitationModel - XS INFO END" << G4endl; 00182 } 00183 00184 } // if water 00185 00186 00187 // return sigma*2*material->GetAtomicNumDensityVector()[1]; 00188 return sigma*2*waterDensity; 00189 // see papers for factor 2 description 00190 00191 }
void G4DNASancheExcitationModel::ExtendLowEnergyLimit | ( | G4double | ) | [inline] |
Definition at line 113 of file G4DNASancheExcitationModel.hh.
References G4Exception(), and JustWarning.
Referenced by G4EmDNAPhysicsChemistry::ConstructProcess().
00114 { 00115 lowEnergyLimit = threshold; 00116 if (lowEnergyLimit < 2*CLHEP::eV) 00117 G4Exception ("*** WARNING : the G4DNASancheExcitationModel class is not validated below 2 eV !","",JustWarning,"") ; 00118 }
void G4DNASancheExcitationModel::Initialise | ( | const G4ParticleDefinition * | , | |
const G4DataVector & | ||||
) | [virtual] |
Implements G4VEmModel.
Definition at line 81 of file G4DNASancheExcitationModel.cc.
References FatalException, fParticleChangeForGamma, G4cout, G4endl, G4Exception(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4VEmModel::HighEnergyLimit(), G4DNAMolecularMaterial::Instance(), G4VEmModel::LowEnergyLimit(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().
00083 { 00084 00085 00086 if (verboseLevel > 3) 00087 G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl; 00088 00089 // Energy limits 00090 00091 if (LowEnergyLimit() < lowEnergyLimit) 00092 { 00093 G4cout << "G4DNASancheExcitationModel: low energy limit increased from " << 00094 LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl; 00095 SetLowEnergyLimit(lowEnergyLimit); 00096 } 00097 00098 if (HighEnergyLimit() > highEnergyLimit) 00099 { 00100 G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " << 00101 HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl; 00102 SetHighEnergyLimit(highEnergyLimit); 00103 } 00104 00105 // 00106 00107 if (verboseLevel > 0) 00108 { 00109 G4cout << "Sanche Excitation model is initialized " << G4endl 00110 << "Energy range: " 00111 << LowEnergyLimit() / eV << " eV - " 00112 << HighEnergyLimit() / eV << " eV" 00113 << G4endl; 00114 } 00115 00116 // Initialize water density pointer 00117 fpWaterDensity = G4DNAMolecularMaterial::Instance()->GetNumMolPerVolTableFor(G4Material::GetMaterial("G4_WATER")); 00118 00119 if (isInitialised) { return; } 00120 fParticleChangeForGamma = GetParticleChangeForGamma(); 00121 isInitialised = true; 00122 00123 char *path = getenv("G4LEDATA"); 00124 std::ostringstream eFullFileName; 00125 eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat"; 00126 std::ifstream input(eFullFileName.str().c_str()); 00127 00128 if (!input) 00129 { 00130 G4Exception("G4DNASancheExcitationModel::Initialise","em0003", 00131 FatalException,"Missing data file:/dna/sigma_excitationvib_e_sanche.dat"); 00132 } 00133 00134 while(!input.eof()) 00135 { 00136 double t; 00137 input>>t; 00138 tdummyVec.push_back(t); 00139 input>>map1[t][0]>>map1[t][1]>>map1[t][2]>>map1[t][3]>>map1[t][4]>>map1[t][5]>>map1[t][6]>>map1[t][7]>>map1[t][8]; 00140 //G4cout<<t<<" "<<map1[t][0]<<map1[t][1]<<map1[t][2]<<map1[t][3]<<map1[t][4]<<map1[t][5]<<map1[t][6]<<map1[t][7]<<map1[t][8]<<G4endl; 00141 } 00142 }
Definition at line 240 of file G4DNASancheExcitationModel.cc.
00241 { 00242 std::vector<double>::iterator t2 = std::upper_bound(tdummyVec.begin(),tdummyVec.end(), t/eV); 00243 std::vector<double>::iterator t1 = t2-1; 00244 00245 double sigma = LinInterpolate((*t1), (*t2), t/eV, map1[*t1][level], map1[*t2][level]); 00246 sigma*=1e-16*cm*cm; 00247 if(sigma==0.)sigma=1e-30; 00248 return (sigma); 00249 }
void G4DNASancheExcitationModel::SampleSecondaries | ( | std::vector< G4DynamicParticle * > * | , | |
const G4MaterialCutsCouple * | , | |||
const G4DynamicParticle * | , | |||
G4double | tmin, | |||
G4double | maxEnergy | |||
) | [virtual] |
Implements G4VEmModel.
Definition at line 195 of file G4DNASancheExcitationModel.cc.
References fParticleChangeForGamma, G4cout, G4endl, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().
00200 { 00201 00202 if (verboseLevel > 3) 00203 G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel" << G4endl; 00204 00205 G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy(); 00206 G4int level = RandomSelect(electronEnergy0); 00207 G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8 00208 G4double newEnergy = electronEnergy0 - excitationEnergy; 00209 00210 /* 00211 if (electronEnergy0 < highEnergyLimit) 00212 { 00213 if (newEnergy >= lowEnergyLimit) 00214 { 00215 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); 00216 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); 00217 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); 00218 } 00219 00220 else 00221 { 00222 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); 00223 fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0); 00224 } 00225 } 00226 */ 00227 00228 if (electronEnergy0 < highEnergyLimit && newEnergy>0.) 00229 { 00230 fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); 00231 fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); 00232 fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); 00233 } 00234 00235 // 00236 }
void G4DNASancheExcitationModel::SetVerboseLevel | ( | int | verbose | ) | [inline] |
Definition at line 79 of file G4DNASancheExcitationModel.hh.
Referenced by G4DNASancheExcitationModel(), Initialise(), and SampleSecondaries().