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