Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GGNuclNuclCrossSection Class Reference

#include <G4GGNuclNuclCrossSection.hh>

+ Inheritance diagram for G4GGNuclNuclCrossSection:

Public Member Functions

 G4GGNuclNuclCrossSection ()
 
virtual ~G4GGNuclNuclCrossSection ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *)
 
G4double GetZandACrossSection (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetCoulombBarier (const G4DynamicParticle *, G4double Z, G4double A, G4double pR, G4double tR)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
G4double GetRatioSD (const G4DynamicParticle *, G4double At, G4double Zt)
 
G4double GetRatioQE (const G4DynamicParticle *, G4double At, G4double Zt)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, const G4Element *)
 
G4double GetHadronNucleonXsc (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double GetHadronNucleonXscPDG (G4ParticleDefinition *, G4double sMand, G4ParticleDefinition *)
 
G4double GetHadronNucleonXscNS (G4ParticleDefinition *, G4double pTkin, G4ParticleDefinition *)
 
G4double GetHNinelasticXscVU (const G4DynamicParticle *, G4int At, G4int Zt)
 
G4double CalculateEcmValue (const G4double, const G4double, const G4double)
 
G4double CalcMandelstamS (const G4double, const G4double, const G4double)
 
G4double GetElasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetInelasticGlauberGribov (const G4DynamicParticle *, G4int Z, G4int A)
 
G4double GetTotalGlauberGribovXsc ()
 
G4double GetElasticGlauberGribovXsc ()
 
G4double GetInelasticGlauberGribovXsc ()
 
G4double GetProductionGlauberGribovXsc ()
 
G4double GetDiffractionGlauberGribovXsc ()
 
G4double GetRadiusConst ()
 
G4double GetNucleusRadius (const G4DynamicParticle *, const G4Element *)
 
G4double GetNucleusRadius (G4double Zt, G4double At)
 
G4double GetNucleusRadiusGG (G4double At)
 
G4double GetNucleusRadiusDE (G4double Z, G4double A)
 
G4double GetNucleusRadiusRMS (G4double Z, G4double A)
 
void SetEnergyLowerLimit (G4double E)
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4bool IsIsoApplicable (const G4DynamicParticle *, G4int Z, G4int A, const G4Element *elm=0, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int Z, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 50 of file G4GGNuclNuclCrossSection.hh.

Constructor & Destructor Documentation

◆ G4GGNuclNuclCrossSection()

G4GGNuclNuclCrossSection::G4GGNuclNuclCrossSection ( )

Definition at line 46 of file G4GGNuclNuclCrossSection.cc.

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}
static const char * Default_Name()
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Proton * Proton()
Definition: G4Proton.cc:93

◆ ~G4GGNuclNuclCrossSection()

G4GGNuclNuclCrossSection::~G4GGNuclNuclCrossSection ( )
virtual

Definition at line 59 of file G4GGNuclNuclCrossSection.cc.

60{}

Member Function Documentation

◆ BuildPhysicsTable()

virtual void G4GGNuclNuclCrossSection::BuildPhysicsTable ( const G4ParticleDefinition )
inlinevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 75 of file G4GGNuclNuclCrossSection.hh.

76 {}

◆ CalcMandelstamS()

G4double G4GGNuclNuclCrossSection::CalcMandelstamS ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 766 of file G4GGNuclNuclCrossSection.cc.

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}
double G4double
Definition: G4Types.hh:64

Referenced by GetHadronNucleonXsc(), and GetHadronNucleonXscNS().

◆ CalculateEcmValue()

G4double G4GGNuclNuclCrossSection::CalculateEcmValue ( const G4double  mp,
const G4double  mt,
const G4double  Plab 
)

Definition at line 749 of file G4GGNuclNuclCrossSection.cc.

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}

◆ CrossSectionDescription()

void G4GGNuclNuclCrossSection::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 63 of file G4GGNuclNuclCrossSection.cc.

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}

◆ Default_Name()

static const char * G4GGNuclNuclCrossSection::Default_Name ( )
inlinestatic

Definition at line 57 of file G4GGNuclNuclCrossSection.hh.

57{return "Glauber-Gribov nucleus nucleus";}

Referenced by G4IonPhysics::ConstructProcess().

◆ DumpPhysicsTable()

virtual void G4GGNuclNuclCrossSection::DumpPhysicsTable ( const G4ParticleDefinition )
inlinevirtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 79 of file G4GGNuclNuclCrossSection.hh.

80 {G4cout << "G4NuclNuclCrossSection: uses Glauber-Gribov formula"<<G4endl;}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout

◆ GetCoulombBarier()

