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