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

#include <G4ChipsNeutronElasticXS.hh>

+ Inheritance diagram for G4ChipsNeutronElasticXS:

Public Member Functions

 G4ChipsNeutronElasticXS ()
 
 ~G4ChipsNeutronElasticXS ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4double GetChipsCrossSection (G4double momentum, G4int Z, G4int N, G4int pdg)
 
G4double GetExchangeT (G4int tZ, G4int tN, G4int pPDG)
 
G4double GetHMaxT ()
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double ComputeCrossSectionPerElement (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, const G4Element *, const G4Material *mat=nullptr)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=nullptr)
 
virtual G4double ComputeIsoCrossSection (G4double kinEnergy, G4double loge, const G4ParticleDefinition *, G4int Z, G4int A, const G4Isotope *iso=nullptr, const G4Element *elm=nullptr, const G4Material *mat=nullptr)
 
virtual const G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy, G4double logE)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
bool ForAllAtomsAndEnergies () const
 
void SetForAllAtomsAndEnergies (G4bool val)
 
const G4StringGetName () const
 
void SetName (const G4String &nam)
 
G4VCrossSectionDataSetoperator= (const G4VCrossSectionDataSet &right)=delete
 
 G4VCrossSectionDataSet (const G4VCrossSectionDataSet &)=delete
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 
G4String name
 

Detailed Description

Definition at line 46 of file G4ChipsNeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsNeutronElasticXS()

G4ChipsNeutronElasticXS::G4ChipsNeutronElasticXS ( )

Definition at line 61 of file G4ChipsNeutronElasticXS.cc.

61 :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
62{
63 lPMin=-8.; // Min tabulated log Momentum (D)
64 lPMax= 8.; // Max tabulated log Momentum (D)
65 dlnP=(lPMax-lPMin)/nLast;// Log step in table (D)
66 onlyCS=true;// Flag to calc only CS (not Si/Bi)(L)
67 lastSIG=0.; // Last calculated cross section (L)
68 lastLP=-10.;// Last log(momOfIncidentHadron) (L)
69 lastTM=0.; // Last t_maximum (L)
70 theSS=0.; // The Last sq.slope of 1st difMax (L)
71 theS1=0.; // The Last mantissa of 1st difMax (L)
72 theB1=0.; // The Last slope of 1st difr. Max (L)
73 theS2=0.; // The Last mantissa of 2nd difMax (L)
74 theB2=0.; // The Last slope of 2nd difr. Max (L)
75 theS3=0.; // The Last mantissa of 3d difrMax (L)
76 theB3=0.; // The Last slope of 3d difructMax (L)
77 theS4=0.; // The Last mantissa of 4th difMax (L)
78 theB4=0.; // The Last slope of 4th difr. Max (L)
79 lastTZ=0; // Last atomic number of the target
80 lastTN=0; // Last # of neutrons in the target
81 lastPIN=0.; // Last initialized max momentum
82 lastCST=0; // Elastic cross-section table
83 lastPAR=0; // Parameters of FunctionalCalculation
84 lastSST=0; // E-dep of sq.slope of the 1st difMax
85 lastS1T=0; // E-dep of mantissa of the 1st difMax
86 lastB1T=0; // E-dep of theSlope of the 1st difMax
87 lastS2T=0; // E-dep of mantissa of the 2nd difMax
88 lastB2T=0; // E-dep of theSlope of the 2nd difMax
89 lastS3T=0; // E-dep of mantissa of the 3d difrMax
90 lastB3T=0; // E-dep of the slope of the 3d difMax
91 lastS4T=0; // E-dep of mantissa of the 4th difMax
92 lastB4T=0; // E-dep of theSlope of the 4th difMax
93 lastN=0; // The last N of calculated nucleus
94 lastZ=0; // The last Z of calculated nucleus
95 lastP=0.; // Last used in cross section Momentum
96 lastTH=0.; // Last threshold momentum
97 lastCS=0.; // Last value of the Cross Section
98 lastI=0; // The last position in the DAMDB
99
100 mNeut= G4Neutron::Neutron()->GetPDGMass()*.001;// MeV to GeV
101 mProt= G4Proton::Proton()->GetPDGMass()*.001;// MeV to GeV
102 mNeut2= mNeut*mNeut;
103}
static const char * Default_Name()
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4Proton * Proton()
Definition G4Proton.cc:90
G4VCrossSectionDataSet(const G4String &nam="")

