// // ******************************************************************** // * 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: G4CoulombScatteringModel.cc,v 1.49 2010/05/27 14:22:05 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-04-beta-cand-01 $ // // ------------------------------------------------------------------- // // GEANT4 Class file // // // File name: G4CoulombScatteringModel // // 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.10.06 V.Ivanchenko use inheritance from G4eCoulombScatteringModel // 09.10.07 V.Ivanchenko reorganized methods, add cut dependence in scattering off e- // 09.06.08 V.Ivanchenko SelectIsotope is moved to the base class // 16.06.09 Consolandi rows 109, 111-112, 183, 185-186 // 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 "G4CoulombScatteringModel.hh" #include "Randomize.hh" #include "G4ParticleChangeForGamma.hh" #include "G4ParticleTable.hh" #include "G4IonTable.hh" #include "G4Proton.hh" #include "G4NucleiProperties.hh" //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... using namespace std; G4CoulombScatteringModel::G4CoulombScatteringModel(const G4String& nam) : G4eCoulombScatteringModel(nam) {} //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... G4CoulombScatteringModel::~G4CoulombScatteringModel() {} //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... G4double G4CoulombScatteringModel::ComputeCrossSectionPerAtom( const G4ParticleDefinition* p, G4double kinEnergy, G4double Z, G4double, G4double cutEnergy, G4double) { //G4cout << "### G4CoulombScatteringModel::ComputeCrossSectionPerAtom for " // << p->GetParticleName()<<" Z= "<size() << G4endl; */ if(trec > tcut) { G4ParticleDefinition* ion = theParticleTable->FindIon(iz, ia, 0, iz); G4double plab = sqrt(finalT*(finalT + 2.0*mass)); G4ThreeVector p2 = (ptot*dir - plab*newDirection).unit(); G4DynamicParticle* newdp = new G4DynamicParticle(ion, p2, trec); fvect->push_back(newdp); } else if(trec > 0.0) { fParticleChange->ProposeLocalEnergyDeposit(trec); fParticleChange->ProposeNonIonizingEnergyDeposit(trec); } } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......