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