G4double G4GGNuclNuclCrossSection::GetCoulombBarier ( const G4DynamicParticle aParticle,
G4double  Z,
G4double  A,
G4double  pR,
G4double  tR 
)

Definition at line 202 of file G4GGNuclNuclCrossSection.cc.

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}
int G4int
Definition: G4Types.hh:66
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4double GetIonMass(G4int Z, G4int A, G4int L=0) const
!! Only ground states are supported now
Definition: G4IonTable.cc:774
G4double GetPDGCharge() const
static G4ParticleTable * GetParticleTable()
G4IonTable * GetIonTable()

Referenced by GetZandACrossSection().

◆ GetDiffractionGlauberGribovXsc()

G4double G4GGNuclNuclCrossSection::GetDiffractionGlauberGribovXsc ( )
inline

Definition at line 105 of file G4GGNuclNuclCrossSection.hh.

105{ return fDiffractionXsc; };

◆ GetElasticGlauberGribov()

G4double G4GGNuclNuclCrossSection::GetElasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 136 of file G4GGNuclNuclCrossSection.hh.

138{
139 GetZandACrossSection(dp, Z, A);
140 return fElasticXsc;
141}
G4double GetZandACrossSection(const G4DynamicParticle *, G4int Z, G4int A)

◆ GetElasticGlauberGribovXsc()

G4double G4GGNuclNuclCrossSection::GetElasticGlauberGribovXsc ( )
inline

Definition at line 102 of file G4GGNuclNuclCrossSection.hh.

102{ return fElasticXsc; };

◆ GetElementCrossSection()

G4double G4GGNuclNuclCrossSection::GetElementCrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 90 of file G4GGNuclNuclCrossSection.cc.

93{
94 G4int A = G4lrint(G4NistManager::Instance()->GetAtomicMassAmu(Z));
95 return GetZandACrossSection(aParticle, Z, A);
96}
static G4NistManager * Instance()
int G4lrint(double ad)
Definition: templates.hh:163

◆ GetHadronNucleonXsc() [1/2]

G4double G4GGNuclNuclCrossSection::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
const G4Element anElement 
)

Definition at line 328 of file G4GGNuclNuclCrossSection.cc.

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}
G4double GetZ() const
Definition: G4Element.hh:131
G4double GetN() const
Definition: G4Element.hh:134
G4double GetHadronNucleonXsc(const G4DynamicParticle *, const G4Element *)

Referenced by GetHadronNucleonXsc().

◆ GetHadronNucleonXsc() [2/2]

G4double G4GGNuclNuclCrossSection::GetHadronNucleonXsc ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 347 of file G4GGNuclNuclCrossSection.cc.

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}
double mag() const
G4double GetMass() const
G4ThreeVector GetMomentum() const
G4double CalcMandelstamS(const G4double, const G4double, const G4double)

◆ GetHadronNucleonXscNS()

G4double G4GGNuclNuclCrossSection::GetHadronNucleonXscNS ( G4ParticleDefinition pParticle,
G4double  pTkin,
G4ParticleDefinition tParticle 
)

Definition at line 443 of file G4GGNuclNuclCrossSection.cc.

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}
G4double GetHadronNucleonXscPDG(G4ParticleDefinition *, G4double sMand, G4ParticleDefinition *)

Referenced by GetRatioQE(), and GetRatioSD().

◆ GetHadronNucleonXscPDG()

G4double G4GGNuclNuclCrossSection::GetHadronNucleonXscPDG ( G4ParticleDefinition pParticle,
G4double  sMand,
G4ParticleDefinition tParticle 
)

Definition at line 386 of file G4GGNuclNuclCrossSection.cc.

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}
@ neutron
bool G4bool
Definition: G4Types.hh:67

Referenced by GetHadronNucleonXscNS().

◆ GetHNinelasticXscVU()

G4double G4GGNuclNuclCrossSection::GetHNinelasticXscVU ( const G4DynamicParticle aParticle,
G4int  At,
G4int  Zt 
)

Definition at line 540 of file G4GGNuclNuclCrossSection.cc.

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}
G4double GetTotalEnergy() const

◆ GetInelasticGlauberGribov()

G4double G4GGNuclNuclCrossSection::GetInelasticGlauberGribov ( const G4DynamicParticle dp,
G4int  Z,
G4int  A 
)
inline

Definition at line 146 of file G4GGNuclNuclCrossSection.hh.

148{
149 GetZandACrossSection(dp, Z, A);
150 return fInelasticXsc;
151}

◆ GetInelasticGlauberGribovXsc()

G4double G4GGNuclNuclCrossSection::GetInelasticGlauberGribovXsc ( )
inline

Definition at line 103 of file G4GGNuclNuclCrossSection.hh.

