Geant4.10
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Data Structures | Typedefs | Functions
G4ScreenedNuclearRecoil.hh File Reference

Definition of the G4ScreenedNuclearRecoil class. More...

#include "globals.hh"
#include "G4VDiscreteProcess.hh"
#include "G4ParticleChange.hh"
#include "c2_function.hh"
#include "CLHEP/Units/SystemOfUnits.h"
#include <map>
#include <vector>

Go to the source code of this file.

Data Structures

struct  G4ScreeningTables
 
class  G4ScreenedCoulombCrossSectionInfo
 
class  G4ScreenedCoulombCrossSection
 
struct  G4CoulombKinematicsInfo
 
class  G4ScreenedCollisionStage
 
class  G4ScreenedCoulombClassicalKinematics
 
class  G4SingleScatter
 
class  G4ScreenedNuclearRecoil
 A process which handles screened Coulomb collisions between nuclei. More...
 
class  G4NativeScreenedCoulombCrossSection
 

Typedefs

typedef c2_const_ptr< G4doubleG4_c2_const_ptr
 
typedef c2_ptr< G4doubleG4_c2_ptr
 
typedef c2_function< G4doubleG4_c2_function
 
typedef struct G4ScreeningTables G4ScreeningTables
 
typedef struct
G4CoulombKinematicsInfo 
G4CoulombKinematicsInfo
 

Functions

