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