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