Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ChipsPionMinusElasticXS.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4ChipsPionMinusElasticXS.cc 93260 2015-10-14 08:37:04Z gcosmo $
28 //
29 //
30 // G4 Physics class: G4ChipsPionMinusElasticXS for pA elastic cross sections
31 // Created: M.V. Kossov, CERN/ITEP(Moscow), 21-Jan-10
32 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 21-Jan-10
33 //
34 // -------------------------------------------------------------------------------
35 // Short description: Interaction cross-sections for the elastic process.
36 // Class extracted from CHIPS and integrated in Geant4 by W.Pokorski
37 // -------------------------------------------------------------------------------
38 //
39 
41 #include "G4SystemOfUnits.hh"
42 #include "G4DynamicParticle.hh"
43 #include "G4ParticleDefinition.hh"
44 #include "G4PionMinus.hh"
45 #include "G4Nucleus.hh"
46 #include "G4ParticleTable.hh"
47 #include "G4NucleiProperties.hh"
48 #include "G4IonTable.hh"
49 
50 #include "G4Log.hh"
51 #include "G4Exp.hh"
52 #include "G4Pow.hh"
53 
54 // factory
55 #include "G4CrossSectionFactory.hh"
56 //
58 
59 G4ChipsPionMinusElasticXS::G4ChipsPionMinusElasticXS():G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
60 {
61  lPMin=-8.; //Min tabulated logarithmMomentum(D)
62  lPMax= 8.; //Max tabulated logarithmMomentum(D)
63  dlnP=(lPMax-lPMin)/nLast;//LogStep inTheTable(D)
64  onlyCS=true;//Flag toCalcul OnlyCS(not Si/Bi)(L)
65  lastSIG=0.; //Last calculated cross section (L)
66  lastLP=-10.;//Last log(mom_of IncidentHadron)(L)
67  lastTM=0.; //Last t_maximum (L)
68  theSS=0.; //TheLastSqSlope of 1st difr.Max(L)
69  theS1=0.; //TheLastMantissa of 1st difr.Max(L)
70  theB1=0.; //TheLastSlope of 1st difruct.Max(L)
71  theS2=0.; //TheLastMantissa of 2nd difr.Max(L)
72  theB2=0.; //TheLastSlope of 2nd difruct.Max(L)
73  theS3=0.; //TheLastMantissa of 3d difr. Max(L)
74  theB3=0.; //TheLastSlope of 3d difruct. Max(L)
75  theS4=0.; //TheLastMantissa of 4th difr.Max(L)
76  theB4=0.; //TheLastSlope of 4th difruct.Max(L)
77  lastTZ=0; // Last atomic number of the target
78  lastTN=0; // Last # of neutrons in the target
79  lastPIN=0.; // Last initialized max momentum
80  lastCST=0; // Elastic cross-section table
81  lastPAR=0; // Parameters ForFunctionCalculation
82  lastSST=0; // E-dep of SqardSlope of 1st difMax
83  lastS1T=0; // E-dep of mantissa of 1st dif.Max
84  lastB1T=0; // E-dep of the slope of 1st difMax
85  lastS2T=0; // E-dep of mantissa of 2nd difrMax
86  lastB2T=0; // E-dep of the slope of 2nd difMax
87  lastS3T=0; // E-dep of mantissa of 3d difr.Max
88  lastB3T=0; // E-dep of the slope of 3d difrMax
89  lastS4T=0; // E-dep of mantissa of 4th difrMax
90  lastB4T=0; // E-dep of the slope of 4th difMax
91  lastN=0; // The last N of calculated nucleus
92  lastZ=0; // The last Z of calculated nucleus
93  lastP=0.; // LastUsed in CrossSection Momentum
94  lastTH=0.; // Last threshold momentum
95  lastCS=0.; // Last value of the Cross Section
96  lastI=0; // The last position in the DAMDB
97 }
98 
100 {
101  std::vector<G4double*>::iterator pos;
102  for (pos=CST.begin(); pos<CST.end(); pos++)
103  { delete [] *pos; }
104  CST.clear();
105  for (pos=PAR.begin(); pos<PAR.end(); pos++)
106  { delete [] *pos; }
107  PAR.clear();
108  for (pos=SST.begin(); pos<SST.end(); pos++)
109  { delete [] *pos; }
110  SST.clear();
111  for (pos=S1T.begin(); pos<S1T.end(); pos++)
112  { delete [] *pos; }
113  S1T.clear();
114  for (pos=B1T.begin(); pos<B1T.end(); pos++)
115  { delete [] *pos; }
116  B1T.clear();
117  for (pos=S2T.begin(); pos<S2T.end(); pos++)
118  { delete [] *pos; }
119  S2T.clear();
120  for (pos=B2T.begin(); pos<B2T.end(); pos++)
121  { delete [] *pos; }
122  B2T.clear();
123  for (pos=S3T.begin(); pos<S3T.end(); pos++)
124  { delete [] *pos; }
125  S3T.clear();
126  for (pos=B3T.begin(); pos<B3T.end(); pos++)
127  { delete [] *pos; }
128  B3T.clear();
129  for (pos=S4T.begin(); pos<S4T.end(); pos++)
130  { delete [] *pos; }
131  S4T.clear();
132  for (pos=B4T.begin(); pos<B4T.end(); pos++)
133  { delete [] *pos; }
134  B4T.clear();
135 }
136 
137 void
139 {
140  outFile << "G4ChipsPionMinusElasticXS provides the elastic cross\n"
141  << "section for pion- nucleus scattering as a function of incident\n"
142  << "momentum. The cross section is calculated using M. Kossov's\n"
143  << "CHIPS parameterization of cross section data.\n";
144 }
145 
146 
148  const G4Element*,
149  const G4Material*)
150 {
151  return true;
152 }
153 
154 // The main member function giving the collision cross section (P is in IU, CS is in mb)
155 // Make pMom in independent units ! (Now it is MeV)
157  const G4Isotope*,
158  const G4Element*,
159  const G4Material*)
160 {
161  G4double pMom=Pt->GetTotalMomentum();
162  G4int tgN = A - tgZ;
163 
164  return GetChipsCrossSection(pMom, tgZ, tgN, -212);
165 }
166 
168 {
169 
170  G4double pEn=pMom;
171  G4bool fCS = false;
172  onlyCS=fCS;
173 
174  G4bool in=false; // By default the isotope must be found in the AMDB
175  lastP = 0.; // New momentum history (nothing to compare with)
176  lastN = tgN; // The last N of the calculated nucleus
177  lastZ = tgZ; // The last Z of the calculated nucleus
178  lastI = colN.size(); // Size of the Associative Memory DB in the heap
179  if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
180  { // The nucleus with projPDG is found in AMDB
181  if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
182  {
183  lastI=i;
184  lastTH =colTH[i]; // Last THreshold (A-dependent)
185  if(pEn<=lastTH)
186  {
187  return 0.; // Energy is below the Threshold value
188  }
189  lastP =colP [i]; // Last Momentum (A-dependent)
190  lastCS =colCS[i]; // Last CrossSect (A-dependent)
191  // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
192  if(lastP == pMom) // Do not recalculate
193  {
194  CalculateCrossSection(fCS,-1,i,-211,lastZ,lastN,pMom); // Update param's only
195  return lastCS*millibarn; // Use theLastCS
196  }
197  in = true; // This is the case when the isotop is found in DB
198  // Momentum pMom is in IU ! @@ Units
199  lastCS=CalculateCrossSection(fCS,-1,i,-211,lastZ,lastN,pMom); // read & update
200  if(lastCS<=0. && pEn>lastTH) // Correct the threshold
201  {
202  lastTH=pEn;
203  }
204  break; // Go out of the LOOP with found lastI
205  }
206  } // End of attampt to find the nucleus in DB
207  if(!in) // This nucleus has not been calculated previously
208  {
210  lastCS=CalculateCrossSection(fCS,0,lastI,-211,lastZ,lastN,pMom);//calculate&create
211  if(lastCS<=0.)
212  {
213  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
214  if(pEn>lastTH)
215  {
216  lastTH=pEn;
217  }
218  }
219  colN.push_back(tgN);
220  colZ.push_back(tgZ);
221  colP.push_back(pMom);
222  colTH.push_back(lastTH);
223  colCS.push_back(lastCS);
224  return lastCS*millibarn;
225  } // End of creation of the new set of parameters
226  else
227  {
228  colP[lastI]=pMom;
229  colCS[lastI]=lastCS;
230  }
231  return lastCS*millibarn;
232 }
233 
234 // Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
235 // F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
236 G4double G4ChipsPionMinusElasticXS::CalculateCrossSection(G4bool CS, G4int F,G4int I,
237  G4int PDG, G4int tgZ, G4int tgN, G4double pIU)
238 {
239  G4double pMom=pIU/GeV; // All calculations are in GeV
240  onlyCS=CS; // Flag to calculate only CS (not Si/Bi)
241  lastLP=G4Log(pMom); // Make a logarithm of the momentum for calculation
242  if(F) // This isotope was found in AMDB =>RETRIEVE/UPDATE
243  {
244  if(F<0) // the AMDB must be loded
245  {
246  lastPIN = PIN[I]; // Max log(P) initialised for this table set
247  lastPAR = PAR[I]; // Pointer to the parameter set
248  lastCST = CST[I]; // Pointer to the total sross-section table
249  lastSST = SST[I]; // Pointer to the first squared slope
250  lastS1T = S1T[I]; // Pointer to the first mantissa
251  lastB1T = B1T[I]; // Pointer to the first slope
252  lastS2T = S2T[I]; // Pointer to the second mantissa
253  lastB2T = B2T[I]; // Pointer to the second slope
254  lastS3T = S3T[I]; // Pointer to the third mantissa
255  lastB3T = B3T[I]; // Pointer to the rhird slope
256  lastS4T = S4T[I]; // Pointer to the 4-th mantissa
257  lastB4T = B4T[I]; // Pointer to the 4-th slope
258  }
259  if(lastLP>lastPIN && lastLP<lPMax)
260  {
261  lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
262  PIN[I]=lastPIN; // Remember the new P-Limit of the tables
263  }
264  }
265  else // This isotope wasn't initialized => CREATE
266  {
267  lastPAR = new G4double[nPoints]; // Allocate memory for parameters of CS function
268  lastPAR[nLast]=0; // Initialization for VALGRIND
269  lastCST = new G4double[nPoints]; // Allocate memory for Tabulated CS function
270  lastSST = new G4double[nPoints]; // Allocate memory for Tabulated first sqaredSlope
271  lastS1T = new G4double[nPoints]; // Allocate memory for Tabulated first mantissa
272  lastB1T = new G4double[nPoints]; // Allocate memory for Tabulated first slope
273  lastS2T = new G4double[nPoints]; // Allocate memory for Tabulated second mantissa
274  lastB2T = new G4double[nPoints]; // Allocate memory for Tabulated second slope
275  lastS3T = new G4double[nPoints]; // Allocate memory for Tabulated third mantissa
276  lastB3T = new G4double[nPoints]; // Allocate memory for Tabulated third slope
277  lastS4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th mantissa
278  lastB4T = new G4double[nPoints]; // Allocate memory for Tabulated 4-th slope
279  lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
280  PIN.push_back(lastPIN); // Fill parameters of CS function to AMDB
281  PAR.push_back(lastPAR); // Fill parameters of CS function to AMDB
282  CST.push_back(lastCST); // Fill Tabulated CS function to AMDB
283  SST.push_back(lastSST); // Fill Tabulated first sq.slope to AMDB
284  S1T.push_back(lastS1T); // Fill Tabulated first mantissa to AMDB
285  B1T.push_back(lastB1T); // Fill Tabulated first slope to AMDB
286  S2T.push_back(lastS2T); // Fill Tabulated second mantissa to AMDB
287  B2T.push_back(lastB2T); // Fill Tabulated second slope to AMDB
288  S3T.push_back(lastS3T); // Fill Tabulated third mantissa to AMDB
289  B3T.push_back(lastB3T); // Fill Tabulated third slope to AMDB
290  S4T.push_back(lastS4T); // Fill Tabulated 4-th mantissa to AMDB
291  B4T.push_back(lastB4T); // Fill Tabulated 4-th slope to AMDB
292  } // End of creation/update of the new set of parameters and tables
293  // =----------= NOW Update (if necessary) and Calculate the Cross Section =-----------=
294  if(lastLP>lastPIN && lastLP<lPMax)
295  {
296  lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
297  }
298  if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
299  if(lastLP>lPMin && lastLP<=lastPIN) // Linear fit is made using precalculated tables
300  {
301  if(lastLP==lastPIN)
302  {
303  G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
304  G4int blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
305  if(blast<0 || blast>=nLast) G4cout<<"G4QEleastCS::CCS:b="<<blast<<","<<nLast<<G4endl;
306  lastSIG = lastCST[blast];
307  if(!onlyCS) // Skip the differential cross-section parameters
308  {
309  theSS = lastSST[blast];
310  theS1 = lastS1T[blast];
311  theB1 = lastB1T[blast];
312  theS2 = lastS2T[blast];
313  theB2 = lastB2T[blast];
314  theS3 = lastS3T[blast];
315  theB3 = lastB3T[blast];
316  theS4 = lastS4T[blast];
317  theB4 = lastB4T[blast];
318  }
319  }
320  else
321  {
322  G4double shift=(lastLP-lPMin)/dlnP; // a shift from the beginning of the table
323  G4int blast=static_cast<int>(shift); // the lower bin number
324  if(blast<0) blast=0;
325  if(blast>=nLast) blast=nLast-1; // low edge of the last bin
326  shift-=blast; // step inside the unit bin
327  G4int lastL=blast+1; // the upper bin number
328  G4double SIGL=lastCST[blast]; // the basic value of the cross-section
329  lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
330  if(!onlyCS) // Skip the differential cross-section parameters
331  {
332  G4double SSTL=lastSST[blast]; // the low bin of the first squared slope
333  theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
334  G4double S1TL=lastS1T[blast]; // the low bin of the first mantissa
335  theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
336  G4double B1TL=lastB1T[blast]; // the low bin of the first slope
337  theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
338  G4double S2TL=lastS2T[blast]; // the low bin of the second mantissa
339  theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
340  G4double B2TL=lastB2T[blast]; // the low bin of the second slope
341  theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
342  G4double S3TL=lastS3T[blast]; // the low bin of the third mantissa
343  theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
344  G4double B3TL=lastB3T[blast]; // the low bin of the third slope
345  theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
346  G4double S4TL=lastS4T[blast]; // the low bin of the 4-th mantissa
347  theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
348  G4double B4TL=lastB4T[blast]; // the low bin of the 4-th slope
349  theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
350  }
351  }
352  }
353  else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
354  if(lastSIG<0.) lastSIG = 0.; // @@ a Warning print can be added
355  return lastSIG;
356 }
357 
358 // It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
359 G4double G4ChipsPionMinusElasticXS::GetPTables(G4double LP, G4double ILP, G4int PDG,
360  G4int tgZ, G4int tgN)
361 {
362  // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
363  static const G4double pwd=2727;
364  const G4int n_pimpel=38; // #of parameters for pp-elastic (<nPoints=128)
365  // -0- -1- -2- -3- -4- -5- -6- -7- -8- -9--10-11-12-
366  G4double pimp_el[n_pimpel]={1.27,1.53,.0676,3.5,.36,.04,.017,.0025,.0557,2.4,7.,.7,.6,
367  .05,5.,74.,3.,3.4,.2,.17,.001,8.,.055,3.64,5.e-5,4000.,1500.,
368  .46,1.2e6,3.5e6,5.e-5,1.e10,8.5e8,1.e10,1.1,3.4e6,6.8e6,0.};
369  // -13-14--15--16--17-18--19--20--21- -22- -23- -24- -25- -26-
370  // -27--28- -29- -30- -31- -32- -33- -34- -35- -36- -37-
371  if(PDG ==-211)
372  {
373  // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
374  //p2=p*p;p3=p2*p;sp=sqrt(p);p2s=p2*sp;lp=log(p);dl1=lp-(3.=par(3));p4=p2*p2; p=|3-mom|
375  //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
376  // par(0) par(7) par(1) par(2) par(4) par(5) par(6)
377  //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
378  // par(8) par(9) par(10) par(11) par(12)par(13) par(14)
379  // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
380  // par(15) par(16) par(17) par(18) par(19) par(20) par(21) par(22) par(23)
381  // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
382  // par(24) par(25) par(26) par(27) par(28) par(29) par(30) par(31)
383  //
384  if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
385  {
386  if ( tgZ == 1 && tgN == 0 )
387  {
388  for (G4int ip=0; ip<n_pimpel; ip++) lastPAR[ip]=pimp_el[ip]; // PiMinus+P
389  }
390  else
391  {
392  G4double a=tgZ+tgN;
393  G4double sa=std::sqrt(a);
394  G4double ssa=std::sqrt(sa);
395  G4double asa=a*sa;
396  G4double a2=a*a;
397  G4double a3=a2*a;
398  G4double a4=a3*a;
399  G4double a5=a4*a;
400  G4double a6=a4*a2;
401  G4double a7=a6*a;
402  G4double a8=a7*a;
403  G4double a9=a8*a;
404  G4double a10=a5*a5;
405  G4double a12=a6*a6;
406  G4double a14=a7*a7;
407  G4double a16=a8*a8;
408  G4double a17=a16*a;
409  //G4double a20=a16*a4;
410  G4double a32=a16*a16;
411  // Reaction cross-section parameters (pel=peh_fit.f)
412  lastPAR[0]=(.95*sa+2.E5/a16)/(1.+17/a); // p1
413  lastPAR[1]=a/(1./4.4+1./a); // p2
414  lastPAR[2]=.22/G4Pow::GetInstance()->powA(a,.33); // p3
415  lastPAR[3]=.5*a/(1.+3./a+1800./a8); // p4
416  lastPAR[4]=3.E-4*G4Pow::GetInstance()->powA(a,.32)/(1.+14./a2); // p5
417  lastPAR[5]=0.; // p6 not used
418  lastPAR[6]=(.55+.001*a2)/(1.+4.E-4*a2); // p7
419  lastPAR[7]=(.0002/asa+4.E-9*a)/(1.+9./a4); // p8
420  lastPAR[8]=0.; // p9 not used
421  // @@ the differential cross-section is parameterized separately for A>6 & A<7
422  if(a<6.5)
423  {
424  G4double a28=a16*a12;
425  // The main pre-exponent (pel_sg)
426  lastPAR[ 9]=4000*a; // p1
427  lastPAR[10]=1.2e7*a8+380*a17; // p2
428  lastPAR[11]=.7/(1.+4.e-12*a16); // p3
429  lastPAR[12]=2.5/a8/(a4+1.e-16*a32); // p4
430  lastPAR[13]=.28*a; // p5
431  lastPAR[14]=1.2*a2+2.3; // p6
432  lastPAR[15]=3.8/a; // p7
433  // The main slope (pel_sl)
434  lastPAR[16]=.01/(1.+.0024*a5); // p1
435  lastPAR[17]=.2*a; // p2
436  lastPAR[18]=9.e-7/(1.+.035*a5); // p3
437  lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a); // p4
438  // The main quadratic (pel_sh)
439  lastPAR[20]=2.25*a3; // p1
440  lastPAR[21]=18.; // p2
441  lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7); // p3
442  lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a); // p4
443  // The 1st max pre-exponent (pel_qq)
444  lastPAR[24]=1.e5/(a8+2.5e12/a16); // p1
445  lastPAR[25]=8.e7/(a12+1.e-27*a28*a28); // p2
446  lastPAR[26]=.0006*a3; // p3
447  // The 1st max slope (pel_qs)
448  lastPAR[27]=10.+4.e-8*a12*a; // p1
449  lastPAR[28]=.114; // p2
450  lastPAR[29]=.003; // p3
451  lastPAR[30]=2.e-23; // p4
452  // The effective pre-exponent (pel_ss)
453  lastPAR[31]=1./(1.+.0001*a8); // p1
454  lastPAR[32]=1.5e-4/(1.+5.e-6*a12); // p2
455  lastPAR[33]=.03; // p3
456  // The effective slope (pel_sb)
457  lastPAR[34]=a/2; // p1
458  lastPAR[35]=2.e-7*a4; // p2
459  lastPAR[36]=4.; // p3
460  lastPAR[37]=64./a3; // p4
461  // The gloria pre-exponent (pel_us)
462  lastPAR[38]=1.e8*G4Exp(.32*asa); // p1
463  lastPAR[39]=20.*G4Exp(.45*asa); // p2
464  lastPAR[40]=7.e3+2.4e6/a5; // p3
465  lastPAR[41]=2.5e5*G4Exp(.085*a3); // p4
466  lastPAR[42]=2.5*a; // p5
467  // The gloria slope (pel_ub)
468  lastPAR[43]=920.+.03*a8*a3; // p1
469  lastPAR[44]=93.+.0023*a12; // p2
470  }
471  else
472  {
473  G4double p1a10=2.2e-28*a10;
474  G4double r4a16=6.e14/a16;
475  G4double s4a16=r4a16*r4a16;
476  // a24
477  // a36
478  // The main pre-exponent (peh_sg)
479  lastPAR[ 9]=4.5*G4Pow::GetInstance()->powA(a,1.15); // p1
480  lastPAR[10]=.06*G4Pow::GetInstance()->powA(a,.6); // p2
481  lastPAR[11]=.6*a/(1.+2.e15/a16); // p3
482  lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32); // p4
483  lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5); // p5
484  lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12); // p6
485  // The main slope (peh_sl)
486  lastPAR[15]=400./a12+2.e-22*a9; // p1
487  lastPAR[16]=1.e-32*a12/(1.+5.e22/a14); // p2
488  lastPAR[17]=1000./a2+9.5*sa*ssa; // p3
489  lastPAR[18]=4.e-6*a*asa+1.e11/a16; // p4
490  lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16); // p5
491  lastPAR[20]=9.+100./a; // p6
492  // The main quadratic (peh_sh)
493  lastPAR[21]=.002*a3+3.e7/a6; // p1
494  lastPAR[22]=7.e-15*a4*asa; // p2
495  lastPAR[23]=9000./a4; // p3
496  // The 1st max pre-exponent (peh_qq)
497  lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4); // p1
498  lastPAR[25]=1.e-5*a2+2.e14/a16; // p2
499  lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12); // p3
500  lastPAR[27]=.016*asa/(1.+5.e16/a16); // p4
501  // The 1st max slope (peh_qs)
502  lastPAR[28]=.002*a4/(1.+7.e7/G4Pow::GetInstance()->powA(a-6.83,14)); // p1
503  lastPAR[29]=2.e6/a6+7.2/G4Pow::GetInstance()->powA(a,.11); // p2
504  lastPAR[30]=11.*a3/(1.+7.e23/a16/a8); // p3
505  lastPAR[31]=100./asa; // p4
506  // The 2nd max pre-exponent (peh_ss)
507  lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4); // p1
508  lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8); // p2
509  lastPAR[34]=1.3+3.e5/a4; // p3
510  lastPAR[35]=500./(a2+50.)+3; // p4
511  lastPAR[36]=1.e-9/a+s4a16*s4a16; // p5
512  // The 2nd max slope (peh_sb)
513  lastPAR[37]=.4*asa+3.e-9*a6; // p1
514  lastPAR[38]=.0005*a5; // p2
515  lastPAR[39]=.002*a5; // p3
516  lastPAR[40]=10.; // p4
517  // The effective pre-exponent (peh_us)
518  lastPAR[41]=.05+.005*a; // p1
519  lastPAR[42]=7.e-8/sa; // p2
520  lastPAR[43]=.8*sa; // p3
521  lastPAR[44]=.02*sa; // p4
522  lastPAR[45]=1.e8/a3; // p5
523  lastPAR[46]=3.e32/(a32+1.e32); // p6
524  // The effective slope (peh_ub)
525  lastPAR[47]=24.; // p1
526  lastPAR[48]=20./sa; // p2
527  lastPAR[49]=7.e3*a/(sa+1.); // p3
528  lastPAR[50]=900.*sa/(1.+500./a3); // p4
529  }
530  // Parameter for lowEnergyNeutrons
531  lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
532  }
533  lastPAR[nLast]=pwd;
534  // and initialize the zero element of the table
535  G4double lp=lPMin; // ln(momentum)
536  G4bool memCS=onlyCS; // ??
537  onlyCS=false;
538  lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
539  onlyCS=memCS;
540  lastSST[0]=theSS;
541  lastS1T[0]=theS1;
542  lastB1T[0]=theB1;
543  lastS2T[0]=theS2;
544  lastB2T[0]=theB2;
545  lastS3T[0]=theS3;
546  lastB3T[0]=theB3;
547  lastS4T[0]=theS4;
548  lastB4T[0]=theB4;
549  }
550  if(LP>ILP)
551  {
552  G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
553  if(ini<0) ini=0;
554  if(ini<nPoints)
555  {
556  G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
557  if(fin>=nPoints) fin=nLast; // Limit of the tabular initialization
558  if(fin>=ini)
559  {
560  G4double lp=0.;
561  for(G4int ip=ini; ip<=fin; ip++) // Calculate tabular CS,S1,B1,S2,B2,S3,B3
562  {
563  lp=lPMin+ip*dlnP; // ln(momentum)
564  G4bool memCS=onlyCS;
565  onlyCS=false;
566  lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
567  onlyCS=memCS;
568  lastSST[ip]=theSS;
569  lastS1T[ip]=theS1;
570  lastB1T[ip]=theB1;
571  lastS2T[ip]=theS2;
572  lastB2T[ip]=theB2;
573  lastS3T[ip]=theS3;
574  lastB3T[ip]=theB3;
575  lastS4T[ip]=theS4;
576  lastB4T[ip]=theB4;
577  }
578  return lp;
579  }
580  else G4cout<<"*Warning*G4ChipsPionMinusElasticXS::GetPTables: PDG="<<PDG
581  <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP
582  <<" > ILP="<<ILP<<" nothing is done!"<<G4endl;
583  }
584  else G4cout<<"*Warning*G4ChipsPionMinusElasticXS::GetPTables: PDG="<<PDG
585  <<", Z="<<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
586  <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
587  }
588  }
589  else
590  {
591  // G4cout<<"*Error*G4ChipsPionMinusElasticXS::GetPTables: PDG="<<PDG<<", Z="<<tgZ
592  // <<", N="<<tgN<<", while it is defined only for PDG=-211"<<G4endl;
593  // throw G4QException("G4ChipsPionMinusElasticXS::GetPTables:onlyPipA implemented");
595  ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
596  << ", while it is defined only for PDG=-211 (pi-)" << G4endl;
597  G4Exception("G4ChipsPionMinusElasticXS::GetPTables()", "HAD_CHPS_0000",
598  FatalException, ed);
599  }
600  return ILP;
601 }
602 
603 // Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
605 {
606  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
607  static const G4double third=1./3.;
608  static const G4double fifth=1./5.;
609  static const G4double sevth=1./7.;
610  if(PDG!=-211)G4cout<<"Warning*G4ChipsPionMinusElasticXS::GetExT:PDG="<<PDG<<G4endl;
611  if(onlyCS)G4cout<<"Warning*G4ChipsPionMinusElasticXS::GetExchanT:onlyCS=1"<<G4endl;
612  if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
613  G4double q2=0.;
614  if(tgZ==1 && tgN==0) // ===> p+p=p+p
615  {
616  G4double E1=lastTM*theB1;
617  G4double R1=(1.-G4Exp(-E1));
618  G4double E2=lastTM*theB2;
619  G4double R2=(1.-G4Exp(-E2*E2*E2));
620  G4double E3=lastTM*theB3;
621  G4double R3=(1.-G4Exp(-E3));
622  G4double I1=R1*theS1/theB1;
623  G4double I2=R2*theS2;
624  G4double I3=R3*theS3;
625  G4double I12=I1+I2;
626  G4double rand=(I12+I3)*G4UniformRand();
627  if (rand<I1 )
628  {
629  G4double ran=R1*G4UniformRand();
630  if(ran>1.) ran=1.;
631  q2=-G4Log(1.-ran)/theB1;
632  }
633  else if(rand<I12)
634  {
635  G4double ran=R2*G4UniformRand();
636  if(ran>1.) ran=1.;
637  q2=-G4Log(1.-ran);
638  if(q2<0.) q2=0.;
639  q2=G4Pow::GetInstance()->powA(q2,third)/theB2;
640  }
641  else
642  {
643  G4double ran=R3*G4UniformRand();
644  if(ran>1.) ran=1.;
645  q2=-G4Log(1.-ran)/theB3;
646  }
647  }
648  else
649  {
650  G4double a=tgZ+tgN;
651  G4double E1=lastTM*(theB1+lastTM*theSS);
652  G4double R1=(1.-G4Exp(-E1));
653  G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
654  G4double tm2=lastTM*lastTM;
655  G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
656  if(a>6.5)E2*=tm2; // for heavy nuclei
657  G4double R2=(1.-G4Exp(-E2));
658  G4double E3=lastTM*theB3;
659  if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
660  G4double R3=(1.-G4Exp(-E3));
661  G4double E4=lastTM*theB4;
662  G4double R4=(1.-G4Exp(-E4));
663  G4double I1=R1*theS1;
664  G4double I2=R2*theS2;
665  G4double I3=R3*theS3;
666  G4double I4=R4*theS4;
667  G4double I12=I1+I2;
668  G4double I13=I12+I3;
669  G4double rand=(I13+I4)*G4UniformRand();
670  if(rand<I1)
671  {
672  G4double ran=R1*G4UniformRand();
673  if(ran>1.) ran=1.;
674  q2=-G4Log(1.-ran)/theB1;
675  if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
676  }
677  else if(rand<I12)
678  {
679  G4double ran=R2*G4UniformRand();
680  if(ran>1.) ran=1.;
681  q2=-G4Log(1.-ran)/theB2;
682  if(q2<0.) q2=0.;
683  if(a<6.5) q2=G4Pow::GetInstance()->powA(q2,third);
684  else q2=G4Pow::GetInstance()->powA(q2,fifth);
685  }
686  else if(rand<I13)
687  {
688  G4double ran=R3*G4UniformRand();
689  if(ran>1.) ran=1.;
690  q2=-G4Log(1.-ran)/theB3;
691  if(q2<0.) q2=0.;
692  if(a>6.5) q2=G4Pow::GetInstance()->powA(q2,sevth);
693  }
694  else
695  {
696  G4double ran=R4*G4UniformRand();
697  if(ran>1.) ran=1.;
698  q2=-G4Log(1.-ran)/theB4;
699  if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
700  }
701  }
702  if(q2<0.) q2=0.;
703  if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
704  if(q2>lastTM)
705  {
706  q2=lastTM;
707  }
708  return q2*GeVSQ;
709 }
710 
711 // Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
712 G4double G4ChipsPionMinusElasticXS::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
713 {
714  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
715  if(onlyCS)G4cout<<"Warning*G4ChipsPionMinusElasticXS::GetSlope:onlCS=true"<<G4endl;
716  if(lastLP<-4.3) return 0.; // S-wave for p<14 MeV/c (kinE<.1MeV)
717  if(PDG !=-211)
718  {
719  // G4cout<<"*Error*G4ChipsPionMinusElasticXS::GetSlope: PDG="<<PDG<<", Z="<<tgZ
720  // <<", N="<<tgN<<", while it is defined only for PDG=-211"<<G4endl;
721  // throw G4QException("G4ChipsPionMinusElasticXS::GetSlope: pipA are implemented");
723  ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
724  << ", while it is defined only for PDG=-211" << G4endl;
725  G4Exception("G4ChipsPionMinusElasticXS::GetSlope()", "HAD_CHPS_0000",
726  FatalException, ed);
727  }
728  if(theB1<0.) theB1=0.;
729  if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QElasticCrossSect::Getslope:"<<theB1<<G4endl;
730  return theB1/GeVSQ;
731 }
732 
733 // Returns half max(Q2=-t) in independent units (MeV^2)
734 G4double G4ChipsPionMinusElasticXS::GetHMaxT()
735 {
736  static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
737  return lastTM*HGeVSQ;
738 }
739 
740 // lastLP is used, so calculating tables, one need to remember and then recover lastLP
741 G4double G4ChipsPionMinusElasticXS::GetTabValues(G4double lp, G4int PDG, G4int tgZ,
742  G4int tgN)
743 {
744  if(PDG!=-211)G4cout<<"*Warn*G4ChipsPionMinusElasticXS::GetTabV: PDG="<<PDG<<G4endl;
745  if(tgZ<0 || tgZ>92)
746  {
747  G4cout<<"*Warning*G4QPionPlusElCS::GetTabValue:(1-92) No isotopes for Z="<<tgZ<<G4endl;
748  return 0.;
749  }
750  G4int iZ=tgZ-1; // Z index
751  if(iZ<0)
752  {
753  iZ=0; // conversion of the neutron target to the proton target
754  tgZ=1;
755  tgN=0;
756  }
757  G4double p=G4Exp(lp); // momentum
758  G4double sp=std::sqrt(p); // sqrt(p)
759  G4double p2=p*p;
760  G4double p3=p2*p;
761  G4double p4=p3*p;
762  if ( tgZ == 1 && tgN == 0 ) // PiMin+P
763  {
764  G4double dl2=lp-lastPAR[14];
765  theSS=lastPAR[37];
766  theS1=(lastPAR[15]+lastPAR[16]*dl2*dl2)/(1.+lastPAR[17]/p4/p)+
767  (lastPAR[18]/p2+lastPAR[19]*p)/(p4+lastPAR[20]*sp);
768  theB1=lastPAR[21]*G4Pow::GetInstance()->powA(p,lastPAR[22])/(1.+lastPAR[23]/p3);
769  theS2=lastPAR[24]+lastPAR[25]/(p4+lastPAR[26]*p);
770  theB2=lastPAR[27]+lastPAR[28]/(p4+lastPAR[29]/sp);
771  theS3=lastPAR[30]+lastPAR[31]/(p4*p4+lastPAR[32]*p2+lastPAR[33]);
772  theB3=lastPAR[34]+lastPAR[35]/(p4+lastPAR[36]);
773  theS4=0.;
774  theB4=0.;
775  // Returns the total elastic pim-p cross-section (to avoid spoiling lastSIG)
776  G4double lr=lp+lastPAR[0]; // lr
777  G4double ld=lp-lastPAR[14];
778  G4double dl3=lp+lastPAR[4]; // lm
779  G4double dl4=lp-lastPAR[6]; // lh
780 //G4cout<<"lastPAR[13] "<<lastPAR[13]<<" lastPAR[6] "<<lastPAR[6]<<" lastPAR[7] "<<lastPAR[7]<<G4endl;
781  return lastPAR[1]/(lr*lr+lastPAR[2])+
782  (lastPAR[8]*ld*ld+lastPAR[9]+lastPAR[10]/sp)/(1.+lastPAR[11]/p4)+
783  lastPAR[12]/(dl3*dl3+lastPAR[5])+lastPAR[13]/(dl4*dl4+lastPAR[7]);
784  }
785  else
786  {
787  G4double p5=p4*p;
788  G4double p6=p5*p;
789  G4double p8=p6*p2;
790  G4double p10=p8*p2;
791  G4double p12=p10*p2;
792  G4double p16=p8*p8;
793  //G4double p24=p16*p8;
794  G4double dl=lp-5.;
795  G4double a=tgZ+tgN;
796  G4double pah=G4Pow::GetInstance()->powA(p,a/2);
797  G4double pa=pah*pah;
798  G4double pa2=pa*pa;
799  if(a<6.5)
800  {
801  theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
802  (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
803  theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
804  theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
805  theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
806  theB2=lastPAR[27]*G4Pow::GetInstance()->powA(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
807  theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
808  theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
809  theS4=p2*(pah*lastPAR[38]*G4Exp(-pah*lastPAR[39])+
810  lastPAR[40]/(1.+lastPAR[41]*G4Pow::GetInstance()->powA(p,lastPAR[42])));
811  theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
812  }
813  else
814  {
815  theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
816  lastPAR[13]/(p5+lastPAR[14]/p16);
817  theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/G4Pow::GetInstance()->powA(p,lastPAR[20]))+
818  lastPAR[17]/(1.+lastPAR[18]/p4);
819  theSS=lastPAR[21]/(p4/G4Pow::GetInstance()->powA(p,lastPAR[23])+lastPAR[22]/p4);
820  theS2=lastPAR[24]/p4/(G4Pow::GetInstance()->powA(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
821  theB2=lastPAR[28]/G4Pow::GetInstance()->powA(p,lastPAR[29])+lastPAR[30]/G4Pow::GetInstance()->powA(p,lastPAR[31]);
822  theS3=lastPAR[32]/G4Pow::GetInstance()->powA(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
823  lastPAR[33]/(1.+lastPAR[34]/p6);
824  theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
825  theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
826  (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
827  theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
828  }
829  // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
830  // p1 p2 p3
831  return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p8)+
832  lastPAR[3]/(p4+lastPAR[4]/p3)+lastPAR[6]/(p4+lastPAR[7]/p4);
833  // p4 p5 p7 p8
834  }
835  return 0.;
836 } // End of GetTableValues
837 
838 // Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
839 G4double G4ChipsPionMinusElasticXS::GetQ2max(G4int PDG, G4int tgZ, G4int tgN,
840  G4double pP)
841 {
842  static const G4double mPi= G4PionMinus::PionMinus()->GetPDGMass()*.001; // MeV to GeV
843  static const G4double mPi2= mPi*mPi;
844 
845  G4double pP2=pP*pP; // squared momentum of the projectile
846  if(tgZ || tgN>-1) // ---> pipA
847  {
848  G4double mt=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(tgZ,tgZ+tgN,0)->GetPDGMass()*.001; // Target mass in GeV
849 
850  G4double dmt=mt+mt;
851  G4double mds=dmt*std::sqrt(pP2+mPi2)+mPi2+mt*mt; // Mondelstam mds
852  return dmt*dmt*pP2/mds;
853  }
854  else
855  {
857  ed << "PDG = " << PDG << ",Z = " << tgZ << ",N = " << tgN
858  << ", while it is defined only for p projectiles & Z_target>0" << G4endl;
859  G4Exception("G4ChipsPionMinusElasticXS::GetQ2max()", "HAD_CHPS_0000",
860  FatalException, ed);
861  return 0;
862  }
863 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
virtual G4double GetIsoCrossSection(const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
virtual void CrossSectionDescription(std::ostream &) const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const char * p
Definition: xmltok.h:285
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
G4double GetExchangeT(G4int tZ, G4int tN, G4int pPDG)
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
G4IonTable * GetIonTable() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
virtual G4bool IsIsoApplicable(const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
bool G4bool
Definition: G4Types.hh:79
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
#define G4_DECLARE_XS_FACTORY(cross_section)
G4double GetPDGMass() const
static G4ParticleTable * GetParticleTable()
static constexpr double gigaelectronvolt
Definition: G4SIunits.hh:209
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
static const G4double pos
static constexpr double millibarn
Definition: G4SIunits.hh:106