| 1 | //
|
|---|
| 2 | // ********************************************************************
|
|---|
| 3 | // * License and Disclaimer *
|
|---|
| 4 | // * *
|
|---|
| 5 | // * The Geant4 software is copyright of the Copyright Holders of *
|
|---|
| 6 | // * the Geant4 Collaboration. It is provided under the terms and *
|
|---|
| 7 | // * conditions of the Geant4 Software License, included in the file *
|
|---|
| 8 | // * LICENSE and available at http://cern.ch/geant4/license . These *
|
|---|
| 9 | // * include a list of copyright holders. *
|
|---|
| 10 | // * *
|
|---|
| 11 | // * Neither the authors of this software system, nor their employing *
|
|---|
| 12 | // * institutes,nor the agencies providing financial support for this *
|
|---|
| 13 | // * work make any representation or warranty, express or implied, *
|
|---|
| 14 | // * regarding this software system or assume any liability for its *
|
|---|
| 15 | // * use. Please see the license in the file LICENSE and URL above *
|
|---|
| 16 | // * for the full disclaimer and the limitation of liability. *
|
|---|
| 17 | // * *
|
|---|
| 18 | // * This code implementation is the result of the scientific and *
|
|---|
| 19 | // * technical work of the GEANT4 collaboration. *
|
|---|
| 20 | // * By using, copying, modifying or distributing the software (or *
|
|---|
| 21 | // * any work based on the software) you agree to acknowledge its *
|
|---|
| 22 | // * use in resulting scientific publications, and indicate your *
|
|---|
| 23 | // * acceptance of all terms of the Geant4 Software license. *
|
|---|
| 24 | // ********************************************************************
|
|---|
| 25 | //
|
|---|
| 26 | // $Id: G4DNAChampionElasticModel.cc,v 1.10 2009/11/03 15:04:25 sincerti Exp $
|
|---|
| 27 | // GEANT4 tag $Name: geant4-09-03-cand-01 $
|
|---|
| 28 | //
|
|---|
| 29 |
|
|---|
| 30 | #include "G4DNAChampionElasticModel.hh"
|
|---|
| 31 |
|
|---|
| 32 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 33 |
|
|---|
| 34 | using namespace std;
|
|---|
| 35 |
|
|---|
| 36 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 37 |
|
|---|
| 38 | G4DNAChampionElasticModel::G4DNAChampionElasticModel(const G4ParticleDefinition*,
|
|---|
| 39 | const G4String& nam)
|
|---|
| 40 | :G4VEmModel(nam),isInitialised(false)
|
|---|
| 41 | {
|
|---|
| 42 |
|
|---|
| 43 | killBelowEnergy = 8.23*eV; // Minimum e- energy for energy loss by excitation
|
|---|
| 44 | lowEnergyLimit = 0 * eV;
|
|---|
| 45 | lowEnergyLimitOfModel = 7.4 * eV; // The model lower energy is 7.4 eV
|
|---|
| 46 | highEnergyLimit = 10 * MeV;
|
|---|
| 47 | SetLowEnergyLimit(lowEnergyLimit);
|
|---|
| 48 | SetHighEnergyLimit(highEnergyLimit);
|
|---|
| 49 |
|
|---|
| 50 | verboseLevel= 0;
|
|---|
| 51 | // Verbosity scale:
|
|---|
| 52 | // 0 = nothing
|
|---|
| 53 | // 1 = warning for energy non-conservation
|
|---|
| 54 | // 2 = details of energy budget
|
|---|
| 55 | // 3 = calculation of cross sections, file openings, sampling of atoms
|
|---|
| 56 | // 4 = entering in methods
|
|---|
| 57 |
|
|---|
| 58 | if( verboseLevel>0 )
|
|---|
| 59 | {
|
|---|
| 60 | G4cout << "Champion Elastic model is constructed " << G4endl
|
|---|
| 61 | << "Energy range: "
|
|---|
| 62 | << lowEnergyLimit / eV << " eV - "
|
|---|
| 63 | << highEnergyLimit / MeV << " MeV"
|
|---|
| 64 | << G4endl;
|
|---|
| 65 | }
|
|---|
| 66 | }
|
|---|
| 67 |
|
|---|
| 68 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 69 |
|
|---|
| 70 | G4DNAChampionElasticModel::~G4DNAChampionElasticModel()
|
|---|
| 71 | {
|
|---|
| 72 | // For total cross section
|
|---|
| 73 |
|
|---|
| 74 | std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
|
|---|
| 75 | for (pos = tableData.begin(); pos != tableData.end(); ++pos)
|
|---|
| 76 | {
|
|---|
| 77 | G4DNACrossSectionDataSet* table = pos->second;
|
|---|
| 78 | delete table;
|
|---|
| 79 | }
|
|---|
| 80 |
|
|---|
| 81 | // For final state
|
|---|
| 82 |
|
|---|
| 83 | eVecm.clear();
|
|---|
| 84 |
|
|---|
| 85 | }
|
|---|
| 86 |
|
|---|
| 87 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 88 |
|
|---|
| 89 | void G4DNAChampionElasticModel::Initialise(const G4ParticleDefinition* /*particle*/,
|
|---|
| 90 | const G4DataVector& /*cuts*/)
|
|---|
| 91 | {
|
|---|
| 92 |
|
|---|
| 93 | if (verboseLevel > 3)
|
|---|
| 94 | G4cout << "Calling G4DNAChampionElasticModel::Initialise()" << G4endl;
|
|---|
| 95 |
|
|---|
| 96 | // Energy limits
|
|---|
| 97 |
|
|---|
| 98 | if (LowEnergyLimit() < lowEnergyLimit)
|
|---|
| 99 | {
|
|---|
| 100 | G4cout << "G4DNAChampionElasticModel: low energy limit increased from " <<
|
|---|
| 101 | LowEnergyLimit()/eV << " eV to " << lowEnergyLimit/eV << " eV" << G4endl;
|
|---|
| 102 | SetLowEnergyLimit(lowEnergyLimit);
|
|---|
| 103 | }
|
|---|
| 104 |
|
|---|
| 105 | if (HighEnergyLimit() > highEnergyLimit)
|
|---|
| 106 | {
|
|---|
| 107 | G4cout << "G4DNAChampionElasticModel: high energy limit decreased from " <<
|
|---|
| 108 | HighEnergyLimit()/MeV << " MeV to " << highEnergyLimit/MeV << " MeV" << G4endl;
|
|---|
| 109 | SetHighEnergyLimit(highEnergyLimit);
|
|---|
| 110 | }
|
|---|
| 111 |
|
|---|
| 112 | // Reading of data files
|
|---|
| 113 |
|
|---|
| 114 | G4double scaleFactor = 1e-16*cm*cm;
|
|---|
| 115 |
|
|---|
| 116 | G4String fileElectron("dna/sigma_elastic_e_champion");
|
|---|
| 117 |
|
|---|
| 118 | G4ParticleDefinition* electronDef = G4Electron::ElectronDefinition();
|
|---|
| 119 | G4String electron;
|
|---|
| 120 |
|
|---|
| 121 | if (electronDef != 0)
|
|---|
| 122 | {
|
|---|
| 123 | // For total cross section
|
|---|
| 124 |
|
|---|
| 125 | electron = electronDef->GetParticleName();
|
|---|
| 126 |
|
|---|
| 127 | tableFile[electron] = fileElectron;
|
|---|
| 128 |
|
|---|
| 129 | G4DNACrossSectionDataSet* tableE = new G4DNACrossSectionDataSet(new G4LogLogInterpolation, eV,scaleFactor );
|
|---|
| 130 | tableE->LoadData(fileElectron);
|
|---|
| 131 | tableData[electron] = tableE;
|
|---|
| 132 |
|
|---|
| 133 | // For final state
|
|---|
| 134 |
|
|---|
| 135 | char *path = getenv("G4LEDATA");
|
|---|
| 136 |
|
|---|
| 137 | if (!path)
|
|---|
| 138 | G4Exception("G4FinalStateElasticChampion::Initialise: G4LEDATA environment variable not set");
|
|---|
| 139 |
|
|---|
| 140 | std::ostringstream eFullFileName;
|
|---|
| 141 | eFullFileName << path << "/dna/sigmadiff_elastic_e_champion.dat";
|
|---|
| 142 | std::ifstream eDiffCrossSection(eFullFileName.str().c_str());
|
|---|
| 143 |
|
|---|
| 144 | if (!eDiffCrossSection) G4Exception("G4DNAChampionElasticModel::Initialise: error opening electron DATA FILE");
|
|---|
| 145 |
|
|---|
| 146 | eTdummyVec.push_back(0.);
|
|---|
| 147 |
|
|---|
| 148 | while(!eDiffCrossSection.eof())
|
|---|
| 149 | {
|
|---|
| 150 | double tDummy;
|
|---|
| 151 | double eDummy;
|
|---|
| 152 | eDiffCrossSection>>tDummy>>eDummy;
|
|---|
| 153 |
|
|---|
| 154 | // SI : mandatory eVecm initialization
|
|---|
| 155 | if (tDummy != eTdummyVec.back())
|
|---|
| 156 | {
|
|---|
| 157 | eTdummyVec.push_back(tDummy);
|
|---|
| 158 | eVecm[tDummy].push_back(0.);
|
|---|
| 159 | }
|
|---|
| 160 |
|
|---|
| 161 | eDiffCrossSection>>eDiffCrossSectionData[tDummy][eDummy];
|
|---|
| 162 |
|
|---|
| 163 | // SI : only if not end of file reached !
|
|---|
| 164 | if (!eDiffCrossSection.eof()) eDiffCrossSectionData[tDummy][eDummy]*=scaleFactor;
|
|---|
| 165 |
|
|---|
| 166 | if (eDummy != eVecm[tDummy].back()) eVecm[tDummy].push_back(eDummy);
|
|---|
| 167 |
|
|---|
| 168 | }
|
|---|
| 169 |
|
|---|
| 170 | // End final state
|
|---|
| 171 |
|
|---|
| 172 | }
|
|---|
| 173 | else G4Exception("G4DNAChampionElasticModel::Initialise: electron is not defined");
|
|---|
| 174 |
|
|---|
| 175 | if (verboseLevel > 2)
|
|---|
| 176 | G4cout << "Loaded cross section files for Champion Elastic model" << G4endl;
|
|---|
| 177 |
|
|---|
| 178 | if( verboseLevel>0 )
|
|---|
| 179 | {
|
|---|
| 180 | G4cout << "Champion Elastic model is initialized " << G4endl
|
|---|
| 181 | << "Energy range: "
|
|---|
| 182 | << LowEnergyLimit() / eV << " eV - "
|
|---|
| 183 | << HighEnergyLimit() / MeV << " MeV"
|
|---|
| 184 | << G4endl;
|
|---|
| 185 | }
|
|---|
| 186 |
|
|---|
| 187 | if(!isInitialised)
|
|---|
| 188 | {
|
|---|
| 189 | isInitialised = true;
|
|---|
| 190 |
|
|---|
| 191 | if(pParticleChange)
|
|---|
| 192 | fParticleChangeForGamma = reinterpret_cast<G4ParticleChangeForGamma*>(pParticleChange);
|
|---|
| 193 | else
|
|---|
| 194 | fParticleChangeForGamma = new G4ParticleChangeForGamma();
|
|---|
| 195 | }
|
|---|
| 196 |
|
|---|
| 197 | // InitialiseElementSelectors(particle,cuts);
|
|---|
| 198 |
|
|---|
| 199 | // Test if water material
|
|---|
| 200 |
|
|---|
| 201 | flagMaterialIsWater= false;
|
|---|
| 202 | densityWater = 0;
|
|---|
| 203 |
|
|---|
| 204 | const G4ProductionCutsTable* theCoupleTable = G4ProductionCutsTable::GetProductionCutsTable();
|
|---|
| 205 |
|
|---|
| 206 | if(theCoupleTable)
|
|---|
| 207 | {
|
|---|
| 208 | G4int numOfCouples = theCoupleTable->GetTableSize();
|
|---|
| 209 |
|
|---|
| 210 | if(numOfCouples>0)
|
|---|
| 211 | {
|
|---|
| 212 | for (G4int i=0; i<numOfCouples; i++)
|
|---|
| 213 | {
|
|---|
| 214 | const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
|
|---|
| 215 | const G4Material* material = couple->GetMaterial();
|
|---|
| 216 |
|
|---|
| 217 | if (material->GetName() == "G4_WATER")
|
|---|
| 218 | {
|
|---|
| 219 | G4double density = material->GetAtomicNumDensityVector()[1];
|
|---|
| 220 | flagMaterialIsWater = true;
|
|---|
| 221 | densityWater = density;
|
|---|
| 222 |
|
|---|
| 223 | if (verboseLevel > 3)
|
|---|
| 224 | G4cout << "****** Water material is found with density(cm^-3)=" << density/(cm*cm*cm) << G4endl;
|
|---|
| 225 | }
|
|---|
| 226 |
|
|---|
| 227 | }
|
|---|
| 228 |
|
|---|
| 229 | } // if(numOfCouples>0)
|
|---|
| 230 |
|
|---|
| 231 | } // if (theCoupleTable)
|
|---|
| 232 |
|
|---|
| 233 | }
|
|---|
| 234 |
|
|---|
| 235 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 236 |
|
|---|
| 237 | G4double G4DNAChampionElasticModel::CrossSectionPerVolume(const G4Material*,
|
|---|
| 238 | const G4ParticleDefinition* p,
|
|---|
| 239 | G4double ekin,
|
|---|
| 240 | G4double,
|
|---|
| 241 | G4double)
|
|---|
| 242 | {
|
|---|
| 243 | if (verboseLevel > 3)
|
|---|
| 244 | G4cout << "Calling CrossSectionPerVolume() of G4DNAChampionElasticModel" << G4endl;
|
|---|
| 245 |
|
|---|
| 246 | // Calculate total cross section for model
|
|---|
| 247 |
|
|---|
| 248 | G4double sigma=0;
|
|---|
| 249 |
|
|---|
| 250 | if (flagMaterialIsWater)
|
|---|
| 251 | {
|
|---|
| 252 | const G4String& particleName = p->GetParticleName();
|
|---|
| 253 |
|
|---|
| 254 | if (ekin < highEnergyLimit)
|
|---|
| 255 | {
|
|---|
| 256 | //SI : XS must not be zero otherwise sampling of secondaries method ignored
|
|---|
| 257 | if (ekin < lowEnergyLimitOfModel) ekin = lowEnergyLimitOfModel;
|
|---|
| 258 | //
|
|---|
| 259 |
|
|---|
| 260 | std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
|
|---|
| 261 | pos = tableData.find(particleName);
|
|---|
| 262 |
|
|---|
| 263 | if (pos != tableData.end())
|
|---|
| 264 | {
|
|---|
| 265 | G4DNACrossSectionDataSet* table = pos->second;
|
|---|
| 266 | if (table != 0)
|
|---|
| 267 | {
|
|---|
| 268 | sigma = table->FindValue(ekin);
|
|---|
| 269 | }
|
|---|
| 270 | }
|
|---|
| 271 | else
|
|---|
| 272 | {
|
|---|
| 273 | G4Exception("G4DNAChampionElasticModel::ComputeCrossSectionPerVolume: attempting to calculate cross section for wrong particle");
|
|---|
| 274 | }
|
|---|
| 275 | }
|
|---|
| 276 |
|
|---|
| 277 | if (verboseLevel > 3)
|
|---|
| 278 | {
|
|---|
| 279 | G4cout << "---> Kinetic energy(eV)=" << ekin/eV << G4endl;
|
|---|
| 280 | G4cout << " - Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
|
|---|
| 281 | G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*densityWater/(1./cm) << G4endl;
|
|---|
| 282 | }
|
|---|
| 283 |
|
|---|
| 284 | } // if (flagMaterialIsWater)
|
|---|
| 285 |
|
|---|
| 286 | return sigma*densityWater;
|
|---|
| 287 | }
|
|---|
| 288 |
|
|---|
| 289 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 290 |
|
|---|
| 291 | void G4DNAChampionElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
|
|---|
| 292 | const G4MaterialCutsCouple* /*couple*/,
|
|---|
| 293 | const G4DynamicParticle* aDynamicElectron,
|
|---|
| 294 | G4double,
|
|---|
| 295 | G4double)
|
|---|
| 296 | {
|
|---|
| 297 |
|
|---|
| 298 | if (verboseLevel > 3)
|
|---|
| 299 | G4cout << "Calling SampleSecondaries() of G4DNAChampionElasticModel" << G4endl;
|
|---|
| 300 |
|
|---|
| 301 | G4double electronEnergy0 = aDynamicElectron->GetKineticEnergy();
|
|---|
| 302 |
|
|---|
| 303 | if (electronEnergy0 < killBelowEnergy)
|
|---|
| 304 | {
|
|---|
| 305 | fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
|
|---|
| 306 | fParticleChangeForGamma->ProposeLocalEnergyDeposit(electronEnergy0);
|
|---|
| 307 | return ;
|
|---|
| 308 | }
|
|---|
| 309 |
|
|---|
| 310 | if (electronEnergy0>= killBelowEnergy && electronEnergy0 < highEnergyLimit)
|
|---|
| 311 | {
|
|---|
| 312 | G4double cosTheta = RandomizeCosTheta(electronEnergy0);
|
|---|
| 313 |
|
|---|
| 314 | G4double phi = 2. * pi * G4UniformRand();
|
|---|
| 315 |
|
|---|
| 316 | G4ThreeVector zVers = aDynamicElectron->GetMomentumDirection();
|
|---|
| 317 | G4ThreeVector xVers = zVers.orthogonal();
|
|---|
| 318 | G4ThreeVector yVers = zVers.cross(xVers);
|
|---|
| 319 |
|
|---|
| 320 | G4double xDir = std::sqrt(1. - cosTheta*cosTheta);
|
|---|
| 321 | G4double yDir = xDir;
|
|---|
| 322 | xDir *= std::cos(phi);
|
|---|
| 323 | yDir *= std::sin(phi);
|
|---|
| 324 |
|
|---|
| 325 | G4ThreeVector zPrimeVers((xDir*xVers + yDir*yVers + cosTheta*zVers));
|
|---|
| 326 |
|
|---|
| 327 | fParticleChangeForGamma->ProposeMomentumDirection(zPrimeVers.unit()) ;
|
|---|
| 328 |
|
|---|
| 329 | fParticleChangeForGamma->SetProposedKineticEnergy(electronEnergy0);
|
|---|
| 330 | }
|
|---|
| 331 |
|
|---|
| 332 | }
|
|---|
| 333 |
|
|---|
| 334 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 335 |
|
|---|
| 336 | G4double G4DNAChampionElasticModel::DifferentialCrossSection
|
|---|
| 337 | (G4ParticleDefinition * particleDefinition, G4double k, G4double theta)
|
|---|
| 338 | {
|
|---|
| 339 |
|
|---|
| 340 | G4double sigma = 0.;
|
|---|
| 341 | G4double valueT1 = 0;
|
|---|
| 342 | G4double valueT2 = 0;
|
|---|
| 343 | G4double valueE21 = 0;
|
|---|
| 344 | G4double valueE22 = 0;
|
|---|
| 345 | G4double valueE12 = 0;
|
|---|
| 346 | G4double valueE11 = 0;
|
|---|
| 347 | G4double xs11 = 0;
|
|---|
| 348 | G4double xs12 = 0;
|
|---|
| 349 | G4double xs21 = 0;
|
|---|
| 350 | G4double xs22 = 0;
|
|---|
| 351 |
|
|---|
| 352 | //SI : ensure the correct computation of cross section at the 180*deg limit
|
|---|
| 353 | if (theta==180.) theta=theta-1e-9;
|
|---|
| 354 |
|
|---|
| 355 | if (particleDefinition == G4Electron::ElectronDefinition())
|
|---|
| 356 | {
|
|---|
| 357 | std::vector<double>::iterator t2 = std::upper_bound(eTdummyVec.begin(),eTdummyVec.end(), k);
|
|---|
| 358 | std::vector<double>::iterator t1 = t2-1;
|
|---|
| 359 |
|
|---|
| 360 | std::vector<double>::iterator e12 = std::upper_bound(eVecm[(*t1)].begin(),eVecm[(*t1)].end(), theta);
|
|---|
| 361 | std::vector<double>::iterator e11 = e12-1;
|
|---|
| 362 |
|
|---|
| 363 | std::vector<double>::iterator e22 = std::upper_bound(eVecm[(*t2)].begin(),eVecm[(*t2)].end(), theta);
|
|---|
| 364 | std::vector<double>::iterator e21 = e22-1;
|
|---|
| 365 |
|
|---|
| 366 | valueT1 =*t1;
|
|---|
| 367 | valueT2 =*t2;
|
|---|
| 368 | valueE21 =*e21;
|
|---|
| 369 | valueE22 =*e22;
|
|---|
| 370 | valueE12 =*e12;
|
|---|
| 371 | valueE11 =*e11;
|
|---|
| 372 |
|
|---|
| 373 | xs11 = eDiffCrossSectionData[valueT1][valueE11];
|
|---|
| 374 | xs12 = eDiffCrossSectionData[valueT1][valueE12];
|
|---|
| 375 | xs21 = eDiffCrossSectionData[valueT2][valueE21];
|
|---|
| 376 | xs22 = eDiffCrossSectionData[valueT2][valueE22];
|
|---|
| 377 |
|
|---|
| 378 | }
|
|---|
| 379 |
|
|---|
| 380 | G4double xsProduct = xs11 * xs12 * xs21 * xs22;
|
|---|
| 381 |
|
|---|
| 382 | if (xs11==0 || xs12==0 ||xs21==0 ||xs22==0) return (0.);
|
|---|
| 383 |
|
|---|
| 384 | if (xsProduct != 0.)
|
|---|
| 385 | {
|
|---|
| 386 | sigma = QuadInterpolator( valueE11, valueE12,
|
|---|
| 387 | valueE21, valueE22,
|
|---|
| 388 | xs11, xs12,
|
|---|
| 389 | xs21, xs22,
|
|---|
| 390 | valueT1, valueT2,
|
|---|
| 391 | k, theta );
|
|---|
| 392 | }
|
|---|
| 393 |
|
|---|
| 394 | return sigma;
|
|---|
| 395 | }
|
|---|
| 396 |
|
|---|
| 397 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 398 |
|
|---|
| 399 | G4double G4DNAChampionElasticModel::LinLogInterpolate(G4double e1,
|
|---|
| 400 | G4double e2,
|
|---|
| 401 | G4double e,
|
|---|
| 402 | G4double xs1,
|
|---|
| 403 | G4double xs2)
|
|---|
| 404 | {
|
|---|
| 405 | G4double d1 = std::log(xs1);
|
|---|
| 406 | G4double d2 = std::log(xs2);
|
|---|
| 407 | G4double value = std::exp(d1 + (d2 - d1)*(e - e1)/ (e2 - e1));
|
|---|
| 408 | return value;
|
|---|
| 409 | }
|
|---|
| 410 |
|
|---|
| 411 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 412 |
|
|---|
| 413 | G4double G4DNAChampionElasticModel::LogLogInterpolate(G4double e1,
|
|---|
| 414 | G4double e2,
|
|---|
| 415 | G4double e,
|
|---|
| 416 | G4double xs1,
|
|---|
| 417 | G4double xs2)
|
|---|
| 418 | {
|
|---|
| 419 | G4double a = (std::log10(xs2)-std::log10(xs1)) / (std::log10(e2)-std::log10(e1));
|
|---|
| 420 | G4double b = std::log10(xs2) - a*std::log10(e2);
|
|---|
| 421 | G4double sigma = a*std::log10(e) + b;
|
|---|
| 422 | G4double value = (std::pow(10.,sigma));
|
|---|
| 423 | return value;
|
|---|
| 424 | }
|
|---|
| 425 |
|
|---|
| 426 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 427 |
|
|---|
| 428 | G4double G4DNAChampionElasticModel::QuadInterpolator(G4double e11, G4double e12,
|
|---|
| 429 | G4double e21, G4double e22,
|
|---|
| 430 | G4double xs11, G4double xs12,
|
|---|
| 431 | G4double xs21, G4double xs22,
|
|---|
| 432 | G4double t1, G4double t2,
|
|---|
| 433 | G4double t, G4double e)
|
|---|
| 434 | {
|
|---|
| 435 | // Log-Log
|
|---|
| 436 | /*
|
|---|
| 437 | G4double interpolatedvalue1 = LogLogInterpolate(e11, e12, e, xs11, xs12);
|
|---|
| 438 | G4double interpolatedvalue2 = LogLogInterpolate(e21, e22, e, xs21, xs22);
|
|---|
| 439 | G4double value = LogLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
|
|---|
| 440 | */
|
|---|
| 441 |
|
|---|
| 442 | // Lin-Log
|
|---|
| 443 | G4double interpolatedvalue1 = LinLogInterpolate(e11, e12, e, xs11, xs12);
|
|---|
| 444 | G4double interpolatedvalue2 = LinLogInterpolate(e21, e22, e, xs21, xs22);
|
|---|
| 445 | G4double value = LinLogInterpolate(t1, t2, t, interpolatedvalue1, interpolatedvalue2);
|
|---|
| 446 | return value;
|
|---|
| 447 | }
|
|---|
| 448 |
|
|---|
| 449 | //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
|
|---|
| 450 |
|
|---|
| 451 | G4double G4DNAChampionElasticModel::RandomizeCosTheta(G4double k)
|
|---|
| 452 | {
|
|---|
| 453 | // ***** Similar method as for screened Rutherford scattering
|
|---|
| 454 |
|
|---|
| 455 | G4int iMax=180;
|
|---|
| 456 | G4double max=0;
|
|---|
| 457 | G4double tmp=0;
|
|---|
| 458 |
|
|---|
| 459 | // Look for maximum :
|
|---|
| 460 | for (G4int i=0; i<iMax; i++)
|
|---|
| 461 | {
|
|---|
| 462 | tmp = DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,G4double(i)*180./(iMax-1));
|
|---|
| 463 | if (tmp>max) max = tmp;
|
|---|
| 464 | }
|
|---|
| 465 |
|
|---|
| 466 | G4double oneOverMax=0;
|
|---|
| 467 | if (max!=0) oneOverMax = 1./max;
|
|---|
| 468 |
|
|---|
| 469 | G4double cosTheta = 0.;
|
|---|
| 470 | G4double fCosTheta = 0.;
|
|---|
| 471 |
|
|---|
| 472 | do
|
|---|
| 473 | {
|
|---|
| 474 | cosTheta = 2. * G4UniformRand() - 1.;
|
|---|
| 475 | fCosTheta = oneOverMax * DifferentialCrossSection(G4Electron::ElectronDefinition(),k/eV,std::acos(cosTheta)*180./pi);
|
|---|
| 476 | }
|
|---|
| 477 | while (fCosTheta < G4UniformRand());
|
|---|
| 478 |
|
|---|
| 479 | if (verboseLevel > 3)
|
|---|
| 480 | {
|
|---|
| 481 | G4cout << "---> Cos(theta)=" << cosTheta << G4endl;
|
|---|
| 482 | }
|
|---|
| 483 |
|
|---|
| 484 | return cosTheta;
|
|---|
| 485 | }
|
|---|