Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4GGNuclNuclCrossSection.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 // 24.11.08 V. Grichine - first implementation
27 // 25.10.12 W.Pokorski - following Vladimir's advice, I removed Z>1 condition
28 //
29 
31 
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 #include "G4ParticleTable.hh"
35 #include "G4IonTable.hh"
36 #include "G4ParticleDefinition.hh"
37 #include "G4HadTmpUtil.hh"
38 #include "G4HadronNucleonXsc.hh"
39 
40 // factory
41 #include "G4CrossSectionFactory.hh"
42 //
44 
45 
47  : G4VCrossSectionDataSet(Default_Name()),
48  fUpperLimit(100000*GeV), fLowerLimit(0.1*MeV),
49  fRadiusConst(1.08*fermi), // 1.1, 1.3 ?
50  fTotalXsc(0.0), fElasticXsc(0.0), fInelasticXsc(0.0), fProductionXsc(0.0),
51  fDiffractionXsc(0.0), fHadronNucleonXsc(0.0)
52 {
53  theProton = G4Proton::Proton();
54  theNeutron = G4Neutron::Neutron();
55  hnXsc = new G4HadronNucleonXsc();
56 }
57 
58 
60 {
61  delete hnXsc;
62 }
63 
64 void
66 {
67  outFile << "G4GGNuclNuclCrossSection calculates total, inelastic and\n"
68  << "elastic cross sections for nucleus-nucleus collisions using\n"
69  << "the Glauber model with Gribov corrections. It is valid for\n"
70  << "all incident energies above 100 keV./n";
71 }
72 
73 G4bool
75  G4int, const G4Material*)
76 {
77  G4bool applicable = true;
78 // G4double kineticEnergy = aDP->GetKineticEnergy();
79 
80 // if (kineticEnergy >= fLowerLimit) applicable = true;
81  return applicable;
82 }
83 
85 //
86 // Calculates total and inelastic Xsc, derives elastic as total - inelastic
87 // accordong to Glauber model with Gribov correction calculated in the dipole
88 // approximation on light cone. Gaussian density helps to calculate rest
89 // integrals of the model. [1] B.Z. Kopeliovich, nucl-th/0306044
90 
91 
94  const G4Material*)
95 {
96  G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
97  return GetZandACrossSection(aParticle, Z, A);
98 }
99 
101 //
102 // Calculates total and inelastic Xsc, derives elastic as total - inelastic
103 // accordong to Glauber model with Gribov correction calculated in the dipole
104 // approximation on light cone. Gaussian density of point-like nucleons helps
105 // to calculate rest integrals of the model. [1] B.Z. Kopeliovich,
106 // nucl-th/0306044 + simplification above
107 
108 
111  G4int tZ, G4int tA)
112 {
113  G4double xsection;
114  G4double sigma;
115  G4double cofInelastic = 2.4;
116  G4double cofTotal = 2.0;
117  G4double nucleusSquare;
118  G4double cB;
119  G4double ratio;
120 
121  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
122  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
123 
124  G4double pTkin = aParticle->GetKineticEnergy();
125  pTkin /= pA;
126 
127  G4double pN = pA - pZ;
128  if( pN < 0. ) pN = 0.;
129 
130  G4double tN = tA - tZ;
131  if( tN < 0. ) tN = 0.;
132 
133  G4double tR = GetNucleusRadius( G4double(tZ),G4double(tA) );
134  G4double pR = GetNucleusRadius(pZ,pA);
135 
136  cB = GetCoulombBarier(aParticle, G4double(tZ), G4double(tA), pR, tR);
137 
138  if ( cB > 0. )
139  {
140  G4DynamicParticle* dProton = new G4DynamicParticle(theProton,
141  G4ParticleMomentum(1.,0.,0.),
142  pTkin);
143 
144  G4DynamicParticle* dNeutron = new G4DynamicParticle(theNeutron,
145  G4ParticleMomentum(1.,0.,0.),
146  pTkin);
147 
148  sigma = (pZ*tZ+pN*tN)*hnXsc->GetHadronNucleonXscNS(dProton, theProton);
149 
150  G4double ppInXsc = hnXsc->GetInelasticHadronNucleonXsc();
151 
152  sigma += (pZ*tN+pN*tZ)*hnXsc->GetHadronNucleonXscNS(dNeutron, theProton);
153 
154  G4double npInXsc = hnXsc->GetInelasticHadronNucleonXsc();
155 
156  delete dProton;
157  delete dNeutron;
158 
159  // G4cout<<"ppInXsc = "<<ppInXsc/millibarn<<"; npInXsc = "<<npInXsc/millibarn<<G4endl;
160  // G4cout<<"npTotXsc = "<<hnXsc->GetTotalHadronNucleonXsc()/millibarn<<"; npElXsc = "
161  // <<hnXsc->GetElasticHadronNucleonXsc()/millibarn<<G4endl;
162 
163  nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
164 
165  ratio = sigma/nucleusSquare;
166  xsection = nucleusSquare*std::log( 1. + ratio );
167  fTotalXsc = xsection;
168  fTotalXsc *= cB;
169 
170  fInelasticXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
171 
172  fInelasticXsc *= cB;
173  fElasticXsc = fTotalXsc - fInelasticXsc;
174 
175  // if (fElasticXsc < DBL_MIN) fElasticXsc = DBL_MIN;
176  /*
177  G4double difratio = ratio/(1.+ratio);
178 
179  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
180  */
181  // production to be checked !!! edit MK xsc
182 
183  //sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscMK(theProton, pTkin, theProton) +
184  // (pZ*tN+pN*tZ)*GetHadronNucleonXscMK(theProton, pTkin, theNeutron);
185 
186  sigma = (pZ*tZ+pN*tN)*ppInXsc + (pZ*tN+pN*tZ)*npInXsc;
187 
188  ratio = sigma/nucleusSquare;
189  fProductionXsc = nucleusSquare*std::log( 1. + cofInelastic*ratio )/cofInelastic;
190 
191  if (fElasticXsc < 0.) fElasticXsc = 0.;
192  }
193  else
194  {
195  fInelasticXsc = 0.;
196  fTotalXsc = 0.;
197  fElasticXsc = 0.;
198  fProductionXsc = 0.;
199  }
200 
201  return fInelasticXsc; // xsection;
202 }
203 
205 //
206 //
207 
210  G4double pR, G4double tR)
211 {
212  G4double ratio;
213  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
214 
215  G4double pTkin = aParticle->GetKineticEnergy();
216  // G4double pPlab = aParticle->GetTotalMomentum();
217  G4double pM = aParticle->GetDefinition()->GetPDGMass();
218  // G4double tM = tZ*proton_mass_c2 + (tA-tZ)*neutron_mass_c2; // ~ 1% accuracy
220  G4double pElab = pTkin + pM;
221  G4double totEcm = std::sqrt(pM*pM + tM*tM + 2.*pElab*tM);
222  // G4double pPcm = pPlab*tM/totEcm;
223  // G4double pTcm = std::sqrt(pM*pM + pPcm*pPcm) - pM;
224  G4double totTcm = totEcm - pM -tM;
225 
227  bC /= pR + tR;
228  bC /= 2.; // 4., 2. parametrisation cof ??? vmg
229 
230  // G4cout<<"pTkin = "<<pTkin/GeV<<"; pPlab = "
231  // <<pPlab/GeV<<"; bC = "<<bC/GeV<<"; pTcm = "<<pTcm/GeV<<G4endl;
232 
233  if( totTcm <= bC ) ratio = 0.;
234  else ratio = 1. - bC/totTcm;
235 
236  // if(ratio < DBL_MIN) ratio = DBL_MIN;
237  if( ratio < 0.) ratio = 0.;
238 
239  // G4cout <<"ratio = "<<ratio<<G4endl;
240  return ratio;
241 }
242 
243 
245 //
246 // Return single-diffraction/inelastic cross-section ratio
247 
250 {
251  G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
252 
253  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
254  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
255 
256  G4double pTkin = aParticle->GetKineticEnergy();
257  pTkin /= pA;
258 
259  G4double pN = pA - pZ;
260  if( pN < 0. ) pN = 0.;
261 
262  G4double tN = tA - tZ;
263  if( tN < 0. ) tN = 0.;
264 
265  G4double tR = GetNucleusRadius(tZ,tA);
266  G4double pR = GetNucleusRadius(pZ,pA);
267 
268  sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
269  (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
270 
271  nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
272  ratio = sigma/nucleusSquare;
273  fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
274  G4double difratio = ratio/(1.+ratio);
275 
276  fDiffractionXsc = 0.5*nucleusSquare*( difratio - std::log( 1. + difratio ) );
277 
278  if (fInelasticXsc > 0.) ratio = fDiffractionXsc/fInelasticXsc;
279  else ratio = 0.;
280 
281  return ratio;
282 }
283 
285 //
286 // Return quasi-elastic/inelastic cross-section ratio
287 
290 {
291  G4double sigma, cofInelastic = 2.4, cofTotal = 2.0, nucleusSquare, ratio;
292 
293  G4double pZ = aParticle->GetDefinition()->GetPDGCharge();
294  G4double pA = aParticle->GetDefinition()->GetBaryonNumber();
295 
296  G4double pTkin = aParticle->GetKineticEnergy();
297  pTkin /= pA;
298 
299  G4double pN = pA - pZ;
300  if( pN < 0. ) pN = 0.;
301 
302  G4double tN = tA - tZ;
303  if( tN < 0. ) tN = 0.;
304 
305  G4double tR = GetNucleusRadius(tZ,tA);
306  G4double pR = GetNucleusRadius(pZ,pA);
307 
308  sigma = (pZ*tZ+pN*tN)*GetHadronNucleonXscNS(theProton, pTkin, theProton) +
309  (pZ*tN+pN*tZ)*GetHadronNucleonXscNS(theProton, pTkin, theNeutron);
310 
311  nucleusSquare = cofTotal*pi*( pR*pR + tR*tR ); // basically 2piRR
312  ratio = sigma/nucleusSquare;
313  fInelasticXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
314 
315  // sigma = GetHNinelasticXsc(aParticle, tA, tZ);
316  ratio = sigma/nucleusSquare;
317  fProductionXsc = nucleusSquare*std::log(1. + cofInelastic*ratio)/cofInelastic;
318 
319  if (fInelasticXsc > fProductionXsc) ratio = (fInelasticXsc-fProductionXsc)/fInelasticXsc;
320  else ratio = 0.;
321  if ( ratio < 0. ) ratio = 0.;
322 
323  return ratio;
324 }
325 
327 //
328 // Returns hadron-nucleon Xsc according to differnt parametrisations:
329 // [2] E. Levin, hep-ph/9710546
330 // [3] U. Dersch, et al, hep-ex/9910052
331 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
332 
333 G4double
335  const G4Element* anElement)
336 {
337  G4int At = G4lrint(anElement->GetN()); // number of nucleons
338  G4int Zt = G4lrint(anElement->GetZ()); // number of protons
339  return GetHadronNucleonXsc(aParticle, At, Zt);
340 }
341 
342 
343 
344 
346 //
347 // Returns hadron-nucleon Xsc according to differnt parametrisations:
348 // [2] E. Levin, hep-ph/9710546
349 // [3] U. Dersch, et al, hep-ex/9910052
350 // [4] M.J. Longo, et al, Phys.Rev.Lett. 33 (1974) 725
351 
352 G4double
354  G4int At, G4int Zt)
355 {
356  G4double xsection = 0.;
357 
359  GetIonTable()->GetIonMass(Zt, At);
360  targ_mass = 0.939*GeV; // ~mean neutron and proton ???
361 
362  G4double proj_mass = aParticle->GetMass();
363  G4double proj_momentum = aParticle->GetMomentum().mag();
364  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
365 
366  sMand /= GeV*GeV; // in GeV for parametrisation
367  proj_momentum /= GeV;
368  const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
369 
370  if(pParticle == theNeutron) // as proton ???
371  {
372  xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
373  }
374  else if(pParticle == theProton)
375  {
376  xsection = G4double(At)*(21.70*std::pow(sMand,0.0808) + 56.08*std::pow(sMand,-0.4525));
377  }
378 
379  xsection *= millibarn;
380  return xsection;
381 }
382 
383 
384 
386 //
387 // Returns hadron-nucleon Xsc according to PDG parametrisation (2005):
388 // http://pdg.lbl.gov/2006/reviews/hadronicrpp.pdf
389 // At = number of nucleons, Zt = number of protons
390 
391 G4double
393  G4double sMand,
394  G4ParticleDefinition* tParticle)
395 {
396  G4double xsection = 0.;
397  // G4bool pORn = (tParticle == theProton || nucleon == theNeutron );
398  G4bool proton = (tParticle == theProton);
399  G4bool neutron = (tParticle == theNeutron);
400 
401  // General PDG fit constants
402 
403  G4double s0 = 5.38*5.38; // in Gev^2
404  G4double eta1 = 0.458;
405  G4double eta2 = 0.458;
406  G4double B = 0.308;
407 
408  // const G4ParticleDefinition* pParticle = aParticle->GetDefinition();
409 
410  if(pParticle == theNeutron) // proton-neutron fit
411  {
412  if ( proton )
413  {
414  xsection = ( 35.80 + B*std::pow(std::log(sMand/s0),2.)
415  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
416  }
417  if ( neutron )
418  {
419  xsection = (35.45 + B*std::pow(std::log(sMand/s0),2.)
420  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2)); // pp for nn
421  }
422  }
423  else if(pParticle == theProton)
424  {
425  if ( proton )
426  {
427  xsection = (35.45 + B*std::pow(std::log(sMand/s0),2.)
428  + 42.53*std::pow(sMand,-eta1) - 33.34*std::pow(sMand,-eta2));
429 
430  }
431  if ( neutron )
432  {
433  xsection = (35.80 + B*std::pow(std::log(sMand/s0),2.)
434  + 40.15*std::pow(sMand,-eta1) - 30.*std::pow(sMand,-eta2));
435  }
436  }
437  xsection *= millibarn; // parametrised in mb
438  return xsection;
439 }
440 
441 
443 //
444 // Returns nucleon-nucleon cross-section based on N. Starkov parametrisation of
445 // data from mainly http://wwwppds.ihep.su:8001/c5-6A.html database
446 // projectile nucleon is pParticle with pTkin shooting target nucleon tParticle
447 
448 G4double
450  G4double pTkin,
451  G4ParticleDefinition* tParticle)
452 {
453  G4double xsection(0);
454  // G4double Delta; DHW 19 May 2011: variable set but not used
455  G4double A0, B0;
456  G4double hpXscv(0);
457  G4double hnXscv(0);
458 
459  G4double targ_mass = tParticle->GetPDGMass();
460  G4double proj_mass = pParticle->GetPDGMass();
461 
462  G4double proj_energy = proj_mass + pTkin;
463  G4double proj_momentum = std::sqrt(pTkin*(pTkin+2*proj_mass));
464 
465  G4double sMand = CalcMandelstamS ( proj_mass , targ_mass , proj_momentum );
466 
467  sMand /= GeV*GeV; // in GeV for parametrisation
468  proj_momentum /= GeV;
469  proj_energy /= GeV;
470  proj_mass /= GeV;
471 
472  // General PDG fit constants
473 
474  // G4double s0 = 5.38*5.38; // in Gev^2
475  // G4double eta1 = 0.458;
476  // G4double eta2 = 0.458;
477  // G4double B = 0.308;
478 
479  if( proj_momentum >= 373.)
480  {
481  return GetHadronNucleonXscPDG(pParticle,sMand,tParticle);
482  }
483  else if( proj_momentum >= 10. ) // high energy: pp = nn = np
484  // if( proj_momentum >= 2.)
485  {
486  // Delta = 1.; // DHW 19 May 2011: variable set but not used
487  // if (proj_energy < 40.) Delta = 0.916+0.0021*proj_energy;
488 
489  if (proj_momentum >= 10.) {
490  B0 = 7.5;
491  A0 = 100. - B0*std::log(3.0e7);
492 
493  xsection = A0 + B0*std::log(proj_energy) - 11
494  + 103*std::pow(2*0.93827*proj_energy + proj_mass*proj_mass+
495  0.93827*0.93827,-0.165); // mb
496  }
497  }
498  else // low energy pp = nn != np
499  {
500  if(pParticle == tParticle) // pp or nn // nn to be pp
501  {
502  if( proj_momentum < 0.73 )
503  {
504  hnXscv = 23 + 50*( std::pow( std::log(0.73/proj_momentum), 3.5 ) );
505  }
506  else if( proj_momentum < 1.05 )
507  {
508  hnXscv = 23 + 40*(std::log(proj_momentum/0.73))*
509  (std::log(proj_momentum/0.73));
510  }
511  else // if( proj_momentum < 10. )
512  {
513  hnXscv = 39.0 +
514  75*(proj_momentum - 1.2)/(std::pow(proj_momentum,3.0) + 0.15);
515  }
516  xsection = hnXscv;
517  }
518  else // pn to be np
519  {
520  if( proj_momentum < 0.8 )
521  {
522  hpXscv = 33+30*std::pow(std::log(proj_momentum/1.3),4.0);
523  }
524  else if( proj_momentum < 1.4 )
525  {
526  hpXscv = 33+30*std::pow(std::log(proj_momentum/0.95),2.0);
527  }
528  else // if( proj_momentum < 10. )
529  {
530  hpXscv = 33.3+
531  20.8*(std::pow(proj_momentum,2.0)-1.35)/
532  (std::pow(proj_momentum,2.50)+0.95);
533  }
534  xsection = hpXscv;
535  }
536  }
537  xsection *= millibarn; // parametrised in mb
538  return xsection;
539 }
540 
542 //
543 // Returns hadron-nucleon inelastic cross-section based on FTF-parametrisation
544 
545 G4double
547  G4int At, G4int Zt)
548 {
549  G4int PDGcode = aParticle->GetDefinition()->GetPDGEncoding();
550  G4int absPDGcode = std::abs(PDGcode);
551  G4double Elab = aParticle->GetTotalEnergy();
552  // (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
553  G4double Plab = aParticle->GetMomentum().mag();
554  // std::sqrt(Elab * Elab - 0.88);
555 
556  Elab /= GeV;
557  Plab /= GeV;
558 
559  G4double LogPlab = std::log( Plab );
560  G4double sqrLogPlab = LogPlab * LogPlab;
561 
562  //G4cout<<"Plab = "<<Plab<<G4endl;
563 
564  G4double NumberOfTargetProtons = Zt;
565  G4double NumberOfTargetNucleons = At;
566  G4double NumberOfTargetNeutrons = NumberOfTargetNucleons - NumberOfTargetProtons;
567 
568  if(NumberOfTargetNeutrons < 0.) NumberOfTargetNeutrons = 0.;
569 
570  G4double Xtotal = 0., Xelastic = 0., Xinelastic =0.;
571 
572  if( absPDGcode > 1000 ) //------Projectile is baryon --------
573  {
574  G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) +
575  0.522*sqrLogPlab - 4.51*LogPlab;
576 
577  G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) +
578  0.513*sqrLogPlab - 4.27*LogPlab;
579 
580  G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) +
581  0.169*sqrLogPlab - 1.85*LogPlab;
582 
583  G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) +
584  0.169*sqrLogPlab - 1.85*LogPlab;
585 
586  Xtotal = ( NumberOfTargetProtons * XtotPP +
587  NumberOfTargetNeutrons * XtotPN );
588 
589  Xelastic = ( NumberOfTargetProtons * XelPP +
590  NumberOfTargetNeutrons * XelPN );
591  }
592 
593  Xinelastic = Xtotal - Xelastic;
594  if(Xinelastic < 0.) Xinelastic = 0.;
595 
596  return Xinelastic*= millibarn;
597 }
598 
600 //
601 //
602 
603 G4double
605  const G4Element* anElement)
606 {
607  G4double At = anElement->GetN();
608  G4double oneThird = 1.0/3.0;
609  G4double cubicrAt = std::pow (At, oneThird);
610 
611  G4double R; // = fRadiusConst*cubicrAt;
612  R = fRadiusConst*cubicrAt;
613 
614  G4double meanA = 21.;
615  G4double tauA1 = 40.;
616  G4double tauA2 = 10.;
617  G4double tauA3 = 5.;
618 
619  G4double a1 = 0.85;
620  G4double b1 = 1. - a1;
621 
622  G4double b2 = 0.3;
623  G4double b3 = 4.;
624 
625  if (At > 20.) // 20.
626  {
627  R *= ( a1 + b1*std::exp( -(At - meanA)/tauA1) );
628  }
629  else if (At > 3.5)
630  {
631  R *= ( 1.0 + b2*( 1. - std::exp( (At - meanA)/tauA2) ) );
632  }
633  else
634  {
635  R *= ( 1.0 + b3*( 1. - std::exp( (At - meanA)/tauA3) ) );
636  }
637 
638  return R;
639 }
640 
642 //
643 //
644 
645 G4double
647 {
648  G4double R;
649  R = GetNucleusRadiusDE(Zt,At);
650  // R = GetNucleusRadiusRMS(Zt,At);
651 
652  return R;
653 }
654 
656 
657 G4double
659 {
660  G4double oneThird = 1.0/3.0;
661  G4double cubicrAt = std::pow (At, oneThird);
662 
663  G4double R; // = fRadiusConst*cubicrAt;
664  R = fRadiusConst*cubicrAt;
665 
666  G4double meanA = 20.;
667  G4double tauA = 20.;
668 
669  if ( At > 20.) // 20.
670  {
671  R *= ( 0.8 + 0.2*std::exp( -(At - meanA)/tauA) );
672  }
673  else
674  {
675  R *= ( 1.0 + 0.1*( 1. - std::exp( (At - meanA)/tauA) ) );
676  }
677 
678  return R;
679 }
680 
682 //
683 //
684 
685 G4double
687 {
688  // algorithm from diffuse-elastic
689 
690  G4double R, r0, a11, a12, a13, a2, a3;
691 
692  a11 = 1.26; // 1.08, 1.16
693  a12 = 1.; // 1.08, 1.16
694  a13 = 1.12; // 1.08, 1.16
695  a2 = 1.1;
696  a3 = 1.;
697 
698  // Special rms radii for light nucleii
699 
700  if (A < 50.)
701  {
702  if (std::abs(A-1.) < 0.5) return 0.89*fermi; // p
703  else if(std::abs(A-2.) < 0.5) return 2.13*fermi; // d
704  else if(std::abs(Z-1.) < 0.5 && std::abs(A-3.) < 0.5) return 1.80*fermi; // t
705 
706  else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96*fermi; // He3
707  else if(std::abs(Z-2.) < 0.5 && std::abs(A-4.) < 0.5) return 1.68*fermi; // He4
708 
709  else if(std::abs(Z-3.) < 0.5) return 2.40*fermi; // Li7
710  else if(std::abs(Z-4.) < 0.5) return 2.51*fermi; // Be9
711 
712  else if( 10. < A && A <= 16. ) r0 = a11*( 1 - std::pow(A, -2./3.) )*fermi; // 1.08*fermi;
713  else if( 15. < A && A <= 20. ) r0 = a12*( 1 - std::pow(A, -2./3.) )*fermi;
714  else if( 20. < A && A <= 30. ) r0 = a13*( 1 - std::pow(A, -2./3.) )*fermi;
715  else r0 = a2*fermi;
716 
717  R = r0*std::pow( A, 1./3. );
718  }
719  else
720  {
721  r0 = a3*fermi;
722 
723  R = r0*std::pow(A, 0.27);
724  }
725  return R;
726 }
727 
728 
730 //
731 // RMS radii from e-A scattering data
732 
733 G4double
735 {
736 
737  if (std::abs(A-1.) < 0.5) return 0.89*fermi; // p
738  else if(std::abs(A-2.) < 0.5) return 2.13*fermi; // d
739  else if(std::abs(Z-1.) < 0.5 && std::abs(A-3.) < 0.5) return 1.80*fermi; // t
740 
741  else if(std::abs(Z-2.) < 0.5 && std::abs(A-3.) < 0.5) return 1.96*fermi; // He3
742  else if(std::abs(Z-2.) < 0.5 && std::abs(A-4.) < 0.5) return 1.68*fermi; // He4
743 
744  else if(std::abs(Z-3.) < 0.5) return 2.40*fermi; // Li7
745  else if(std::abs(Z-4.) < 0.5) return 2.51*fermi; // Be9
746 
747  else return 1.24*std::pow(A, 0.28 )*fermi; // A > 9
748 }
749 
750 
752 //
753 //
754 
756  const G4double mt,
757  const G4double Plab)
758 {
759  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
760  G4double Ecm = std::sqrt ( mp * mp + mt * mt + 2 * Elab * mt );
761  // G4double Pcm = Plab * mt / Ecm;
762  // G4double KEcm = std::sqrt ( Pcm * Pcm + mp * mp ) - mp;
763 
764  return Ecm ; // KEcm;
765 }
766 
767 
769 //
770 //
771 
773  const G4double mt,
774  const G4double Plab)
775 {
776  G4double Elab = std::sqrt ( mp * mp + Plab * Plab );
777  G4double sMand = mp*mp + mt*mt + 2*Elab*mt ;
778 
779  return sMand;
780 }
781 
782 //
783 //