◆ ~G4ChipsNeutronElasticXS()

G4ChipsNeutronElasticXS::~G4ChipsNeutronElasticXS ( )

Definition at line 105 of file G4ChipsNeutronElasticXS.cc.

106{
107 std::vector<G4double*>::iterator pos;
108 for (pos=CST.begin(); pos<CST.end(); pos++)
109 { delete [] *pos; }
110 CST.clear();
111 for (pos=PAR.begin(); pos<PAR.end(); pos++)
112 { delete [] *pos; }
113 PAR.clear();
114 for (pos=SST.begin(); pos<SST.end(); pos++)
115 { delete [] *pos; }
116 SST.clear();
117 for (pos=S1T.begin(); pos<S1T.end(); pos++)
118 { delete [] *pos; }
119 S1T.clear();
120 for (pos=B1T.begin(); pos<B1T.end(); pos++)
121 { delete [] *pos; }
122 B1T.clear();
123 for (pos=S2T.begin(); pos<S2T.end(); pos++)
124 { delete [] *pos; }
125 S2T.clear();
126 for (pos=B2T.begin(); pos<B2T.end(); pos++)
127 { delete [] *pos; }
128 B2T.clear();
129 for (pos=S3T.begin(); pos<S3T.end(); pos++)
130 { delete [] *pos; }
131 S3T.clear();
132 for (pos=B3T.begin(); pos<B3T.end(); pos++)
133 { delete [] *pos; }
134 B3T.clear();
135 for (pos=S4T.begin(); pos<S4T.end(); pos++)
136 { delete [] *pos; }
137 S4T.clear();
138 for (pos=B4T.begin(); pos<B4T.end(); pos++)
139 { delete [] *pos; }
140 B4T.clear();
141}

Member Function Documentation

◆ CrossSectionDescription()

void G4ChipsNeutronElasticXS::CrossSectionDescription ( std::ostream & outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 144 of file G4ChipsNeutronElasticXS.cc.

145{
146 outFile << "G4ChipsNeutronElasticXS provides the elastic cross\n"
147 << "section for neutron nucleus scattering as a function of incident\n"
148 << "momentum. The cross section is calculated using M. Kossov's\n"
149 << "CHIPS parameterization of cross section data.\n";
150}

◆ Default_Name()

static const char * G4ChipsNeutronElasticXS::Default_Name ( )
inlinestatic

◆ GetChipsCrossSection()

G4double G4ChipsNeutronElasticXS::GetChipsCrossSection ( G4double momentum,
G4int Z,
G4int N,
G4int pdg )
virtual

!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)

Definition at line 172 of file G4ChipsNeutronElasticXS.cc.

