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