// // ******************************************************************** // * 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. * // ******************************************************************** // // // // original by H.P. Wellisch // modified by J.L. Chuma, TRIUMF, 19-Nov-1996 // last modified: 27-Mar-1997 // J.P.Wellisch: 23-Apr-97: minor simplifications // modified by J.L.Chuma 24-Jul-97 to set the total momentum in Cinema and // EvaporationEffects // modified by J.L.Chuma 21-Oct-97 put std::abs() around the totalE^2-mass^2 // in calculation of total momentum in // Cinema and EvaporationEffects // Chr. Volcker, 10-Nov-1997: new methods and class variables. // HPW added utilities for low energy neutron transport. (12.04.1998) // M.G. Pia, 2 Oct 1998: modified GetFermiMomentum to avoid memory leaks #include "G4Nucleus.hh" #include "G4NucleiProperties.hh" #include "Randomize.hh" #include "G4HadronicException.hh" G4Nucleus::G4Nucleus() { pnBlackTrackEnergy = 0.0; dtaBlackTrackEnergy = 0.0; pnBlackTrackEnergyfromAnnihilation = 0.0; dtaBlackTrackEnergyfromAnnihilation = 0.0; excitationEnergy = 0.0; momentum = G4ThreeVector(0.,0.,0.); fermiMomentum = 1.52*hbarc/fermi; theTemp = 293.16*kelvin; } G4Nucleus::G4Nucleus( const G4double A, const G4double Z ) { SetParameters( A, Z ); pnBlackTrackEnergy = 0.0; dtaBlackTrackEnergy = 0.0; pnBlackTrackEnergyfromAnnihilation = 0.0; dtaBlackTrackEnergyfromAnnihilation = 0.0; excitationEnergy = 0.0; momentum = G4ThreeVector(0.,0.,0.); fermiMomentum = 1.52*hbarc/fermi; theTemp = 293.16*kelvin; } G4Nucleus::G4Nucleus( const G4Material *aMaterial ) { ChooseParameters( aMaterial ); pnBlackTrackEnergy = 0.0; dtaBlackTrackEnergy = 0.0; pnBlackTrackEnergyfromAnnihilation = 0.0; dtaBlackTrackEnergyfromAnnihilation = 0.0; excitationEnergy = 0.0; momentum = G4ThreeVector(0.,0.,0.); fermiMomentum = 1.52*hbarc/fermi; theTemp = aMaterial->GetTemperature(); } G4Nucleus::~G4Nucleus() {} G4ReactionProduct G4Nucleus:: GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp) const { G4double velMag = aVelocity.mag(); G4ReactionProduct result; G4double value = 0; G4double random = 1; G4double norm = 3.*std::sqrt(k_Boltzmann*temp*aMass*G4Neutron::Neutron()->GetPDGMass()); norm /= G4Neutron::Neutron()->GetPDGMass(); norm *= 5.; norm += velMag; norm /= velMag; while(value/normGetPDGMass()); G4double px, py, pz; px = GetThermalPz(theTarget.GetMass(), currentTemp); py = GetThermalPz(theTarget.GetMass(), currentTemp); pz = GetThermalPz(theTarget.GetMass(), currentTemp); theTarget.SetMomentum(px, py, pz); G4double tMom = std::sqrt(px*px+py*py+pz*pz); G4double tEtot = std::sqrt((tMom+theTarget.GetMass())* (tMom+theTarget.GetMass())- 2.*tMom*theTarget.GetMass()); if(1-tEtot/theTarget.GetMass()>0.001) { theTarget.SetTotalEnergy(tEtot); } else { theTarget.SetKineticEnergy(tMom*tMom/(2.*theTarget.GetMass())); } return theTarget; } void G4Nucleus::ChooseParameters( const G4Material *aMaterial ) { G4double random = G4UniformRand(); G4double sum = 0; const G4ElementVector *theElementVector = aMaterial->GetElementVector(); unsigned int i; for(i=0; iGetNumberOfElements(); ++i ) { sum += aMaterial->GetAtomicNumDensityVector()[i]; } G4double running = 0; for(i=0; iGetNumberOfElements(); ++i ) { running += aMaterial->GetAtomicNumDensityVector()[i]; if( running/sum > random ) { aEff = (*theElementVector)[i]->GetA()*mole/g; zEff = (*theElementVector)[i]->GetZ(); break; } } } void G4Nucleus::SetParameters( const G4double A, const G4double Z ) { G4int myZ = G4int(Z + 0.5); G4int myA = G4int(A + 0.5); if( myA<1 || myZ<0 || myZ>myA ) { throw G4HadronicException(__FILE__, __LINE__, "G4Nucleus::SetParameters called with non-physical parameters"); } aEff = A; // atomic weight zEff = Z; // atomic number } G4DynamicParticle * G4Nucleus::ReturnTargetParticle() const { // choose a proton or a neutron as the target particle G4DynamicParticle *targetParticle = new G4DynamicParticle; if( G4UniformRand() < zEff/aEff ) targetParticle->SetDefinition( G4Proton::Proton() ); else targetParticle->SetDefinition( G4Neutron::Neutron() ); return targetParticle; } G4double G4Nucleus::AtomicMass( const G4double A, const G4double Z ) const { // Now returns (atomic mass - electron masses) return G4NucleiProperties::GetNuclearMass(A, Z); } G4double G4Nucleus::GetThermalPz( const G4double mass, const G4double temp ) const { G4double result = G4RandGauss::shoot(); result *= std::sqrt(k_Boltzmann*temp*mass); // Das ist impuls (Pz), // nichtrelativistische rechnung // Maxwell verteilung angenommen return result; } G4double G4Nucleus::EvaporationEffects( G4double kineticEnergy ) { // derived from original FORTRAN code EXNU by H. Fesefeldt (10-Dec-1986) // // Nuclear evaporation as function of atomic number // and kinetic energy (MeV) of primary particle // // returns kinetic energy (MeV) // if( aEff < 1.5 ) { pnBlackTrackEnergy = dtaBlackTrackEnergy = 0.0; return 0.0; } G4double ek = kineticEnergy/GeV; G4float ekin = std::min( 4.0, std::max( 0.1, ek ) ); const G4float atno = std::min( 120., aEff ); const G4float gfa = 2.0*((aEff-1.0)/70.)*std::exp(-(aEff-1.0)/70.); // // 0.35 value at 1 GeV // 0.05 value at 0.1 GeV // G4float cfa = std::max( 0.15, 0.35 + ((0.35-0.05)/2.3)*std::log(ekin) ); G4float exnu = 7.716 * cfa * std::exp(-cfa) * ((atno-1.0)/120.)*std::exp(-(atno-1.0)/120.); G4float fpdiv = std::max( 0.5, 1.0-0.25*ekin*ekin ); // // pnBlackTrackEnergy is the kinetic energy (in GeV) available for // proton/neutron black track particles // dtaBlackTrackEnergy is the kinetic energy (in GeV) available for // deuteron/triton/alpha black track particles // pnBlackTrackEnergy = exnu*fpdiv; dtaBlackTrackEnergy = exnu*(1.0-fpdiv); if( G4int(zEff+0.1) != 82 ) { G4double ran1 = -6.0; G4double ran2 = -6.0; for( G4int i=0; i<12; ++i ) { ran1 += G4UniformRand(); ran2 += G4UniformRand(); } pnBlackTrackEnergy *= 1.0 + ran1*gfa; dtaBlackTrackEnergy *= 1.0 + ran2*gfa; } pnBlackTrackEnergy = std::max( 0.0, pnBlackTrackEnergy ); dtaBlackTrackEnergy = std::max( 0.0, dtaBlackTrackEnergy ); while( pnBlackTrackEnergy+dtaBlackTrackEnergy >= ek ) { pnBlackTrackEnergy *= 1.0 - 0.5*G4UniformRand(); dtaBlackTrackEnergy *= 1.0 - 0.5*G4UniformRand(); } // G4cout << "EvaporationEffects "<= ekOrg/GeV) { pnBlackTrackEnergyfromAnnihilation *= ekOrg/GeV/blackSum; dtaBlackTrackEnergyfromAnnihilation *= ekOrg/GeV/blackSum; } return (pnBlackTrackEnergyfromAnnihilation+dtaBlackTrackEnergyfromAnnihilation)*GeV; } G4double G4Nucleus::Cinema( G4double kineticEnergy ) { // derived from original FORTRAN code CINEMA by H. Fesefeldt (14-Oct-1987) // // input: kineticEnergy (MeV) // returns modified kinetic energy (MeV) // static const G4double expxu = 82.; // upper bound for arg. of exp static const G4double expxl = -expxu; // lower bound for arg. of exp G4double ek = kineticEnergy/GeV; G4double ekLog = std::log( ek ); G4double aLog = std::log( aEff ); G4double em = std::min( 1.0, 0.2390 + 0.0408*aLog*aLog ); G4double temp1 = -ek * std::min( 0.15, 0.0019*aLog*aLog*aLog ); G4double temp2 = std::exp( std::max( expxl, std::min( expxu, -(ekLog-em)*(ekLog-em)*2.0 ) ) ); G4double result = 0.0; if( std::abs( temp1 ) < 1.0 ) { if( temp2 > 1.0e-10 )result = temp1*temp2; } else result = temp1*temp2; if( result < -ek )result = -ek; return result*GeV; } // // methods for class G4Nucleus ... by Christian Volcker // G4ThreeVector G4Nucleus::GetFermiMomentum() { // chv: .. we assume zero temperature! // momentum is equally distributed in each phasespace volume dpx, dpy, dpz. G4double ranflat1= CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum); G4double ranflat2= CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum); G4double ranflat3= CLHEP::RandFlat::shoot((G4double)0.,(G4double)fermiMomentum); G4double ranmax = (ranflat1>ranflat2? ranflat1: ranflat2); ranmax = (ranmax>ranflat3? ranmax : ranflat3); // Isotropic momentum distribution G4double costheta = 2.*G4UniformRand() - 1.0; G4double sintheta = std::sqrt(1.0 - costheta*costheta); G4double phi = 2.0*pi*G4UniformRand(); G4double pz=costheta*ranmax; G4double px=sintheta*std::cos(phi)*ranmax; G4double py=sintheta*std::sin(phi)*ranmax; G4ThreeVector p(px,py,pz); return p; } G4ReactionProductVector* G4Nucleus::Fragmentate() { // needs implementation! return NULL; } void G4Nucleus::AddMomentum(const G4ThreeVector aMomentum) { momentum+=(aMomentum); } void G4Nucleus::AddExcitationEnergy( G4double anEnergy ) { excitationEnergy+=anEnergy; } /* end of file */