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