173{
174
175 G4double pEn=pMom;
176 onlyCS=false;
177
178 G4bool in=false; // By default the isotope must be found in the AMDB
179 lastP = 0.; // New momentum history (nothing to compare with)
180 lastN = tgN; // The last N of the calculated nucleus
181 lastZ = tgZ; // The last Z of the calculated nucleus
182 lastI = (G4int)colN.size(); // Size of the Associative Memory DB in the heap
183 if(lastI) for(G4int i=0; i<lastI; ++i) // Loop over proj/tgZ/tgN lines of DB
184 { // The nucleus with projPDG is found in AMDB
185 if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
186 {
187 lastI=i;
188 lastTH =colTH[i]; // Last THreshold (A-dependent)
189 if(pEn<=lastTH)
190 {
191 return 0.; // Energy is below the Threshold value
192 }
193 lastP =colP [i]; // Last Momentum (A-dependent)
194 lastCS =colCS[i]; // Last CrossSect (A-dependent)
195 // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
196 if(lastP == pMom) // Do not recalculate
197 {
198 CalculateCrossSection(false,-1,i,2112,lastZ,lastN,pMom); // Update param's only
199 return lastCS*millibarn; // Use theLastCS
200 }
201 in = true; // This is the case when the isotop is found in DB
202 // Momentum pMom is in IU ! @@ Units
203 lastCS=CalculateCrossSection(false,-1,i,2112,lastZ,lastN,pMom); // read & update
204 if(lastCS<=0. && pEn>lastTH) // Correct the threshold
205 {
206 lastTH=pEn;
207 }
208 break; // Go out of the LOOP with found lastI
209 }
210 }
211 if(!in) // This nucleus has not been calculated previously
212 {
213 //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
214 lastCS=CalculateCrossSection(false,0,lastI,2112,lastZ,lastN,pMom);//calculate&create
215 if(lastCS<=0.)
216 {
217 lastTH = 0; // ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
218 if(pEn>lastTH)
219 {
220 lastTH=pEn;
221 }
222 }
223 colN.push_back(tgN);
224 colZ.push_back(tgZ);
225 colP.push_back(pMom);
226 colTH.push_back(lastTH);
227 colCS.push_back(lastCS);
228 return lastCS*millibarn;
229 } // End of creation of the new set of parameters
230 else
231 {
232 colP[lastI]=pMom;
233 colCS[lastI]=lastCS;
234 }
235 return lastCS*millibarn;
236}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85

Referenced by G4QuasiElRatios::ChExer(), G4ChipsComponentXS::GetElasticElementCrossSection(), GetIsoCrossSection(), G4ChipsComponentXS::GetTotalElementCrossSection(), G4ChipsElasticModel::SampleInvariantT(), and G4QuasiElRatios::Scatter().

◆ GetExchangeT()

G4double G4ChipsNeutronElasticXS::GetExchangeT ( G4int tZ,
G4int tN,
G4int pPDG )

Definition at line 1868 of file G4ChipsNeutronElasticXS.cc.

