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