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