// // ******************************************************************** // * 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: G4eCoulombScatteringModel.cc,v 1.89 2010/05/27 14:22:05 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ // // ------------------------------------------------------------------- // // GEANT4 Class file // // // File name: G4eCoulombScatteringModel // // Author: Vladimir Ivanchenko // // Creation date: 22.08.2005 // // Modifications: // // 01.08.06 V.Ivanchenko extend upper limit of table to TeV and review the // logic of building - only elements from G4ElementTable // 08.08.06 V.Ivanchenko build internal table in ekin scale, introduce faclim // 19.08.06 V.Ivanchenko add inline function ScreeningParameter // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e- // 09.06.08 V.Ivanchenko add SelectIsotope and sampling of the recoil ion // 16.06.09 C.Consolandi fixed computation of effective mass // 27.05.10 V.Ivanchenko added G4WentzelOKandVIxSection class to // compute cross sections and sample scattering angle // // // Class Description: // // ------------------------------------------------------------------- // //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... #include "G4eCoulombScatteringModel.hh" #include "Randomize.hh" #include "G4DataVector.hh" #include "G4ElementTable.hh" #include "G4ParticleChangeForGamma.hh" #include "G4Proton.hh" #include "G4ParticleTable.hh" #include "G4ProductionCutsTable.hh" #include "G4NucleiProperties.hh" #include "G4Pow.hh" #include "G4LossTableManager.hh" #include "G4NistManager.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... using namespace std; G4eCoulombScatteringModel::G4eCoulombScatteringModel(const G4String& nam) : G4VEmModel(nam), cosThetaMin(1.0), cosThetaMax(-1.0), isInitialised(false) { fNistManager = G4NistManager::Instance(); theParticleTable = G4ParticleTable::GetParticleTable(); theProton = G4Proton::Proton(); currentMaterial = 0; currentElement = 0; lowEnergyLimit = 1*eV; recoilThreshold = 0.*keV; particle = 0; currentCouple = 0; wokvi = new G4WentzelOKandVIxSection(); } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... G4eCoulombScatteringModel::~G4eCoulombScatteringModel() { delete wokvi; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... void G4eCoulombScatteringModel::Initialise(const G4ParticleDefinition* p, const G4DataVector& cuts) { SetupParticle(p); currentCouple = 0; cosThetaMin = cos(PolarAngleLimit()); wokvi->Initialise(p, cosThetaMin); /* G4cout << "G4eCoulombScatteringModel: factorA2(GeV^2) = " << factorA2/(GeV*GeV) << " 1-cos(ThetaLimit)= " << 1 - cosThetaMin << " cos(thetaMax)= " << cosThetaMax << G4endl; */ pCuts = G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3); //G4cout << "!!! G4eCoulombScatteringModel::Initialise for " // << p->GetParticleName() << " cos(TetMin)= " << cosThetaMin // << " cos(TetMax)= " << cosThetaMax <GetParticleName()<<" Z= "<GetParticleName() // << " cut= " << cutEnergy<< G4endl; // Choose nucleus currentElement = SelectRandomAtom(couple,particle, kinEnergy,cutEnergy,kinEnergy); G4double Z = currentElement->GetZ(); if(ComputeCrossSectionPerAtom(particle,kinEnergy, Z, kinEnergy, cutEnergy, kinEnergy) == 0.0) { return; } G4int iz = G4int(Z); G4int ia = SelectIsotopeNumber(currentElement); G4double targetMass = G4NucleiProperties::GetNuclearMass(ia, iz); G4ThreeVector newDirection = wokvi->SampleSingleScattering(cosTetMinNuc, cosThetaMax, elecRatio); G4double cost = newDirection.z(); G4ThreeVector direction = dp->GetMomentumDirection(); newDirection.rotateUz(direction); fParticleChange->ProposeMomentumDirection(newDirection); // recoil sampling assuming a small recoil // and first order correction to primary 4-momentum G4double mom2 = wokvi->GetMomentumSquare(); G4double trec = mom2*(1.0 - cost)/(targetMass + (mass + kinEnergy)*(1.0 + cost)); G4double finalT = kinEnergy - trec; //G4cout<<"G4eCoulombScatteringModel: finalT= "<