103{ return fInelasticXsc; };

◆ GetNucleusRadius() [1/2]

G4double G4GGNuclNuclCrossSection::GetNucleusRadius ( const G4DynamicParticle ,
const G4Element anElement 
)

Definition at line 598 of file G4GGNuclNuclCrossSection.cc.

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}
const G4double oneThird

Referenced by GetRatioQE(), GetRatioSD(), and GetZandACrossSection().

◆ GetNucleusRadius() [2/2]

G4double G4GGNuclNuclCrossSection::GetNucleusRadius ( G4double  Zt,
G4double  At 
)

Definition at line 640 of file G4GGNuclNuclCrossSection.cc.

641{
642 G4double R;
643 R = GetNucleusRadiusDE(Zt,At);
644 // R = GetNucleusRadiusRMS(Zt,At);
645
646 return R;
647}
G4double GetNucleusRadiusDE(G4double Z, G4double A)

◆ GetNucleusRadiusDE()

G4double G4GGNuclNuclCrossSection::GetNucleusRadiusDE ( G4double  Z,
G4double  A 
)

Definition at line 680 of file G4GGNuclNuclCrossSection.cc.

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}

Referenced by GetNucleusRadius().

◆ GetNucleusRadiusGG()

G4double G4GGNuclNuclCrossSection::GetNucleusRadiusGG ( G4double  At)

Definition at line 652 of file G4GGNuclNuclCrossSection.cc.

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}

◆ GetNucleusRadiusRMS()

G4double G4GGNuclNuclCrossSection::GetNucleusRadiusRMS ( G4double  Z,
G4double  A 
)

Definition at line 728 of file G4GGNuclNuclCrossSection.cc.

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}

◆ GetProductionGlauberGribovXsc()

G4double G4GGNuclNuclCrossSection::GetProductionGlauberGribovXsc ( )
inline

Definition at line 104 of file G4GGNuclNuclCrossSection.hh.

104{ return fProductionXsc; };

◆ GetRadiusConst()

G4double G4GGNuclNuclCrossSection::GetRadiusConst ( )
inline

Definition at line 106 of file G4GGNuclNuclCrossSection.hh.

106{ return fRadiusConst; };

◆ GetRatioQE()

G4double G4GGNuclNuclCrossSection::GetRatioQE ( const G4DynamicParticle aParticle,
G4double  At,
G4double  Zt 
)

Definition at line 282 of file G4GGNuclNuclCrossSection.cc.

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}
G4double GetHadronNucleonXscNS(G4ParticleDefinition *, G4double pTkin, G4ParticleDefinition *)
G4double GetNucleusRadius(const G4DynamicParticle *, const G4Element *)
const G4double pi

◆ GetRatioSD()

G4double G4GGNuclNuclCrossSection::GetRatioSD ( const G4DynamicParticle aParticle,
G4double  At,
G4double  Zt 
)

Definition at line 242 of file G4GGNuclNuclCrossSection.cc.

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}

◆ GetTotalGlauberGribovXsc()

G4double G4GGNuclNuclCrossSection::GetTotalGlauberGribovXsc ( )
inline

Definition at line 101 of file G4GGNuclNuclCrossSection.hh.

101{ return fTotalXsc; };

◆ GetZandACrossSection()

G4double G4GGNuclNuclCrossSection::GetZandACrossSection ( const G4DynamicParticle aParticle,
G4int  Z,
G4int  A 
)

Definition at line 107 of file G4GGNuclNuclCrossSection.cc.

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}
G4ThreeVector G4ParticleMomentum
G4double GetCoulombBarier(const G4DynamicParticle *, G4double Z, G4double A, G4double pR, G4double tR)
G4double GetHadronNucleonXscNS(const G4DynamicParticle *, const G4ParticleDefinition *)
G4double GetInelasticHadronNucleonXsc()

Referenced by GetElasticGlauberGribov(), GetElementCrossSection(), and GetInelasticGlauberGribov().

◆ IsElementApplicable()

G4bool G4GGNuclNuclCrossSection::IsElementApplicable ( const G4DynamicParticle ,
G4int  Z,
const G4Material  
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 72 of file G4GGNuclNuclCrossSection.cc.

74{
75 G4bool applicable = true;
76// G4double kineticEnergy = aDP->GetKineticEnergy();
77
78// if (kineticEnergy >= fLowerLimit) applicable = true;
79 return applicable;
80}

◆ SetEnergyLowerLimit()

void G4GGNuclNuclCrossSection::SetEnergyLowerLimit ( G4double  E)
inline

Definition at line 115 of file G4GGNuclNuclCrossSection.hh.

115{fLowerLimit=E;};

The documentation for this class was generated from the following files: