Geant4 Cross Reference |
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 // 27 // ------------------------------------------------------------------- 28 // GEANT 4 class implementation file 29 // 30 // CERN, Geneva, Switzerland 31 // 32 // File name: G4ProtonField.cc 33 // 34 // Author: Alessandro Brunengo (Alessandro.Brunengo@ge.infn.it) 35 // 36 // Creation date: 5 June 2000 37 // ------------------------------------------------------------------- 38 39 #include "G4ProtonField.hh" 40 #include "G4PhysicalConstants.hh" 41 #include "G4SystemOfUnits.hh" 42 #include "G4VNuclearDensity.hh" 43 #include "G4FermiMomentum.hh" 44 #include "G4V3DNucleus.hh" 45 #include "G4Pow.hh" 46 47 G4ProtonField::G4ProtonField(G4V3DNucleus * aNucleus) : 48 G4VNuclearField(aNucleus), theDensity(theNucleus->GetNuclearDensity()) 49 { 50 theA = theNucleus->GetMassNumber(); 51 theZ = theNucleus->GetCharge(); 52 theBarrier = GetBarrier(); 53 theRadius = 2.*theNucleus->GetOuterRadius(); 54 theFermi.Init(theA, theZ); 55 for (G4double aR=0.;aR<theRadius; aR+=0.3*fermi) 56 { 57 G4ThreeVector aPosition(0,0,aR); 58 G4double density = GetDensity(aPosition); 59 G4double fermiMom = GetFermiMomentum(density); 60 theFermiMomBuffer.push_back(fermiMom); 61 } 62 { 63 G4ThreeVector aPosition(0,0,theRadius); 64 G4double density = GetDensity(aPosition); 65 G4double fermiMom = GetFermiMomentum(density); 66 theFermiMomBuffer.push_back(fermiMom); 67 } 68 { 69 G4ThreeVector aPosition(0,0,theRadius+0.001*fermi); 70 theFermiMomBuffer.push_back(0); 71 } 72 { 73 G4ThreeVector aPosition(0,0,1.*m); 74 theFermiMomBuffer.push_back(0); 75 } 76 } 77 78 79 G4ProtonField::~G4ProtonField() 80 { } 81 82 G4double G4ProtonField::GetField(const G4ThreeVector & aPosition) 83 { 84 //G4cout << " Fermi Potential " << (fermiMom*fermiMom)/(2*proton_mass_c2) <<G4endl; 85 G4double x = aPosition.mag(); 86 unsigned int index = static_cast<unsigned int>(x/(0.3*fermi)); 87 if((index+2) > theFermiMomBuffer.size()) return theFermiMomBuffer.back(); 88 G4double y1 = theFermiMomBuffer[index]; 89 G4double y2 = theFermiMomBuffer[index+1]; 90 G4double x1 = (0.3*fermi)*index; 91 G4double x2 = (0.3*fermi)*(index+1); 92 G4double fermiMom = y1 + (x-x1)*(y2-y1)/(x2-x1); 93 G4double y = -1*(fermiMom*fermiMom)/(2*proton_mass_c2)+theBarrier; 94 // G4cout <<" Protonfield test "<<index<<" "<< x1<<" "<<y1<<" "<<x2<<" "<<y2<<" "<<x<<" "<<y<<" "<<theBarrier<<G4endl; 95 return y; 96 } 97 98 G4double G4ProtonField::GetBarrier() 99 { 100 G4double coulombBarrier = (1.44/1.14) * MeV * theZ / (1.0 + G4Pow::GetInstance()->Z13(theA)); 101 //GF G4double bindingEnergy = G4NucleiPropertiesTable::GetBindingEnergy(Z, A); 102 G4double bindingEnergy =0; 103 /* 104 * G4cout << " coulombBarrier/bindingEnergy : " 105 * << coulombBarrier << " /" << bindingEnergy << G4endl; 106 */ 107 return bindingEnergy/theA+coulombBarrier; 108 } 109