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