// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // $Id: G4DNASancheExcitationModel.cc,v 1.4 2010/11/11 22:32:22 sincerti Exp $ // GEANT4 tag $Name: geant4-09-04-ref-00 $ // // Created by Z. Francis #include "G4DNASancheExcitationModel.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... using namespace std; //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... G4DNASancheExcitationModel::G4DNASancheExcitationModel(const G4ParticleDefinition*, const G4String& nam) :G4VEmModel(nam),isInitialised(false) { lowEnergyLimit = 2 * eV; highEnergyLimit = 100 * eV; SetLowEnergyLimit(lowEnergyLimit); SetHighEnergyLimit(highEnergyLimit); nLevels = 9; verboseLevel= 0; // Verbosity scale: // 0 = nothing // 1 = warning for energy non-conservation // 2 = details of energy budget // 3 = calculation of cross sections, file openings, sampling of atoms // 4 = entering in methods if (verboseLevel > 0) { G4cout << "Sanche Excitation model is constructed " << G4endl << "Energy range: " << lowEnergyLimit / eV << " eV - " << highEnergyLimit / eV << " eV" << G4endl; } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... G4DNASancheExcitationModel::~G4DNASancheExcitationModel() {} //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... void G4DNASancheExcitationModel::Initialise(const G4ParticleDefinition* /*particle*/, const G4DataVector& /*cuts*/) { if (verboseLevel > 3) G4cout << "Calling G4DNASancheExcitationModel::Initialise()" << G4endl; // Energy limits if (LowEnergyLimit() < lowEnergyLimit) { G4cout << "G4DNASancheExcitationModel: low energy limit increased from " << LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl; SetLowEnergyLimit(lowEnergyLimit); } if (HighEnergyLimit() > highEnergyLimit) { G4cout << "G4DNASancheExcitationModel: high energy limit decreased from " << HighEnergyLimit()/eV << " eV to " << highEnergyLimit/eV << " eV" << G4endl; SetHighEnergyLimit(highEnergyLimit); } // if (verboseLevel > 0) G4cout << "Sanche Excitation model is initialized " << G4endl << "Energy range: " << LowEnergyLimit() / eV << " eV - " << HighEnergyLimit() / eV << " eV" << G4endl; if(!isInitialised) { isInitialised = true; if(pParticleChange) fParticleChangeForGamma = reinterpret_cast(pParticleChange); else fParticleChangeForGamma = new G4ParticleChangeForGamma(); } // InitialiseElementSelectors(particle,cuts); char *path = getenv("G4LEDATA"); std::ostringstream eFullFileName; eFullFileName << path << "/dna/sigma_excitationvib_e_sanche.dat"; std::ifstream input(eFullFileName.str().c_str()); if (!input) { G4Exception("G4DNASancheExcitationModel:::ERROR OPENING XS DATA FILE"); } while(!input.eof()) { double t; input>>t; tdummyVec.push_back(t); 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]; //G4cout< 3) G4cout << "Calling CrossSectionPerVolume() of G4DNASancheExcitationModel" << G4endl; // Calculate total cross section for model G4double sigma=0; if (material->GetName() == "G4_WATER") { if (particleDefinition == G4Electron::ElectronDefinition()) { if (ekin >= lowEnergyLimit && ekin < highEnergyLimit) { sigma = Sum(ekin); } } if (verboseLevel > 3) { G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl; G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl; G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl; } } // if water return sigma*2*material->GetAtomicNumDensityVector()[1]; // see papers for factor 2 description } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... void G4DNASancheExcitationModel::SampleSecondaries(std::vector*, const G4MaterialCutsCouple*, const G4DynamicParticle* aDynamicElectron, G4double, G4double) { if (verboseLevel > 3) G4cout << "Calling SampleSecondaries() of G4DNASancheExcitationModel" << G4endl; G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy(); G4int level = RandomSelect(electronEnergy0); G4double excitationEnergy = VibrationEnergy(level); // levels go from 0 to 8 G4double newEnergy = electronEnergy0 - excitationEnergy; /* if (electronEnergy0 < highEnergyLimit) { if (newEnergy >= lowEnergyLimit) { fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); } else { fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill); fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0); } } */ if (electronEnergy0 < highEnergyLimit && newEnergy>0.) { fParticleChangeForGamma->ProposeMomentumDirection(aDynamicElectron->GetMomentumDirection()); fParticleChangeForGamma->SetProposedKineticEnergy(newEnergy); fParticleChangeForGamma->ProposeLocalEnergyDeposit(excitationEnergy); } // } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4double G4DNASancheExcitationModel::PartialCrossSection(G4double t, G4int level) { std::vector::iterator t2 = std::upper_bound(tdummyVec.begin(),tdummyVec.end(), t/eV); std::vector::iterator t1 = t2-1; double sigma = LinInterpolate((*t1), (*t2), t/eV, map1[*t1][level], map1[*t2][level]); sigma*=1e-16*cm*cm; if(sigma==0.)sigma=1e-30; return (sigma); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4double G4DNASancheExcitationModel::VibrationEnergy(G4int level) { G4double energies[9] = {0.01, 0.024, 0.061, 0.092, 0.204, 0.417, 0.460, 0.500, 0.835}; return(energies[level]*eV); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... G4int G4DNASancheExcitationModel::RandomSelect(G4double k) { // Level Selection Counting can be done here ! G4int i = nLevels; G4double value = 0.; std::deque values; while (i > 0) { i--; G4double partial = PartialCrossSection(k,i); values.push_front(partial); value += partial; } value *= G4UniformRand(); i = nLevels; while (i > 0) { i--; if (values[i] > value) { //outcount<> "<