G4_c2_functionZBLScreening (G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
 
G4_c2_functionMoliereScreening (G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
 
G4_c2_functionLJScreening (G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
 
G4_c2_functionLJZBLScreening (G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
 

Detailed Description

Definition of the G4ScreenedNuclearRecoil class.

Definition in file G4ScreenedNuclearRecoil.hh.

Typedef Documentation

Definition at line 70 of file G4ScreenedNuclearRecoil.hh.

Definition at line 74 of file G4ScreenedNuclearRecoil.hh.

Definition at line 73 of file G4ScreenedNuclearRecoil.hh.

Function Documentation

G4_c2_function& LJScreening ( G4int  z1,
G4int  z2,
size_t  npoints,
G4double  rMax,
G4double auval 
)

Definition at line 787 of file G4ScreenedNuclearRecoil.cc.

References python.hepunit::angstrom, and c2_factory< float_type >::lin_log_interpolating_function().

Referenced by G4NativeScreenedCoulombCrossSection::G4NativeScreenedCoulombCrossSection(), and LJZBLScreening().

788 {
789 //from Loftager, Besenbacher, Jensen & Sorensen
790 //PhysRev A20, 1443++, 1979
791  G4double au=0.8853*0.529*angstrom/std::sqrt(std::pow(z1, 0.6667)+std::pow(z2,0.6667));
792  std::vector<G4double> r(npoints), phi(npoints);
793 
794  for(size_t i=0; i<npoints; i++) {
795  G4double rr=(float)i/(float)(npoints-1);
796  r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
797 
798  G4double y=std::sqrt(9.67*r[i]/au);
799  G4double ysq=y*y;
800  G4double phipoly=1+y+0.3344*ysq+0.0485*y*ysq+0.002647*ysq*ysq;
801  phi[i]=phipoly*std::exp(-y);
802  // G4cout << r[i] << " " << phi[i] << G4endl;
803  }
804 
805  // compute the derivative at the origin for the spline
806  G4double logphiprime0=(9.67/2.0)*(2*0.3344-1.0); // #avoid 0/0 on first element
807  logphiprime0 *= (1.0/au); // #put back in natural units
808 
809  *auval=au;
810  return c2.lin_log_interpolating_function().load(r, phi, false, logphiprime0*phi[0],true,0);
811 }
int angstrom
Definition: hepunit.py:36
double G4double
Definition: G4Types.hh:76
G4_c2_function& LJZBLScreening ( G4int  z1,
G4int  z2,
size_t  npoints,
G4double  rMax,
G4double auval 
)

connector in between. These numbers are selected so the switchover

Definition at line 813 of file G4ScreenedNuclearRecoil.cc.

References c2_piecewise_function_p< float_type >::append_function(), c2_factory< float_type >::connector_function(), LJScreening(), c2_factory< float_type >::piecewise_function(), c2_const_ptr< float_type >::release_for_return(), c2_function< float_type >::set_domain(), c2_function< float_type >::xmax(), c2_function< float_type >::xmin(), and ZBLScreening().

Referenced by G4NativeScreenedCoulombCrossSection::G4NativeScreenedCoulombCrossSection().

814 {
815 // hybrid of LJ and ZBL, uses LJ if x < 0.25*auniv, ZBL if x > 1.5*auniv, and
816 /// connector in between. These numbers are selected so the switchover
817 // is very near the point where the functions naturally cross.
818  G4double auzbl, aulj;
819 
820  c2p zbl=ZBLScreening(z1, z2, npoints, rMax, &auzbl);
821  c2p lj=LJScreening(z1, z2, npoints, rMax, &aulj);
822 
823  G4double au=(auzbl+aulj)*0.5;
824  lj->set_domain(lj->xmin(), 0.25*au);
825  zbl->set_domain(1.5*au,zbl->xmax());
826 
827  c2p conn=c2.connector_function(lj->xmax(), lj, zbl->xmin(), zbl, true,0);
828  c2_piecewise_function_p<G4double> &pw=c2.piecewise_function();
829  c2p keepit(pw);
830  pw.append_function(lj);
831  pw.append_function(conn);
832  pw.append_function(zbl);
833 
834  *auval=au;
835  keepit.release_for_return();
836  return pw;
837 }
G4_c2_function & LJScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
void append_function(const c2_function< float_type > &func)
append a new function to the sequence
void set_domain(float_type amin, float_type amax)
set the domain for this function.
Definition: c2_function.hh:301
G4_c2_function & ZBLScreening(G4int z1, G4int z2, size_t npoints, G4double rMax, G4double *auval)
double G4double
Definition: G4Types.hh:76
float_type xmax() const
return the upper bound of the domain for this function as set by set_domain()
Definition: c2_function.hh:299
float_type xmin() const
return the lower bound of the domain for this function as set by set_domain()
Definition: c2_function.hh:297
create a c2_function which is a piecewise assembly of other c2_functions.The functions must have incr...
Definition: c2_function.hh:84
G4_c2_function& MoliereScreening ( G4int  z1,
G4int  z2,
size_t  npoints,
G4double  rMax,
G4double auval 
)

Definition at line 761 of file G4ScreenedNuclearRecoil.cc.

References python.hepunit::angstrom, and c2_factory< float_type >::lin_log_interpolating_function().

Referenced by G4NativeScreenedCoulombCrossSection::G4NativeScreenedCoulombCrossSection().

762 {
763  static const size_t ncoef=3;
764  static G4double scales[ncoef]={-6.0, -1.2, -0.3};
765  static G4double coefs[ncoef]={0.10, 0.55, 0.35};
766 
767  G4double au=0.8853*0.529*angstrom/std::sqrt(std::pow(z1, 0.6667)+std::pow(z2,0.6667));
768  std::vector<G4double> r(npoints), phi(npoints);
769 
770  for(size_t i=0; i<npoints; i++) {
771  G4double rr=(float)i/(float)(npoints-1);
772  r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
773  G4double sum=0.0;
774  for(size_t j=0; j<ncoef; j++) sum+=coefs[j]*std::exp(scales[j]*r[i]/au);
775  phi[i]=sum;
776  }
777 
778  // compute the derivative at the origin for the spline
779  G4double phiprime0=0.0;
780  for(size_t j=0; j<ncoef; j++) phiprime0+=scales[j]*coefs[j]*std::exp(scales[j]*r[0]/au);
781  phiprime0*=(1.0/au); // put back in natural units;
782 
783  *auval=au;
784  return c2.lin_log_interpolating_function().load(r, phi, false, phiprime0,true,0);
785 }
int angstrom
Definition: hepunit.py:36
double G4double
Definition: G4Types.hh:76
G4_c2_function& ZBLScreening ( G4int  z1,
G4int  z2,
size_t  npoints,
G4double  rMax,
G4double auval 
)

Definition at line 735 of file G4ScreenedNuclearRecoil.cc.

References python.hepunit::angstrom, and c2_factory< float_type >::lin_log_interpolating_function().

Referenced by G4NativeScreenedCoulombCrossSection::G4NativeScreenedCoulombCrossSection(), and LJZBLScreening().

736 {
737  static const size_t ncoef=4;
738  static G4double scales[ncoef]={-3.2, -0.9432, -0.4028, -0.2016};
739  static G4double coefs[ncoef]={0.1818,0.5099,0.2802,0.0281};
740 
741  G4double au=0.8854*angstrom*0.529/(std::pow(z1, 0.23)+std::pow(z2,0.23));
742  std::vector<G4double> r(npoints), phi(npoints);
743 
744  for(size_t i=0; i<npoints; i++) {
745  G4double rr=(float)i/(float)(npoints-1);
746  r[i]=rr*rr*rMax; // use quadratic r scale to make sampling fine near the center
747  G4double sum=0.0;
748  for(size_t j=0; j<ncoef; j++) sum+=coefs[j]*std::exp(scales[j]*r[i]/au);
749  phi[i]=sum;
750  }
751 
752  // compute the derivative at the origin for the spline
753  G4double phiprime0=0.0;
754  for(size_t j=0; j<ncoef; j++) phiprime0+=scales[j]*coefs[j]*std::exp(scales[j]*r[0]/au);
755  phiprime0*=(1.0/au); // put back in natural units;
756 
757  *auval=au;
758  return c2.lin_log_interpolating_function().load(r, phi, false, phiprime0,true,0);
759 }
int angstrom
Definition: hepunit.py:36
double G4double
Definition: G4Types.hh:76