Geant4 9.6.0
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 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=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
 

Static Public Member Functions

static const char * Default_Name ()
 

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 47 of file G4ChipsNeutronElasticXS.hh.

Constructor & Destructor Documentation

◆ G4ChipsNeutronElasticXS()

G4ChipsNeutronElasticXS::G4ChipsNeutronElasticXS ( )

Definition at line 54 of file G4ChipsNeutronElasticXS.cc.

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

◆ ~G4ChipsNeutronElasticXS()

G4ChipsNeutronElasticXS::~G4ChipsNeutronElasticXS ( )

Definition at line 94 of file G4ChipsNeutronElasticXS.cc.

95{
96 std::vector<G4double*>::iterator pos;
97 for (pos=CST.begin(); pos<CST.end(); pos++)
98 { delete [] *pos; }
99 CST.clear();
100 for (pos=PAR.begin(); pos<PAR.end(); pos++)
101 { delete [] *pos; }
102 PAR.clear();
103 for (pos=SST.begin(); pos<SST.end(); pos++)
104 { delete [] *pos; }
105 SST.clear();
106 for (pos=S1T.begin(); pos<S1T.end(); pos++)
107 { delete [] *pos; }
108 S1T.clear();
109 for (pos=B1T.begin(); pos<B1T.end(); pos++)
110 { delete [] *pos; }
111 B1T.clear();
112 for (pos=S2T.begin(); pos<S2T.end(); pos++)
113 { delete [] *pos; }
114 S2T.clear();
115 for (pos=B2T.begin(); pos<B2T.end(); pos++)
116 { delete [] *pos; }
117 B2T.clear();
118 for (pos=S3T.begin(); pos<S3T.end(); pos++)
119 { delete [] *pos; }
120 S3T.clear();
121 for (pos=B3T.begin(); pos<B3T.end(); pos++)
122 { delete [] *pos; }
123 B3T.clear();
124 for (pos=S4T.begin(); pos<S4T.end(); pos++)
125 { delete [] *pos; }
126 S4T.clear();
127 for (pos=B4T.begin(); pos<B4T.end(); pos++)
128 { delete [] *pos; }
129 B4T.clear();
130}

Member Function Documentation

◆ 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 154 of file G4ChipsNeutronElasticXS.cc.

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

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 1853 of file G4ChipsNeutronElasticXS.cc.

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

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

◆ GetHMaxT()

G4double G4ChipsNeutronElasticXS::GetHMaxT ( )

Definition at line 1975 of file G4ChipsNeutronElasticXS.cc.

1976{
1977 static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
1978 return lastTM*HGeVSQ;
1979}

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 141 of file G4ChipsNeutronElasticXS.cc.

145{
146 G4double pMom=Pt->GetTotalMomentum();
147 G4int tgN = A - tgZ;
148
149 return GetChipsCrossSection(pMom, tgZ, tgN, 2112);
150}
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 132 of file G4ChipsNeutronElasticXS.cc.

135{
136 G4ParticleDefinition* particle = Pt->GetDefinition();
137 if (particle == G4Proton::Proton() ) return true;
138 return false;
139}
G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93

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