1869{
1870 static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
1871 static const G4double third=1./3.;
1872 static const G4double fifth=1./5.;
1873 static const G4double sevth=1./7.;
1874 if(PDG!=2112) G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetExT:PDG="<<PDG<<G4endl;
1875 if(onlyCS) G4cout<<"*Warning*G4ChipsNeutronElasticXS::GetExchangeT:onCS=1"<<G4endl;
1876 if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
1877 G4double q2=0.;
1878 if(tgZ==1 && tgN==0) // ===> n+p=n+p
1879 {
1880 G4double E1=lastTM*theB1;
1881 G4double R1=(1.-std::exp(-E1));
1882 G4double E2=lastTM*theB2;
1883 G4double R2=(1.-std::exp(-E2));
1884 G4double I1=R1*theS1;
1885 G4double I2=R2*theS2/theB2;
1886 //G4double I3=R3*theS3/theB3;
1887 G4double I12=I1+I2;
1888 //G4double rand=(I12+I3)*G4UniformRand();
1889 G4double rand=I12*G4UniformRand();
1890 if (rand<I1 )
1891 {
1892 G4double ran=R1*G4UniformRand();
1893 if(ran>1.) ran=1.;
1894 q2=-std::log(1.-ran)/theB1; // t-chan
1895 }
1896 else
1897 {
1898 G4double ran=R2*G4UniformRand();
1899 if(ran>1.) ran=1.;
1900 q2=lastTM+std::log(1.-ran)/theB2; // u-chan (ChEx)
1901 }
1902 }
1903 else
1904 {
1905 G4double a=tgZ+tgN;
1906 G4double E1=lastTM*(theB1+lastTM*theSS);
1907 G4double R1=(1.-std::exp(-E1));
1908 G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
1909 G4double tm2=lastTM*lastTM;
1910 G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
1911 if(a>6.5)E2*=tm2; // for heavy nuclei
1912 G4double R2=(1.-std::exp(-E2));
1913 G4double E3=lastTM*theB3;
1914 if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
1915 G4double R3=(1.-std::exp(-E3));
1916 G4double E4=lastTM*theB4;
1917 G4double R4=(1.-std::exp(-E4));
1918 G4double I1=R1*theS1;
1919 G4double I2=R2*theS2;
1920 G4double I3=R3*theS3;
1921 G4double I4=R4*theS4;
1922 G4double I12=I1+I2;
1923 G4double I13=I12+I3;
1924 G4double rand=(I13+I4)*G4UniformRand();
1925 if(rand<I1)
1926 {
1927 G4double ran=R1*G4UniformRand();
1928 if(ran>1.) ran=1.;
1929 q2=-std::log(1.-ran)/theB1;
1930 if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
1931 }
1932 else if(rand<I12)
1933 {
1934 G4double ran=R2*G4UniformRand();
1935 if(ran>1.) ran=1.;
1936 q2=-std::log(1.-ran)/theB2;
1937 if(q2<0.) q2=0.;
1938 if(a<6.5) q2=std::pow(q2,third);
1939 else q2=std::pow(q2,fifth);
1940 }
1941 else if(rand<I13)
1942 {
1943 G4double ran=R3*G4UniformRand();
1944 if(ran>1.) ran=1.;
1945 q2=-std::log(1.-ran)/theB3;
1946 if(q2<0.) q2=0.;
1947 if(a>6.5) q2=std::pow(q2,sevth);
1948 }
1949 else
1950 {
1951 G4double ran=R4*G4UniformRand();
1952 if(ran>1.) ran=1.;
1953 q2=-std::log(1.-ran)/theB4;
1954 if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
1955 }
1956 }
1957 if(q2<0.) q2=0.;
1958 if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QNeutronElCroSect::GetExchangeT: -t="<<q2<<G4endl;
1959 if(q2>lastTM)
1960 {
1961 q2=lastTM;
1962 }
1963 return q2*GeVSQ;
1964}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52

Referenced by G4QuasiElRatios::ChExer(), G4ChipsElasticModel::SampleInvariantT(), and G4QuasiElRatios::Scatter().

◆ GetHMaxT()

G4double G4ChipsNeutronElasticXS::GetHMaxT ( )

Definition at line 1990 of file G4ChipsNeutronElasticXS.cc.

1991{
1992 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
1993 return lastTM*HGeVSQ;
1994}

Referenced by G4QuasiElRatios::ChExer(), and G4QuasiElRatios::Scatter().

◆ GetIsoCrossSection()

G4double G4ChipsNeutronElasticXS::GetIsoCrossSection ( const G4DynamicParticle * Pt,
G4int tgZ,
G4int A,
const G4Isotope * iso = 0,
const G4Element * elm = 0,
const G4Material * mat = 0 )
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 159 of file G4ChipsNeutronElasticXS.cc.

163{
164 G4double pMom=Pt->GetTotalMomentum();
165 G4int tgN = A - tgZ;
166
167 return GetChipsCrossSection(pMom, tgZ, tgN, 2112);
168}
const G4double A[17]
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetTotalMomentum() const

◆ IsIsoApplicable()

G4bool G4ChipsNeutronElasticXS::IsIsoApplicable ( const G4DynamicParticle * Pt,
G4int Z,
G4int A,
const G4Element * elm,
const G4Material * mat )
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 152 of file G4ChipsNeutronElasticXS.cc.

155{
156 return true;
157}

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