Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TripathiCrossSection Class Reference

#include <G4TripathiCrossSection.hh>

+ Inheritance diagram for G4TripathiCrossSection:

Public Member Functions

 G4TripathiCrossSection ()
 
 ~G4TripathiCrossSection ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *aPart, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 42 of file G4TripathiCrossSection.hh.

Constructor & Destructor Documentation

◆ G4TripathiCrossSection()

G4TripathiCrossSection::G4TripathiCrossSection ( )

Definition at line 43 of file G4TripathiCrossSection.cc.

◆ ~G4TripathiCrossSection()

G4TripathiCrossSection::~G4TripathiCrossSection ( )

Definition at line 47 of file G4TripathiCrossSection.cc.

48{}

Member Function Documentation

◆ GetElementCrossSection()

G4double G4TripathiCrossSection::GetElementCrossSection ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 62 of file G4TripathiCrossSection.cc.

65{
66 G4double result = 0.;
67 G4double targetAtomicNumber = G4NistManager::Instance()->GetAtomicMassAmu(ZZ);
68 G4double nTargetProtons = ZZ;
69
70 G4double kineticEnergy = aPart->GetKineticEnergy()/MeV;
71 G4double nProjProtons = aPart->GetDefinition()->GetPDGCharge();
72 G4double projectileAtomicNumber =
74
75 const G4double nuleonRadius=1.1E-15;
76 const G4double myNuleonRadius=1.36E-15;
77
78 // needs target mass
79 G4double targetMass =
81 ->GetIonMass(G4lrint(nTargetProtons), G4lrint(targetAtomicNumber));
82 G4LorentzVector pTarget(0,0,0,targetMass);
83 G4LorentzVector pProjectile(aPart->Get4Momentum());
84 pTarget = pTarget+pProjectile;
85 G4double E_cm = (pTarget.mag()-targetMass-pProjectile.m())/MeV;
86 if(E_cm <= DBL_MIN) { return result; }
87 // done
88 G4double r_rms_p = 0.6 * myNuleonRadius *
89 std::pow(projectileAtomicNumber, 1./3.);
90 G4double r_rms_t = 0.6 * myNuleonRadius *
91 std::pow(targetAtomicNumber, 1./3.);
92
93 // done
94 G4double r_p = 1.29*r_rms_p/nuleonRadius ;
95 G4double r_t = 1.29*r_rms_t/nuleonRadius;
96
97 // done
98 G4double Radius = r_p + r_t +
99 1.2*(std::pow(targetAtomicNumber, 1./3.) +
100 std::pow(projectileAtomicNumber, 1./3.))/std::pow(E_cm, 1./3.);
101
102 //done
103 G4double B = 1.44*nProjProtons*nTargetProtons/Radius;
104 if(E_cm <= B) return result;
105 // done
106 G4double Energy = kineticEnergy/projectileAtomicNumber;
107
108 // done
109 //
110 // Note that this correction to G4TripathiCrossSection is just to accurately
111 // reflect Tripathi's algorithm. However, if you're using alpha
112 // particles/protons consider using the more accurate
113 // G4TripathiLightCrossSection, which Tripathi developed specifically for
114 // light systems.
115 //
116
117 G4double D;
118 if (nProjProtons==1 && projectileAtomicNumber==1)
119 {
120 D = 2.05;
121 }
122 else if (nProjProtons==2 && projectileAtomicNumber==4)
123 {
124 D = 2.77-(8.0E-3*targetAtomicNumber)+
125 (1.8E-5*targetAtomicNumber*targetAtomicNumber)
126 - 0.8/(1+std::exp((250.-Energy)/75.));
127 }
128 else
129 {
130 //
131 // This is the original value used in the G4TripathiCrossSection
132 // implementation, and was used for all projectile/target conditions.
133 // I'm not touching this, although judging from Tripathi's paper, this is
134 // valid for cases where the nucleon density changes little with A.
135 //
136 D = 1.75;
137 }
138 // done
139 G4double C_E = D * (1-std::exp(-Energy/40.)) -
140 0.292*std::exp(-Energy/792.)*std::cos(0.229*std::pow(Energy, 0.453));
141
142 // done
143 G4double S = std::pow(projectileAtomicNumber, 1./3.)*
144 std::pow(targetAtomicNumber, 1./3.)/
145 (std::pow(projectileAtomicNumber, 1./3.) +
146 std::pow(targetAtomicNumber, 1./3.));
147
148 // done
149 G4double deltaE = 1.85*S + 0.16*S/std::pow(E_cm,1./3.) - C_E +
150 0.91*(targetAtomicNumber-2.*nTargetProtons)*nProjProtons/
151 (targetAtomicNumber*projectileAtomicNumber);
152
153 // done
154 result = pi * nuleonRadius*nuleonRadius *
155 std::pow(( std::pow(targetAtomicNumber, 1./3.) +
156 std::pow(projectileAtomicNumber, 1./3.) + deltaE),2.) *
157 (1-B/E_cm);
158
159 if(result < 0.) { result = 0.; }
160 return result*m2;
161
162}
double G4double
Definition: G4Types.hh:64
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
G4double GetPDGCharge() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()
const G4double pi
int G4lrint(double ad)
Definition: templates.hh:163
#define DBL_MIN
Definition: templates.hh:75

◆ IsElementApplicable()

G4bool G4TripathiCrossSection::IsElementApplicable ( const G4DynamicParticle aPart,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 51 of file G4TripathiCrossSection.cc.

53{
54 G4bool result = false;
55 if ( (aPart->GetDefinition()->GetBaryonNumber()>2.5) &&
56 ( aPart->GetKineticEnergy()/aPart->GetDefinition()->GetBaryonNumber()<1*GeV) ) {
57 result = true;
58 }
59 return result;
60}
bool G4bool
Definition: G4Types.hh:67

Referenced by G4GeneralSpaceNNCrossSection::GetElementCrossSection().


The documentation for this class was generated from the following files: