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

#include <G4FTFParameters.hh>

Public Member Functions

 G4FTFParameters ()
 
 ~G4FTFParameters ()
 
void InitForInteraction (const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
 
void SethNcmsEnergy (const G4double s)
 
void SetTotalCrossSection (const G4double Xtotal)
 
void SetElastisCrossSection (const G4double Xelastic)
 
void SetInelasticCrossSection (const G4double Xinelastic)
 
void SetProbabilityOfElasticScatt (const G4double Xtotal, const G4double Xelastic)
 
void SetProbabilityOfElasticScatt (const G4double aValue)
 
void SetProbabilityOfAnnihilation (const G4double aValue)
 
void SetRadiusOfHNinteractions2 (const G4double Radius2)
 
void SetSlope (const G4double Slope)
 
void SetGamma0 (const G4double Gamma0)
 
G4double GammaElastic (const G4double impactsquare)
 
void SetAvaragePt2ofElasticScattering (const G4double aPt2)
 
void SetParams (const G4int ProcN, const G4double A1, const G4double B1, const G4double A2, const G4double B2, const G4double A3, const G4double Atop, const G4double Ymin)
 
void SetDeltaProbAtQuarkExchange (const G4double aValue)
 
void SetProbOfSameQuarkExchange (const G4double aValue)
 
void SetProjMinDiffMass (const G4double aValue)
 
void SetProjMinNonDiffMass (const G4double aValue)
 
void SetProbLogDistrPrD (const G4double aValue)
 
void SetTarMinDiffMass (const G4double aValue)
 
void SetTarMinNonDiffMass (const G4double aValue)
 
void SetAveragePt2 (const G4double aValue)
 
void SetProbLogDistr (const G4double aValue)
 
void SetPt2Kink (const G4double aValue)
 
void SetQuarkProbabilitiesAtGluonSplitUp (const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
 
void SetMaxNumberOfCollisions (const G4double aValue, const G4double bValue)
 
void SetProbOfInteraction (const G4double aValue)
 
void SetCofNuclearDestructionPr (const G4double aValue)
 
void SetCofNuclearDestruction (const G4double aValue)
 
void SetR2ofNuclearDestruction (const G4double aValue)
 
void SetExcitationEnergyPerWoundedNucleon (const G4double aValue)
 
void SetDofNuclearDestruction (const G4double aValue)
 
void SetPt2ofNuclearDestruction (const G4double aValue)
 
void SetMaxPt2ofNuclearDestruction (const G4double aValue)
 
G4double GetTotalCrossSection ()
 
G4double GetElasticCrossSection ()
 
G4double GetInelasticCrossSection ()
 
G4double GetProbabilityOfInteraction (const G4double impactsquare)
 
G4double GetInelasticProbability (const G4double impactsquare)
 
G4double GetProbabilityOfElasticScatt ()
 
G4double GetSlope ()
 
G4double GetProbabilityOfAnnihilation ()
 
G4double GetAvaragePt2ofElasticScattering ()
 
G4double GetProcProb (const G4int ProcN, const G4double y)
 
G4double GetDeltaProbAtQuarkExchange ()
 
G4double GetProbOfSameQuarkExchange ()
 
G4double GetProjMinDiffMass ()
 
G4double GetProjMinNonDiffMass ()
 
G4double GetProbLogDistrPrD ()
 
G4double GetTarMinDiffMass ()
 
G4double GetTarMinNonDiffMass ()
 
G4double GetAveragePt2 ()
 
G4double GetProbLogDistr ()
 
G4double GetPt2Kink ()
 
std::vector< G4doubleGetQuarkProbabilitiesAtGluonSplitUp ()
 
G4double GetMaxNumberOfCollisions ()
 
G4double GetProbOfInteraction ()
 
G4double GetCofNuclearDestructionPr ()
 
G4double GetCofNuclearDestruction ()
 
G4double GetR2ofNuclearDestruction ()
 
G4double GetExcitationEnergyPerWoundedNucleon ()
 
G4double GetDofNuclearDestruction ()
 
G4double GetPt2ofNuclearDestruction ()
 
G4double GetMaxPt2ofNuclearDestruction ()
 

Public Attributes

G4double FTFhNcmsEnergy
 
G4double FTFXtotal
 
G4double FTFXelastic
 
G4double FTFXinelastic
 
G4double FTFXannihilation
 
G4double ProbabilityOfAnnihilation
 
G4double ProbabilityOfElasticScatt
 
G4double RadiusOfHNinteractions2
 
G4double FTFSlope
 
G4double AvaragePt2ofElasticScattering
 
G4double FTFGamma0
 
G4double ProcParams [5][7]
 
G4double DeltaProbAtQuarkExchange
 
G4double ProbOfSameQuarkExchange
 
G4double ProjMinDiffMass
 
G4double ProjMinNonDiffMass
 
G4double ProbLogDistrPrD
 
G4double TarMinDiffMass
 
G4double TarMinNonDiffMass
 
G4double AveragePt2
 
G4double ProbLogDistr
 
G4double Pt2kink
 
std::vector< G4doubleQuarkProbabilitiesAtGluonSplitUp
 
G4double MaxNumberOfCollisions
 
G4double ProbOfInelInteraction
 
G4double CofNuclearDestructionPr
 
G4double CofNuclearDestruction
 
G4double R2ofNuclearDestruction
 
G4double ExcitationEnergyPerWoundedNucleon
 
G4double DofNuclearDestruction
 
G4double Pt2ofNuclearDestruction
 
G4double MaxPt2ofNuclearDestruction
 

Detailed Description

Definition at line 305 of file G4FTFParameters.hh.

Constructor & Destructor Documentation

◆ G4FTFParameters()

G4FTFParameters::G4FTFParameters ( )

Definition at line 62 of file G4FTFParameters.cc.

63{
64 StringMass = new G4LundStringFragmentation; // for estimation of min. mass of diffr. states
65 Reset();
66 csGGinstance =
68 if (!csGGinstance) {
69 csGGinstance = new G4ComponentGGHadronNucleusXsc();
70 }
71
72 // Set parameters of a string kink
73 SetPt2Kink( 0.0*GeV*GeV ); // To switch off kinky strings (bad results obtained with 6.0*GeV*GeV)
74 G4double Puubar( 1.0/3.0 ), Pddbar( 1.0/3.0 ), Pssbar( 1.0/3.0 ); // SU(3) symmetry
75 //G4double Puubar( 0.41 ), Pddbar( 0.41 ), Pssbar( 0.18 ); // Broken SU(3) symmetry
76 SetQuarkProbabilitiesAtGluonSplitUp( Puubar, Pddbar, Pssbar );
77}
double G4double
Definition: G4Types.hh:83
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
void SetQuarkProbabilitiesAtGluonSplitUp(const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
void SetPt2Kink(const G4double aValue)

◆ ~G4FTFParameters()

G4FTFParameters::~G4FTFParameters ( )

Definition at line 1515 of file G4FTFParameters.cc.

1515 {
1516 if ( StringMass ) delete StringMass;
1517}

Member Function Documentation

◆ GammaElastic()

G4double G4FTFParameters::GammaElastic ( const G4double  impactsquare)
inline

Definition at line 488 of file G4FTFParameters.hh.

488 {
489 return ( FTFGamma0 * G4Exp( -FTFSlope * impactsquare ) );
490}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179

Referenced by GetInelasticProbability().

◆ GetAvaragePt2ofElasticScattering()

G4double G4FTFParameters::GetAvaragePt2ofElasticScattering ( )
inline

Definition at line 691 of file G4FTFParameters.hh.

691 {
693}
G4double AvaragePt2ofElasticScattering

Referenced by G4ElasticHNScattering::ElasticScattering(), and InitForInteraction().

◆ GetAveragePt2()

G4double G4FTFParameters::GetAveragePt2 ( )
inline

Definition at line 721 of file G4FTFParameters.hh.

721 {
722 return AveragePt2;
723}

Referenced by InitForInteraction().

◆ GetCofNuclearDestruction()

G4double G4FTFParameters::GetCofNuclearDestruction ( )
inline

Definition at line 757 of file G4FTFParameters.hh.

757 {
759}
G4double CofNuclearDestruction

Referenced by InitForInteraction().

◆ GetCofNuclearDestructionPr()

G4double G4FTFParameters::GetCofNuclearDestructionPr ( )
inline

Definition at line 753 of file G4FTFParameters.hh.

753 {
755}
G4double CofNuclearDestructionPr

Referenced by InitForInteraction().

◆ GetDeltaProbAtQuarkExchange()

G4double G4FTFParameters::GetDeltaProbAtQuarkExchange ( )
inline

Definition at line 697 of file G4FTFParameters.hh.

697 {
699}
G4double DeltaProbAtQuarkExchange

Referenced by InitForInteraction().

◆ GetDofNuclearDestruction()

G4double G4FTFParameters::GetDofNuclearDestruction ( )
inline

Definition at line 769 of file G4FTFParameters.hh.

769 {
771}
G4double DofNuclearDestruction

Referenced by InitForInteraction().

◆ GetElasticCrossSection()

G4double G4FTFParameters::GetElasticCrossSection ( )
inline

Definition at line 657 of file G4FTFParameters.hh.

657 {
658 return FTFXelastic;
659}

◆ GetExcitationEnergyPerWoundedNucleon()

G4double G4FTFParameters::GetExcitationEnergyPerWoundedNucleon ( )
inline

Definition at line 765 of file G4FTFParameters.hh.

765 {
767}
G4double ExcitationEnergyPerWoundedNucleon

Referenced by InitForInteraction().

◆ GetInelasticCrossSection()

G4double G4FTFParameters::GetInelasticCrossSection ( )
inline

Definition at line 661 of file G4FTFParameters.hh.

661 {
662 return FTFXinelastic;
663}

◆ GetInelasticProbability()

G4double G4FTFParameters::GetInelasticProbability ( const G4double  impactsquare)
inline

Definition at line 681 of file G4FTFParameters.hh.

681 {
682 G4double Gamma = GammaElastic( impactsquare );
683 return 2*Gamma - Gamma*Gamma;
684}
G4double GammaElastic(const G4double impactsquare)

◆ GetMaxNumberOfCollisions()

G4double G4FTFParameters::GetMaxNumberOfCollisions ( )
inline

Definition at line 745 of file G4FTFParameters.hh.

745 {
747}
G4double MaxNumberOfCollisions

◆ GetMaxPt2ofNuclearDestruction()

G4double G4FTFParameters::GetMaxPt2ofNuclearDestruction ( )
inline

Definition at line 777 of file G4FTFParameters.hh.

777 {
779}
G4double MaxPt2ofNuclearDestruction

◆ GetProbabilityOfAnnihilation()

G4double G4FTFParameters::GetProbabilityOfAnnihilation ( )
inline

Definition at line 686 of file G4FTFParameters.hh.

686 {
688}
G4double ProbabilityOfAnnihilation

◆ GetProbabilityOfElasticScatt()

G4double G4FTFParameters::GetProbabilityOfElasticScatt ( )
inline

Definition at line 677 of file G4FTFParameters.hh.

677 {
679}
G4double ProbabilityOfElasticScatt

◆ GetProbabilityOfInteraction()

G4double G4FTFParameters::GetProbabilityOfInteraction ( const G4double  impactsquare)
inline

Definition at line 669 of file G4FTFParameters.hh.

669 {
670 if ( RadiusOfHNinteractions2 > impactsquare ) {
671 return 1.0;
672 } else {
673 return 0.0;
674 }
675}
G4double RadiusOfHNinteractions2

Referenced by G4FTFParticipants::GetList().

◆ GetProbLogDistr()

G4double G4FTFParameters::GetProbLogDistr ( )
inline

Definition at line 729 of file G4FTFParameters.hh.

729 {
730 return ProbLogDistr;
731}

Referenced by InitForInteraction().

◆ GetProbLogDistrPrD()

G4double G4FTFParameters::GetProbLogDistrPrD ( )
inline

Definition at line 725 of file G4FTFParameters.hh.

725 {
726 return ProbLogDistrPrD;
727}
G4double ProbLogDistrPrD

Referenced by InitForInteraction().

◆ GetProbOfInteraction()

G4double G4FTFParameters::GetProbOfInteraction ( )
inline

Definition at line 749 of file G4FTFParameters.hh.

749 {
751}
G4double ProbOfInelInteraction

◆ GetProbOfSameQuarkExchange()

G4double G4FTFParameters::GetProbOfSameQuarkExchange ( )
inline

Definition at line 701 of file G4FTFParameters.hh.

701 {
703}
G4double ProbOfSameQuarkExchange

Referenced by InitForInteraction().

◆ GetProcProb()

G4double G4FTFParameters::GetProcProb ( const G4int  ProcN,
const G4double  y 
)

Definition at line 848 of file G4FTFParameters.cc.

848 {
849 G4double Prob( 0.0 );
850 if ( y < ProcParams[ProcN][6] ) {
851 Prob = ProcParams[ProcN][5];
852 if (Prob < 0.) Prob=0.;
853 return Prob;
854 }
855 Prob = ProcParams[ProcN][0] * G4Exp( -ProcParams[ProcN][1]*y ) +
856 ProcParams[ProcN][2] * G4Exp( -ProcParams[ProcN][3]*y ) +
857 ProcParams[ProcN][4];
858 if (Prob < 0.) Prob=0.;
859 return Prob;
860}
G4double ProcParams[5][7]

Referenced by G4DiffractiveExcitation::ExciteParticipants().

◆ GetProjMinDiffMass()

G4double G4FTFParameters::GetProjMinDiffMass ( )
inline

◆ GetProjMinNonDiffMass()

G4double G4FTFParameters::GetProjMinNonDiffMass ( )
inline

Definition at line 709 of file G4FTFParameters.hh.

709 {
710 return ProjMinNonDiffMass;
711}
G4double ProjMinNonDiffMass

Referenced by G4DiffractiveExcitation::ExciteParticipants(), and InitForInteraction().

◆ GetPt2Kink()

G4double G4FTFParameters::GetPt2Kink ( )
inline

Definition at line 735 of file G4FTFParameters.hh.

735 {
736 return Pt2kink;
737}

Referenced by G4DiffractiveExcitation::CreateStrings().

◆ GetPt2ofNuclearDestruction()

G4double G4FTFParameters::GetPt2ofNuclearDestruction ( )
inline

Definition at line 773 of file G4FTFParameters.hh.

773 {
775}
G4double Pt2ofNuclearDestruction

Referenced by InitForInteraction().

◆ GetQuarkProbabilitiesAtGluonSplitUp()

std::vector< G4double > G4FTFParameters::GetQuarkProbabilitiesAtGluonSplitUp ( )
inline

Definition at line 739 of file G4FTFParameters.hh.

739 {
741}
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp

Referenced by G4DiffractiveExcitation::CreateStrings().

◆ GetR2ofNuclearDestruction()

G4double G4FTFParameters::GetR2ofNuclearDestruction ( )
inline

Definition at line 761 of file G4FTFParameters.hh.

761 {
763}
G4double R2ofNuclearDestruction

Referenced by InitForInteraction().

◆ GetSlope()

G4double G4FTFParameters::GetSlope ( )
inline

Definition at line 665 of file G4FTFParameters.hh.

665 {
666 return FTFSlope;
667}

Referenced by InitForInteraction().

◆ GetTarMinDiffMass()

G4double G4FTFParameters::GetTarMinDiffMass ( )
inline

◆ GetTarMinNonDiffMass()

G4double G4FTFParameters::GetTarMinNonDiffMass ( )
inline

Definition at line 717 of file G4FTFParameters.hh.

717 {
718 return TarMinNonDiffMass;
719}
G4double TarMinNonDiffMass

Referenced by G4DiffractiveExcitation::ExciteParticipants(), and InitForInteraction().

◆ GetTotalCrossSection()

G4double G4FTFParameters::GetTotalCrossSection ( )
inline

Definition at line 653 of file G4FTFParameters.hh.

653 {
654 return FTFXtotal;
655}

◆ InitForInteraction()

void G4FTFParameters::InitForInteraction ( const G4ParticleDefinition particle,
G4int  theA,
G4int  theZ,
G4double  s 
)

Definition at line 81 of file G4FTFParameters.cc.

83{
84 Reset();
85
86 G4int ProjectilePDGcode = particle->GetPDGEncoding();
87 G4int ProjectileabsPDGcode = std::abs( ProjectilePDGcode );
88 G4double ProjectileMass = particle->GetPDGMass();
89 G4double ProjectileMass2 = ProjectileMass * ProjectileMass;
90
91 G4int ProjectileBaryonNumber( 0 ), AbsProjectileBaryonNumber( 0 ), AbsProjectileCharge( 0 );
92 G4bool ProjectileIsNucleus = false;
93
94 if ( std::abs( particle->GetBaryonNumber() ) > 1 ) { // The projectile is a nucleus
95 ProjectileIsNucleus = true;
96 ProjectileBaryonNumber = particle->GetBaryonNumber();
97 AbsProjectileBaryonNumber = std::abs( ProjectileBaryonNumber );
98 AbsProjectileCharge = std::abs( G4int( particle->GetPDGCharge() ) );
99 if ( ProjectileBaryonNumber > 1 ) {
100 ProjectilePDGcode = 2212; ProjectileabsPDGcode = 2212; // Proton
101 } else {
102 ProjectilePDGcode = -2212; ProjectileabsPDGcode = 2212; // Anti-Proton
103 }
104 ProjectileMass = G4Proton::Proton()->GetPDGMass();
105 ProjectileMass2 = sqr( ProjectileMass );
106 }
107
108 G4double TargetMass = G4Proton::Proton()->GetPDGMass();
109 G4double TargetMass2 = TargetMass * TargetMass;
110
111 G4double Plab = PlabPerParticle;
112 G4double Elab = std::sqrt( Plab*Plab + ProjectileMass2 );
113 G4double KineticEnergy = Elab - ProjectileMass;
114
115 G4double S = ProjectileMass2 + TargetMass2 + 2.0*TargetMass*Elab;
116
117 #ifdef debugFTFparams
118 G4cout << "--------- FTF Parameters --------------" << G4endl << "Proj Plab "
119 << ProjectilePDGcode << " " << Plab << G4endl << "Mass KinE " << ProjectileMass
120 << " " << KineticEnergy << G4endl << " A Z " << theA << " " << theZ << G4endl;
121 #endif
122
123 G4double Ylab, Xtotal( 0.0 ), Xelastic( 0.0 ), Xannihilation( 0.0 );
124 G4int NumberOfTargetNucleons;
125
126 Ylab = 0.5 * G4Log( (Elab + Plab)/(Elab - Plab) );
127
128 G4double ECMSsqr = S/GeV/GeV;
129 G4double SqrtS = std::sqrt( S )/GeV;
130
131 #ifdef debugFTFparams
132 G4cout << "Sqrt(s) " << SqrtS << G4endl;
133 #endif
134
135 TargetMass /= GeV; TargetMass2 /= (GeV*GeV);
136 ProjectileMass /= GeV; ProjectileMass2 /= (GeV*GeV);
137
138 Plab /= GeV;
139 G4double Xftf = 0.0;
140
141 G4int NumberOfTargetProtons = theZ;
142 G4int NumberOfTargetNeutrons = theA - theZ;
143 NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
144
145 // ---------- hadron projectile ----------------
146 if ( AbsProjectileBaryonNumber <= 1 ) { // Projectile is hadron or baryon
147
148 // Interaction on P
149 G4double xTtP = csGGinstance->GetTotalIsotopeCrossSection( particle, KineticEnergy, 1, 1);
150 G4double xElP = csGGinstance->GetElasticIsotopeCrossSection(particle, KineticEnergy, 1, 1);
151
152 // Interaction on N
153 G4double xTtN = csGGinstance->GetTotalIsotopeCrossSection( particle, KineticEnergy, 0, 1);
154 G4double xElN = csGGinstance->GetElasticIsotopeCrossSection(particle, KineticEnergy, 0, 1);
155
156 // Average properties of h+N interactions
157 Xtotal = ( NumberOfTargetProtons * xTtP + NumberOfTargetNeutrons * xTtN ) / NumberOfTargetNucleons;
158 Xelastic = ( NumberOfTargetProtons * xElP + NumberOfTargetNeutrons * xElN ) / NumberOfTargetNucleons;
159 Xannihilation = 0.0;
160
161 Xtotal /= millibarn;
162 Xelastic /= millibarn;
163
164 #ifdef debugFTFparams
165 G4cout<<"Estimated cross sections (total and elastic) of h+N interactions "<<Xtotal<<" "<<Xelastic<<" (mb)"<<G4endl;
166 #endif
167 }
168
169 // ---------- nucleus projectile ----------------
170 if ( ProjectileIsNucleus && ProjectileBaryonNumber > 1 ) {
171
172 #ifdef debugFTFparams
173 G4cout<<"Projectile is a nucleus: A and Z - "<<ProjectileBaryonNumber<<" "<<ProjectileCharge<<G4endl;
174 #endif
175
177 // Interaction on P
178 G4double XtotPP = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
179 G4double XelPP = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
180
182 // Interaction on N
183 G4double XtotPN = csGGinstance->GetTotalIsotopeCrossSection(Neutron, KineticEnergy, 0, 1);
184 G4double XelPN = csGGinstance->GetElasticIsotopeCrossSection(Neutron, KineticEnergy, 0, 1);
185
186 #ifdef debugFTFparams
187 G4cout << "XsPP (total and elastic) " << XtotPP/millibarn << " " << XelPP/millibarn <<" (mb)"<< G4endl
188 << "XsPN (total and elastic) " << XtotPN/millibarn << " " << XelPN/millibarn <<" (mb)"<< G4endl;
189 #endif
190
191 Xtotal = (
192 AbsProjectileCharge * NumberOfTargetProtons * XtotPP +
193 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
194 NumberOfTargetNeutrons * XtotPP
195 +
196 ( AbsProjectileCharge * NumberOfTargetNeutrons +
197 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
198 NumberOfTargetProtons ) * XtotPN
199 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
200 Xelastic= (
201 AbsProjectileCharge * NumberOfTargetProtons * XelPP +
202 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
203 NumberOfTargetNeutrons * XelPP
204 +
205 ( AbsProjectileCharge * NumberOfTargetNeutrons +
206 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
207 NumberOfTargetProtons ) * XelPN
208 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
209
210 Xannihilation = 0.0;
211 Xtotal /= millibarn;
212 Xelastic /= millibarn;
213 }
214
215 // ---------- The projectile is anti-baryon or anti-nucleus ----------------
216 // anti Sigma^0_c anti Delta^-
217 if ( ProjectilePDGcode >= -4112 && ProjectilePDGcode <= -1114 ) {
218 // Only non-strange and strange baryons are considered
219
220 #ifdef debugFTFparams
221 G4cout<<"Projectile is a anti-baryon or anti-nucleus - "<<ProjectileBaryonNumber<<" "<<ProjectileCharge<<G4endl;
222 G4cout<<"(Only non-strange and strange baryons are considered)"<<G4endl;
223 #endif
224
225 G4double X_a( 0.0 ), X_b( 0.0 ), X_c( 0.0 ), X_d( 0.0 );
226 G4double MesonProdThreshold = ProjectileMass + TargetMass +
227 ( 2.0 * 0.14 + 0.016 ); // 2 Mpi + DeltaE;
228
229 if ( PlabPerParticle < 40.0*MeV ) { // Low energy limits. Projectile at rest.
230 Xtotal = 1512.9; // mb
231 Xelastic = 473.2; // mb
232 X_a = 625.1; // mb
233 X_b = 9.780; // mb
234 X_c = 49.989; // mb
235 X_d = 6.614; // mb
236 } else { // Total and elastic cross section of PbarP interactions a'la Arkhipov
237 G4double LogS = G4Log( ECMSsqr / 33.0625 );
238 G4double Xasmpt = 36.04 + 0.304*LogS*LogS; // mb
239 LogS = G4Log( SqrtS / 20.74 );
240 G4double Basmpt = 11.92 + 0.3036*LogS*LogS; // GeV^(-2)
241 G4double R0 = std::sqrt( 0.40874044*Xasmpt - Basmpt ); // GeV^(-1)
242
243 G4double FlowF = SqrtS / std::sqrt( ECMSsqr*ECMSsqr + ProjectileMass2*ProjectileMass2 +
244 TargetMass2*TargetMass2 - 2.0*ECMSsqr*ProjectileMass2
245 - 2.0*ECMSsqr*TargetMass2
246 - 2.0*ProjectileMass2*TargetMass2 );
247
248 Xtotal = Xasmpt * ( 1.0 + 13.55*FlowF/R0/R0/R0*
249 (1.0 - 4.47/SqrtS + 12.38/ECMSsqr - 12.43/SqrtS/ECMSsqr) ); // mb
250
251 Xasmpt = 4.4 + 0.101*LogS*LogS; // mb
252 Xelastic = Xasmpt * ( 1.0 + 59.27*FlowF/R0/R0/R0*
253 (1.0 - 6.95/SqrtS + 23.54/ECMSsqr - 25.34/SqrtS/ECMSsqr ) ); // mb
254
255 X_a = 25.0*FlowF; // mb, 3-shirts diagram
256
257 if ( SqrtS < MesonProdThreshold ) {
258 X_b = 3.13 + 140.0*G4Pow::GetInstance()->powA( MesonProdThreshold - SqrtS, 2.5 ); // mb anti-quark-quark annihilation
259 Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
260 } else {
261 X_b = 6.8/SqrtS; // mb anti-quark-quark annihilation
262 Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
263 }
264
265 X_c = 2.0*FlowF*sqr( ProjectileMass + TargetMass )/ECMSsqr; // mb rearrangement
266
267 X_d = 23.3/ECMSsqr; // mb anti-quark-quark string creation
268 }
269
270 G4double Xann_on_P( 0.0), Xann_on_N( 0.0 );
271
272 if ( ProjectilePDGcode == -2212 ) { // Pbar+P/N
273 Xann_on_P = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
274 Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
275 } else if ( ProjectilePDGcode == -2112 ) { // NeutrBar+P/N
276 Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
277 Xann_on_N = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
278 } else if ( ProjectilePDGcode == -3122 ) { // LambdaBar+P/N
279 Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
280 Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
281 } else if ( ProjectilePDGcode == -3112 ) { // Sigma-Bar+P/N
282 Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
283 Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
284 } else if ( ProjectilePDGcode == -3212 ) { // Sigma0Bar+P/N
285 Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
286 Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
287 } else if ( ProjectilePDGcode == -3222 ) { // Sigma+Bar+P/N
288 Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
289 Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
290 } else if ( ProjectilePDGcode == -3312 ) { // Xi-Bar+P/N
291 Xann_on_P = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
292 Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
293 } else if ( ProjectilePDGcode == -3322 ) { // Xi0Bar+P/N
294 Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
295 Xann_on_N = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
296 } else if ( ProjectilePDGcode == -3334 ) { // Omega-Bar+P/N
297 Xann_on_P = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
298 Xann_on_N = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
299 } else {
300 G4cout << "Unknown anti-baryon for FTF annihilation" << G4endl;
301 }
302
303 //G4cout << "Sum " << Xann_on_P << G4endl;
304
305 if ( ! ProjectileIsNucleus ) { // Projectile is anti-baryon
306 Xannihilation = ( NumberOfTargetProtons * Xann_on_P + NumberOfTargetNeutrons * Xann_on_N )
307 / NumberOfTargetNucleons;
308 } else { // Projectile is a nucleus
309 Xannihilation = (
310 ( AbsProjectileCharge * NumberOfTargetProtons +
311 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
312 NumberOfTargetNeutrons ) * Xann_on_P
313 +
314 ( AbsProjectileCharge * NumberOfTargetNeutrons +
315 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
316 NumberOfTargetProtons ) * Xann_on_N
317 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
318 }
319
320 //G4double Xftf = 0.0;
321 MesonProdThreshold = ProjectileMass + TargetMass + (0.14 + 0.08); // Mpi + DeltaE
322 if ( SqrtS > MesonProdThreshold ) {
323 Xftf = 36.0 * ( 1.0 - MesonProdThreshold/SqrtS );
324 }
325
326 Xtotal = Xelastic + Xannihilation + Xftf;
327
328 #ifdef debugFTFparams
329 G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
330 << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << " (mb)"<< G4endl
331 << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal << " "
332 << Xannihilation/(Xtotal - Xelastic) << G4endl;
333 #endif
334
335 }
336
337 if ( Xtotal == 0.0 ) { // Projectile is undefined, Nucleon assumed
338
340 // Interaction on P
341 G4double XtotPP = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
342 G4double XelPP = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
343
344 // Interaction on N
345 G4double XtotPN = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 0, 1);
346 G4double XelPN = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 0, 1);
347
348 Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN )
349 / NumberOfTargetNucleons;
350 Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN )
351 / NumberOfTargetNucleons;
352 Xannihilation = 0.0;
353 Xtotal /= millibarn;
354 Xelastic /= millibarn;
355 };
356
357 // Geometrical parameters
358 SetTotalCrossSection( Xtotal );
359 SetElastisCrossSection( Xelastic );
360 SetInelasticCrossSection( Xtotal - Xelastic );
361
362 // Interactions with elastic and inelastic collisions
363 SetProbabilityOfElasticScatt( Xtotal, Xelastic );
364
365 SetRadiusOfHNinteractions2( Xtotal/pi/10.0 );
366
367 if ( ( Xtotal - Xelastic ) == 0.0 ) {
369 } else {
370 SetProbabilityOfAnnihilation( Xannihilation / (Xtotal - Xelastic) );
371 }
372
373 if(Xelastic > 0.0) {
374 SetSlope( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 );// Slope parameter of elastic scattering
375 // (GeV/c)^(-2))
376 // Parameters of elastic scattering
377 // Gaussian parametrization of elastic scattering amplitude assumed
378 SetAvaragePt2ofElasticScattering( 1.0/( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 )*GeV*GeV );
379 } else {
380 SetSlope(1.0);
382 }
383 SetGamma0( GetSlope()*Xtotal/10.0/2.0/pi );
384
385 G4double Xinel = Xtotal - Xelastic;
386
387 #ifdef debugFTFparams
388 G4cout<< "Slope of hN elastic scattering" << GetSlope() << G4endl;
389 G4cout << "AvaragePt2ofElasticScattering " << GetAvaragePt2ofElasticScattering() << G4endl;
390 G4cout<<"Parameters of excitation for projectile "<<ProjectilePDGcode<< G4endl;
391 #endif
392
393 if ( (ProjectilePDGcode == 2212) || (ProjectilePDGcode == 2112) ) { // Projectile is proton or neutron
394
395 // A process probability is parameterized as Prob = A_1*exp(-A_2*y) + A_3*exp(-A_4*y) + A_top
396 // y is a rapidity of a partcle in the target nucleus. Ymin is a minimal rapidity below it X=0
397
398 // Proc# A1 B1 A2 B2 A3 Atop Ymin
399 /* original hadr-string-diff-V10-03-07 (similar to 10.3.x)
400 SetParams( 0, 13.71, 1.75, -214.5, 4.25, 0.0, 0.5 , 1.1 ); // Qexchange without Exc.
401 SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
402 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
403 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
404 SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiplier
405 */
406
407 // Proc#
408 SetParams( 0, fParCollBaryonProj.GetProc0A1(), fParCollBaryonProj.GetProc0B1(),
409 fParCollBaryonProj.GetProc0A2(), fParCollBaryonProj.GetProc0B2(),
410 fParCollBaryonProj.GetProc0A3(), fParCollBaryonProj.GetProc0Atop(),
411 fParCollBaryonProj.GetProc0Ymin() ); // Qexchange without Exc.
412 SetParams( 1, fParCollBaryonProj.GetProc1A1(), fParCollBaryonProj.GetProc1B1(),
413 fParCollBaryonProj.GetProc1A2(), fParCollBaryonProj.GetProc1B2(),
414 fParCollBaryonProj.GetProc1A3(), fParCollBaryonProj.GetProc1Atop(),
415 fParCollBaryonProj.GetProc1Ymin() ); // Qexchange with Exc.
416 if ( Xinel > 0.0 ) {
417 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Projectile diffraction
418 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Target diffraction
419
420 SetParams( 4, fParCollBaryonProj.GetProc4A1(), fParCollBaryonProj.GetProc4B1(),
421 fParCollBaryonProj.GetProc4A2(), fParCollBaryonProj.GetProc4B2(),
422 fParCollBaryonProj.GetProc4A3(), fParCollBaryonProj.GetProc4Atop(),
423 fParCollBaryonProj.GetProc4Ymin() ); // Qexchange with Exc. Additional multiplier
424 } else { // if Xinel=0., zero everything out (obviously)
425 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
426 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
427 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
428 }
429
430 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
431 //
432 // It is not decided what to do with diffraction dissociation in Had-Nucl and Nucl-Nucl interactions
433 // For the moment both ProjDiffDisso & TgtDiffDisso for A > 10 are set to false,
434 // so both projectile and target diffraction are turned OFF
435 //
436 if ( ! fParCollBaryonProj.IsProjDiffDissociation() )
437 SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
438 if ( ! fParCollBaryonProj.IsTgtDiffDissociation() )
439 SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
440 }
441
443
444 if ( NumberOfTargetNucleons > 26 ) {
446 } else {
448 }
449
450 SetProjMinDiffMass( fParCollBaryonProj.GetProjMinDiffMass() ); // GeV
451 SetProjMinNonDiffMass( fParCollBaryonProj.GetProjMinNonDiffMass() ); // GeV
452
453 SetTarMinDiffMass( fParCollBaryonProj.GetTgtMinDiffMass() ); // GeV
454 SetTarMinNonDiffMass( fParCollBaryonProj.GetTgtMinNonDiffMass() ); // GeV
455
456 SetAveragePt2( fParCollBaryonProj.GetAveragePt2() ); // GeV^2
457 SetProbLogDistrPrD( fParCollBaryonProj.GetProbLogDistrPrD() );
458 SetProbLogDistr( fParCollBaryonProj.GetProbLogDistr() );
459
460 } else if ( ProjectilePDGcode == -2212 || ProjectilePDGcode == -2112 ) { // Projectile is anti_proton or anti_neutron
461
462 // Proc# A1 B1 A2 B2 A3 Atop Ymin
463 SetParams( 0, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange without Exc.
464 SetParams( 1, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange with Exc.
465 if ( Xinel > 0.) {
466 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Projectile diffraction
467 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Target diffraction
468 SetParams( 4, 1.0, 0.0 , 0.0, 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
469 } else {
470 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
471 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
472 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
473 }
474
475 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
476 //
477 // It is not decided what to do with diffraction dissociation in Had-Nucl and Nucl-Nucl interactions
478 // For the moment both ProjDiffDisso & TgtDiffDisso are set to false,
479 // so both projectile and target diffraction are turned OFF
480 //
481 if ( ! fParCollBaryonProj.IsProjDiffDissociation() )
482 SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
483 if ( ! fParCollBaryonProj.IsTgtDiffDissociation() )
484 SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
485 }
486
489 SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
490 SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
491 SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
492 SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
493 SetAveragePt2( 0.3 ); // GeV^2
494 SetProbLogDistrPrD( 0.55 );
495 SetProbLogDistr( 0.55 );
496
497 } else if ( ProjectileabsPDGcode == 211 || ProjectilePDGcode == 111 ) { // Projectile is Pion
498
499 // Proc# A1 B1 A2 B2 A3 Atop Ymin
500 /* --> original code
501 SetParams( 0, 150.0, 1.8 , -247.3, 2.3, 0., 1. , 2.3 );
502 SetParams( 1, 5.77, 0.6 , -5.77, 0.8, 0., 0. , 0.0 );
503 SetParams( 2, 2.27, 0.5 , -98052.0, 4.0, 0., 0. , 3.0 );
504 SetParams( 3, 7.0, 0.9, -85.28, 1.9, 0.08, 0. , 2.2 );
505 SetParams( 4, 1.0, 0.0 , -11.02, 1.0, 0.0, 0. , 2.4 ); // Qexchange with Exc. Additional multiply
506 */
507 // Proc#
508 SetParams( 0, fParCollPionProj.GetProc0A1(), fParCollPionProj.GetProc0B1(),
509 fParCollPionProj.GetProc0A2(), fParCollPionProj.GetProc0B2(),
510 fParCollPionProj.GetProc0A3(), fParCollPionProj.GetProc0Atop(),
511 fParCollPionProj.GetProc0Ymin() ); // Qexchange without Exc.
512 SetParams( 1, fParCollPionProj.GetProc1A1(), fParCollPionProj.GetProc1B1(),
513 fParCollPionProj.GetProc1A2(), fParCollPionProj.GetProc1B2(),
514 fParCollPionProj.GetProc1A3(), fParCollPionProj.GetProc1Atop(),
515 fParCollPionProj.GetProc1Ymin() ); // Qexchange with Exc.
516 SetParams( 2, fParCollPionProj.GetProc2A1(), fParCollPionProj.GetProc2B1(),
517 fParCollPionProj.GetProc2A2(), fParCollPionProj.GetProc2B2(),
518 fParCollPionProj.GetProc2A3(), fParCollPionProj.GetProc2Atop(),
519 fParCollPionProj.GetProc2Ymin() ); // Projectile diffraction
520 SetParams( 3, fParCollPionProj.GetProc3A1(), fParCollPionProj.GetProc3B1(),
521 fParCollPionProj.GetProc3A2(), fParCollPionProj.GetProc3B2(),
522 fParCollPionProj.GetProc3A3(), fParCollPionProj.GetProc3Atop(),
523 fParCollPionProj.GetProc3Ymin() ); // Target diffraction
524 SetParams( 4, fParCollPionProj.GetProc4A1(), fParCollPionProj.GetProc4B1(),
525 fParCollPionProj.GetProc4A2(), fParCollPionProj.GetProc4B2(),
526 fParCollPionProj.GetProc4A3(), fParCollPionProj.GetProc4Atop(),
527 fParCollPionProj.GetProc4Ymin() ); // Qexchange with Exc. Additional multiply
528
529 // NOTE: how can it be |ProjectileBaryonNumber| > 10 if projectile is a pion ???
530 //
531 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
532 if ( ! fParCollPionProj.IsProjDiffDissociation() )
533 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
534 if ( ! fParCollPionProj.IsTgtDiffDissociation() )
535 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
536 }
537
538 /* original code -->
539 SetDeltaProbAtQuarkExchange( 0.56 );
540 SetProjMinDiffMass( 1.0 ); // GeV
541 SetProjMinNonDiffMass( 1.0 ); // GeV
542 SetTarMinDiffMass( 1.16 ); // GeV
543 SetTarMinNonDiffMass( 1.16 ); // GeV
544 SetAveragePt2( 0.3 ); // GeV^2
545 SetProbLogDistrPrD( 0.55 );
546 SetProbLogDistr( 0.55 );
547 */
548
549 // JVY update, Aug.8, 2018 --> Feb.14, 2019
550 //
552 SetProjMinDiffMass( fParCollPionProj.GetProjMinDiffMass() ); // GeV
553 SetProjMinNonDiffMass( fParCollPionProj.GetProjMinNonDiffMass() ); // GeV
554 SetTarMinDiffMass( fParCollPionProj.GetTgtMinDiffMass() ); // GeV
555 SetTarMinNonDiffMass( fParCollPionProj.GetTgtMinNonDiffMass() ); // GeV
556 SetAveragePt2( fParCollPionProj.GetAveragePt2() ); // GeV^2
557 SetProbLogDistrPrD( fParCollPionProj.GetProbLogDistrPrD() );
558 SetProbLogDistr( fParCollPionProj.GetProbLogDistr() );
559
560 // ---> end update
561
562 } else if ( ProjectileabsPDGcode == 321 || ProjectileabsPDGcode == 311 ||
563 ProjectilePDGcode == 130 || ProjectilePDGcode == 310 ) { // Projectile is Kaon
564
565 // Proc# A1 B1 A2 B2 A3 Atop Ymin
566 SetParams( 0, 60.0 , 2.5 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
567 SetParams( 1, 6.0 , 1.0 , -24.33 , 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
568 SetParams( 2, 2.76, 1.2 , -22.5 , 2.7 ,0.04, 0.0 , 1.40 ); // Projectile diffraction
569 SetParams( 3, 1.09, 0.5 , -8.88 , 2. ,0.05, 0.0 , 1.40 ); // Target diffraction
570 SetParams( 4, 1.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
571 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
572 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
573 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
574 }
575
577 SetProjMinDiffMass( 0.7 ); // GeV
578 SetProjMinNonDiffMass( 0.7 ); // GeV
579 SetTarMinDiffMass( 1.16 ); // GeV
580 SetTarMinNonDiffMass( 1.16 ); // GeV
581 SetAveragePt2( 0.3 ); // GeV^2
582 SetProbLogDistrPrD( 0.55 );
583 SetProbLogDistr( 0.55 );
584
585 } else { // Projectile is not p, n, Pi0, Pi+, Pi-, K+, K-, K0 or their anti-particles
586
587 if ( ProjectileabsPDGcode > 1000 ) { // The projectile is a baryon as P or N
588 // Proc# A1 B1 A2 B2 A3 Atop Ymin
589 SetParams( 0, 13.71, 1.75, -30.69, 3.0 , 0.0, 1.0 , 0.93 ); // Qexchange without Exc.
590 SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
591 if ( Xinel > 0.) {
592 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
593 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
594 SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiply
595 } else {
596 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
597 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
598 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
599 }
600
601 } else { // The projectile is a meson as K+-0
602 // Proc# A1 B1 A2 B2 A3 Atop Ymin
603 SetParams( 0, 60.0 , 2.5 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
604 SetParams( 1, 6.0 , 1.0 , -24.33 , 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
605 SetParams( 2, 2.76, 1.2 , -22.5 , 2.7 ,0.04, 0.0 , 1.40 ); // Projectile diffraction
606 SetParams( 3, 1.09, 0.5 , -8.88 , 2. ,0.05, 0.0 , 1.40 ); // Target diffraction
607 SetParams( 4, 1.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
608 }
609
610 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
611 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
612 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
613 }
614
617
618 SetProjMinDiffMass( GetMinMass(particle)/GeV );
619 SetProjMinNonDiffMass( GetMinMass(particle)/GeV );
620
622 SetTarMinDiffMass( GetMinMass(Neutron)/GeV );
623 SetTarMinNonDiffMass( GetMinMass(Neutron)/GeV );
624
625 SetAveragePt2( 0.3 ); // GeV^2
626 SetProbLogDistrPrD( 0.55 );
627 SetProbLogDistr( 0.55 );
628
629 }
630
631 #ifdef debugFTFparams
632 G4cout<<"DeltaProbAtQuarkExchange "<< GetDeltaProbAtQuarkExchange() << G4endl;
633 G4cout<<"ProbOfSameQuarkExchange "<< GetProbOfSameQuarkExchange() << G4endl;
634 G4cout<<"ProjMinDiffMass "<< GetProjMinDiffMass()/GeV <<" GeV"<< G4endl;
635 G4cout<<"ProjMinNonDiffMass "<< GetProjMinNonDiffMass() <<" GeV"<< G4endl;
636 G4cout<<"TarMinDiffMass "<< GetTarMinDiffMass() <<" GeV"<< G4endl;
637 G4cout<<"TarMinNonDiffMass "<< GetTarMinNonDiffMass() <<" GeV"<< G4endl;
638 G4cout<<"AveragePt2 "<< GetAveragePt2() <<" GeV^2"<< G4endl;
639 G4cout<<"ProbLogDistrPrD "<< GetProbLogDistrPrD() << G4endl;
640 G4cout<<"ProbLogDistrTrD "<< GetProbLogDistr() << G4endl;
641 #endif
642
643 // Set parameters of nuclear destruction
644
645 if ( ProjectileabsPDGcode < 1000 ) { // Meson projectile
646
647 SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
648 //
649 // target destruction
650 //
651 /* original code --->
652 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)*
653 G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
654
655 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
656 SetDofNuclearDestruction( 0.3 );
657 SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
658 ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
659 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
660 SetExcitationEnergyPerWoundedNucleon( 40.0*MeV );
661 */
662 double coeff = fParCollMesonProj.GetNuclearTgtDestructP1();
663 //
664 // NOTE (JVY): Set this switch to false/true on line 138
665 //
666 if ( fParCollMesonProj.IsNuclearTgtDestructP1_ADEP() )
667 {
668 coeff *= G4double(NumberOfTargetNucleons);
669 }
670 double exfactor = G4Exp( fParCollMesonProj.GetNuclearTgtDestructP2()
671 * (Ylab-fParCollMesonProj.GetNuclearTgtDestructP3()) );
672 coeff *= exfactor;
673 coeff /= ( 1.+ exfactor );
674
676
678 SetDofNuclearDestruction( fParCollMesonProj.GetDofNuclearDestruct() );
679 coeff = fParCollMesonProj.GetPt2NuclearDestructP2();
680 exfactor = G4Exp( fParCollMesonProj.GetPt2NuclearDestructP3()
681 * (Ylab-fParCollMesonProj.GetPt2NuclearDestructP4()) );
682 coeff *= exfactor;
683 coeff /= ( 1. + exfactor );
684 SetPt2ofNuclearDestruction( (fParCollMesonProj.GetPt2NuclearDestructP1()+coeff)*CLHEP::GeV*CLHEP::GeV );
685
688
689 } else if ( ProjectilePDGcode == -2212 || ProjectilePDGcode == -2112 ) { // for anti-baryon projectile
690
691 SetMaxNumberOfCollisions( Plab, 2.0 );
692
693 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)*
694 G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
695 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
697 SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
698 ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
699 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
701 if ( Plab < 2.0 ) { // 2 GeV/c
702 // For slow anti-baryon we have to garanty putting on mass-shell
704 SetR2ofNuclearDestruction( 1.5*fermi*fermi ); // this is equivalent to setting a few line above
705 // is it even necessary ?
707 SetPt2ofNuclearDestruction( 0.035*GeV*GeV );
708 SetMaxPt2ofNuclearDestruction( 0.04*GeV*GeV );
709 }
710
711 } else { // Projectile baryon assumed
712
713 // NOTE (JVY) FIXME !!! Will decide later how/if to make this one configurable...
714 //
715 SetMaxNumberOfCollisions( Plab, 2.0 );
716
717 // projectile destruction - does NOT really matter for particle projectile, only for a nucleus projectile
718 //
719 double coeff = 0.;
720 coeff = fParCollBaryonProj.GetNuclearProjDestructP1();
721 //
722 // NOTE (JVY): Set this switch to false/true on line 136
723 //
724 if ( fParCollBaryonProj.IsNuclearProjDestructP1_NBRNDEP() )
725 {
726 coeff *= G4double(AbsProjectileBaryonNumber);
727 }
728 double exfactor = G4Exp( fParCollBaryonProj.GetNuclearProjDestructP2()*(Ylab-fParCollBaryonProj.GetNuclearProjDestructP3()) );
729 coeff *= exfactor;
730 coeff /= ( 1.+ exfactor );
732
733 // target desctruction
734 //
735 coeff = fParCollBaryonProj.GetNuclearTgtDestructP1();
736 //
737 // NOTE (JVY): Set this switch to false/true on line 138
738 //
739 if ( fParCollBaryonProj.IsNuclearTgtDestructP1_ADEP() )
740 {
741 coeff *= G4double(NumberOfTargetNucleons);
742 }
743 exfactor = G4Exp( fParCollBaryonProj.GetNuclearTgtDestructP2()*(Ylab-fParCollBaryonProj.GetNuclearTgtDestructP3()) );
744 coeff *= exfactor;
745 coeff /= ( 1.+ exfactor );
747
749 SetDofNuclearDestruction( fParCollBaryonProj.GetDofNuclearDestruct() );
750
751 coeff = fParCollBaryonProj.GetPt2NuclearDestructP2();
752 exfactor = G4Exp( fParCollBaryonProj.GetPt2NuclearDestructP3()*(Ylab-fParCollBaryonProj.GetPt2NuclearDestructP4()) );
753 coeff *= exfactor;
754 coeff /= ( 1. + exfactor );
755 SetPt2ofNuclearDestruction( (fParCollBaryonProj.GetPt2NuclearDestructP1()+coeff)*CLHEP::GeV*CLHEP::GeV );
756
759
760 }
761
762 #ifdef debugFTFparams
763 G4cout<<"CofNuclearDestructionPr "<< GetCofNuclearDestructionPr() << G4endl;
764 G4cout<<"CofNuclearDestructionTr "<< GetCofNuclearDestruction() << G4endl;
765 G4cout<<"R2ofNuclearDestruction "<< GetR2ofNuclearDestruction()/fermi/fermi <<" fermi^2"<< G4endl;
766 G4cout<<"DofNuclearDestruction "<< GetDofNuclearDestruction() << G4endl;
767 G4cout<<"Pt2ofNuclearDestruction "<< GetPt2ofNuclearDestruction()/GeV/GeV <<" GeV^2"<< G4endl;
768 G4cout<<"ExcitationEnergyPerWoundedNucleon "<< GetExcitationEnergyPerWoundedNucleon() <<" MeV"<< G4endl;
769 #endif
770
771 //SetCofNuclearDestruction( 0.47*G4Exp( 2.0*(Ylab - 2.5) )/( 1.0 + G4Exp( 2.0*(Ylab - 2.5) ) ) );
772 //SetPt2ofNuclearDestruction( ( 0.035 + 0.1*G4Exp( 4.0*(Ylab - 3.0) )/( 1.0 + G4Exp( 4.0*(Ylab - 3.0) ) ) )*GeV*GeV );
773
774 //SetMagQuarkExchange( 120.0 ); // 210.0 PipP
775 //SetSlopeQuarkExchange( 2.0 );
776 //SetDeltaProbAtQuarkExchange( 0.6 );
777 //SetProjMinDiffMass( 0.7 ); // GeV 1.1
778 //SetProjMinNonDiffMass( 0.7 ); // GeV
779 //SetProbabilityOfProjDiff( 0.0); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
780 //SetTarMinDiffMass( 1.1 ); // GeV
781 //SetTarMinNonDiffMass( 1.1 ); // GeV
782 //SetProbabilityOfTarDiff( 0.0 ); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
783
784 //SetAveragePt2( 0.0 ); // GeV^2 0.3
785 //------------------------------------
786 //SetProbabilityOfElasticScatt( 1.0, 1.0); //(Xtotal, Xelastic);
787 //SetProbabilityOfProjDiff( 1.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 0->1
788 //SetProbabilityOfTarDiff( 4.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 2->4
789 //SetAveragePt2( 0.3 ); // (0.15)
790 //SetAvaragePt2ofElasticScattering( 0.0 );
791
792 //SetMaxNumberOfCollisions( Plab, 6.0 ); //(4.0*(Plab + 0.01), Plab); // 6.0 );
793 //SetAveragePt2( 0.15 );
794 //SetCofNuclearDestruction(-1.);//( 0.75 ); // (0.25)
795 //SetExcitationEnergyPerWoundedNucleon(0.);//( 30.0*MeV ); // (75.0*MeV)
796 //SetDofNuclearDestruction(0.);//( 0.2 ); //0.4 // 0.3 0.5
797
798 /*
799 SetAveragePt2(0.3);
800 SetCofNuclearDestructionPr(0.);
801 SetCofNuclearDestruction(0.); //( 0.5 ); (0.25)
802 SetExcitationEnergyPerWoundedNucleon(0.); // 30.0*MeV; (75.0*MeV)
803 SetDofNuclearDestruction(0.); // 0.2; 0.4; 0.3; 0.5
804 SetPt2ofNuclearDestruction(0.); //(2.*0.075*GeV*GeV); ( 0.3*GeV*GeV ); (0.168*GeV*GeV)
805 */
806
807 //SetExcitationEnergyPerWoundedNucleon(0.001);
808 //SetPt2Kink( 0.0*GeV*GeV );
809
810 //SetRadiusOfHNinteractions2( Xtotal/pi/10.0 /2.);
811 //SetRadiusOfHNinteractions2( (Xtotal - Xelastic)/pi/10.0 );
812 //SetProbabilityOfElasticScatt( 1.0, 0.0);
813 /*
814 G4cout << "Pt2 " << GetAveragePt2()<<" "<<GetAveragePt2()/GeV/GeV<<G4endl;
815 G4cout << "Cnd " << GetCofNuclearDestruction() << G4endl;
816 G4cout << "Dnd " << GetDofNuclearDestruction() << G4endl;
817 G4cout << "Pt2 " << GetPt2ofNuclearDestruction()/GeV/GeV << G4endl;
818 */
819
820}
double S(double temp)
G4double G4Log(G4double x)
Definition: G4Log.hh:226
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double GetProc3B2() const
double GetProc4A3() const
double GetNuclearTgtDestructP3() const
double GetProc3B1() const
bool IsNuclearProjDestructP1_NBRNDEP() const
double GetProc0A2() const
double GetExciEnergyPerWoundedNucleon() const
double GetProc1Atop() const
double GetProc2A3() const
double GetNuclearProjDestructP2() const
double GetProbLogDistrPrD() const
double GetPt2NuclearDestructP4() const
double GetPt2NuclearDestructP3() const
double GetNuclearProjDestructP3() const
double GetProc4A2() const
double GetPt2NuclearDestructP1() const
double GetProc2Ymin() const
double GetProjMinNonDiffMass() const
double GetProc3Ymin() const
double GetProc2A2() const
double GetProc4B1() const
double GetProc1A1() const
double GetProc0B2() const
double GetProbLogDistr() const
double GetProc1B2() const
double GetProc3Atop() const
double GetProc1A3() const
double GetProc1B1() const
double GetProc1Ymin() const
double GetProc2B1() const
double GetNuclearTgtDestructP2() const
double GetPt2NuclearDestructP2() const
double GetR2ofNuclearDestruct() const
double GetDeltaProbAtQuarkExchange() const
double GetDofNuclearDestruct() const
double GetProc0A3() const
bool IsProjDiffDissociation() const
double GetProc0Atop() const
double GetProc2B2() const
double GetProc2A1() const
double GetProc2Atop() const
double GetProc1A2() const
double GetProc4Atop() const
double GetProc3A1() const
double GetAveragePt2() const
double GetProc4Ymin() const
double GetProc0B1() const
double GetTgtMinNonDiffMass() const
bool IsTgtDiffDissociation() const
double GetProbOfSameQuarkExchange() const
double GetProc4B2() const
double GetProc0Ymin() const
double GetProc3A2() const
double GetProc3A3() const
bool IsNuclearTgtDestructP1_ADEP() const
double GetMaxPt2ofNuclearDestruct() const
double GetProc4A1() const
double GetTgtMinDiffMass() const
double GetNuclearTgtDestructP1() const
double GetProjMinDiffMass() const
double GetProc0A1() const
double GetNuclearProjDestructP1() const
void SetTarMinNonDiffMass(const G4double aValue)
void SetProjMinDiffMass(const G4double aValue)
G4double GetProbLogDistrPrD()
G4double GetCofNuclearDestructionPr()
void SetTotalCrossSection(const G4double Xtotal)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
G4double GetAveragePt2()
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void SetElastisCrossSection(const G4double Xelastic)
G4double GetPt2ofNuclearDestruction()
G4double GetProbLogDistr()
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
G4double GetProjMinNonDiffMass()
G4double GetAvaragePt2ofElasticScattering()
void SetSlope(const G4double Slope)
G4double GetTarMinDiffMass()
void SetDeltaProbAtQuarkExchange(const G4double aValue)
G4double GetTarMinNonDiffMass()
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
void SetTarMinDiffMass(const G4double aValue)
void SetProjMinNonDiffMass(const G4double aValue)
void SetCofNuclearDestructionPr(const G4double aValue)
void SetR2ofNuclearDestruction(const G4double aValue)
void SetProbLogDistrPrD(const G4double aValue)
void SetGamma0(const G4double Gamma0)
void SetAveragePt2(const G4double aValue)
G4double GetDeltaProbAtQuarkExchange()
G4double GetExcitationEnergyPerWoundedNucleon()
void SetProbabilityOfAnnihilation(const G4double aValue)
void SetDofNuclearDestruction(const G4double aValue)
G4double GetProjMinDiffMass()
G4double GetProbOfSameQuarkExchange()
G4double GetDofNuclearDestruction()
void SetRadiusOfHNinteractions2(const G4double Radius2)
void SetParams(const G4int ProcN, const G4double A1, const G4double B1, const G4double A2, const G4double B2, const G4double A3, const G4double Atop, const G4double Ymin)
void SetPt2ofNuclearDestruction(const G4double aValue)
G4double GetR2ofNuclearDestruction()
G4double GetCofNuclearDestruction()
void SetProbLogDistr(const G4double aValue)
void SetCofNuclearDestruction(const G4double aValue)
void SetMaxNumberOfCollisions(const G4double aValue, const G4double bValue)
void SetProbOfSameQuarkExchange(const G4double aValue)
void SetInelasticCrossSection(const G4double Xinelastic)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)=0
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)=0
T sqr(const T &x)
Definition: templates.hh:128

Referenced by G4FTFModel::Init().

◆ SetAvaragePt2ofElasticScattering()

void G4FTFParameters::SetAvaragePt2ofElasticScattering ( const G4double  aPt2)
inline

Definition at line 540 of file G4FTFParameters.hh.

540 {
542}

Referenced by InitForInteraction().

◆ SetAveragePt2()

void G4FTFParameters::SetAveragePt2 ( const G4double  aValue)
inline

Definition at line 580 of file G4FTFParameters.hh.

580 {
581 AveragePt2 = aValue*CLHEP::GeV*CLHEP::GeV;
582}

Referenced by InitForInteraction().

◆ SetCofNuclearDestruction()

void G4FTFParameters::SetCofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 628 of file G4FTFParameters.hh.

628 {
629 CofNuclearDestruction = aValue;
630}

Referenced by InitForInteraction().

◆ SetCofNuclearDestructionPr()

void G4FTFParameters::SetCofNuclearDestructionPr ( const G4double  aValue)
inline

Definition at line 624 of file G4FTFParameters.hh.

624 {
626}

Referenced by InitForInteraction().

◆ SetDeltaProbAtQuarkExchange()

void G4FTFParameters::SetDeltaProbAtQuarkExchange ( const G4double  aValue)
inline

Definition at line 556 of file G4FTFParameters.hh.

556 {
558}

Referenced by InitForInteraction().

◆ SetDofNuclearDestruction()

void G4FTFParameters::SetDofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 640 of file G4FTFParameters.hh.

640 {
641 DofNuclearDestruction = aValue;
642}

Referenced by InitForInteraction().

◆ SetElastisCrossSection()

void G4FTFParameters::SetElastisCrossSection ( const G4double  Xelastic)
inline

Definition at line 502 of file G4FTFParameters.hh.

502 {
503 FTFXelastic = Xelastic;
504}

Referenced by InitForInteraction().

◆ SetExcitationEnergyPerWoundedNucleon()

void G4FTFParameters::SetExcitationEnergyPerWoundedNucleon ( const G4double  aValue)
inline

Definition at line 636 of file G4FTFParameters.hh.

636 {
638}

Referenced by InitForInteraction().

◆ SetGamma0()

void G4FTFParameters::SetGamma0 ( const G4double  Gamma0)
inline

Definition at line 535 of file G4FTFParameters.hh.

535 {
536 FTFGamma0 = Gamma0;
537}

Referenced by InitForInteraction().

◆ SethNcmsEnergy()

void G4FTFParameters::SethNcmsEnergy ( const G4double  s)
inline

Definition at line 492 of file G4FTFParameters.hh.

492 {
494}

◆ SetInelasticCrossSection()

void G4FTFParameters::SetInelasticCrossSection ( const G4double  Xinelastic)
inline

Definition at line 506 of file G4FTFParameters.hh.

506 {
507 FTFXinelastic = Xinelastic;
508}

Referenced by InitForInteraction().

◆ SetMaxNumberOfCollisions()

void G4FTFParameters::SetMaxNumberOfCollisions ( const G4double  aValue,
const G4double  bValue 
)
inline

Definition at line 607 of file G4FTFParameters.hh.

608 {
609 if ( Plab > Pbound ) {
610 MaxNumberOfCollisions = Plab/Pbound;
611 SetProbOfInteraction( -1.0 );
612 } else {
613 //MaxNumberOfCollisions = -1.0;
614 //SetProbOfInteraction( G4Exp( 0.25*(Plab-Pbound) ) );
616 SetProbOfInteraction( -1.0 );
617 }
618}
void SetProbOfInteraction(const G4double aValue)

Referenced by InitForInteraction().

◆ SetMaxPt2ofNuclearDestruction()

void G4FTFParameters::SetMaxPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 648 of file G4FTFParameters.hh.

648 {
650}

Referenced by InitForInteraction().

◆ SetParams()

void G4FTFParameters::SetParams ( const G4int  ProcN,
const G4double  A1,
const G4double  B1,
const G4double  A2,
const G4double  B2,
const G4double  A3,
const G4double  Atop,
const G4double  Ymin 
)
inline

Definition at line 546 of file G4FTFParameters.hh.

549 {
550 ProcParams[ProcN][0] = A1; ProcParams[ProcN][1] = B1;
551 ProcParams[ProcN][2] = A2; ProcParams[ProcN][3] = B2;
552 ProcParams[ProcN][4] = A3;
553 ProcParams[ProcN][5] = Atop; ProcParams[ProcN][6] = Ymin;
554}

Referenced by InitForInteraction().

◆ SetProbabilityOfAnnihilation()

void G4FTFParameters::SetProbabilityOfAnnihilation ( const G4double  aValue)
inline

Definition at line 523 of file G4FTFParameters.hh.

523 {
525}

Referenced by InitForInteraction().

◆ SetProbabilityOfElasticScatt() [1/2]

void G4FTFParameters::SetProbabilityOfElasticScatt ( const G4double  aValue)
inline

Definition at line 519 of file G4FTFParameters.hh.

519 {
521}

◆ SetProbabilityOfElasticScatt() [2/2]

void G4FTFParameters::SetProbabilityOfElasticScatt ( const G4double  Xtotal,
const G4double  Xelastic 
)
inline

Definition at line 510 of file G4FTFParameters.hh.

511 {
512 if ( Xtotal == 0.0 ) {
514 } else {
515 ProbabilityOfElasticScatt = Xelastic / Xtotal;
516 }
517}

Referenced by G4FTFModel::Init(), and InitForInteraction().

◆ SetProbLogDistr()

void G4FTFParameters::SetProbLogDistr ( const G4double  aValue)
inline

Definition at line 588 of file G4FTFParameters.hh.

588 {
589 ProbLogDistr = aValue;
590}

Referenced by InitForInteraction().

◆ SetProbLogDistrPrD()

void G4FTFParameters::SetProbLogDistrPrD ( const G4double  aValue)
inline

Definition at line 584 of file G4FTFParameters.hh.

584 {
585 ProbLogDistrPrD = aValue;
586}

Referenced by InitForInteraction().

◆ SetProbOfInteraction()

void G4FTFParameters::SetProbOfInteraction ( const G4double  aValue)
inline

Definition at line 620 of file G4FTFParameters.hh.

620 {
621 ProbOfInelInteraction = aValue;
622}

Referenced by SetMaxNumberOfCollisions().

◆ SetProbOfSameQuarkExchange()

void G4FTFParameters::SetProbOfSameQuarkExchange ( const G4double  aValue)
inline

Definition at line 560 of file G4FTFParameters.hh.

560 {
562}

Referenced by InitForInteraction().

◆ SetProjMinDiffMass()

void G4FTFParameters::SetProjMinDiffMass ( const G4double  aValue)
inline

Definition at line 564 of file G4FTFParameters.hh.

564 {
565 ProjMinDiffMass = aValue*CLHEP::GeV;
566}

Referenced by InitForInteraction().

◆ SetProjMinNonDiffMass()

void G4FTFParameters::SetProjMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 568 of file G4FTFParameters.hh.

568 {
569 ProjMinNonDiffMass = aValue*CLHEP::GeV;
570}

Referenced by InitForInteraction().

◆ SetPt2Kink()

void G4FTFParameters::SetPt2Kink ( const G4double  aValue)
inline

Definition at line 594 of file G4FTFParameters.hh.

594 {
595 Pt2kink = aValue;
596}

Referenced by G4FTFParameters().

◆ SetPt2ofNuclearDestruction()

void G4FTFParameters::SetPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 644 of file G4FTFParameters.hh.

644 {
646}

Referenced by InitForInteraction().

◆ SetQuarkProbabilitiesAtGluonSplitUp()

void G4FTFParameters::SetQuarkProbabilitiesAtGluonSplitUp ( const G4double  Puubar,
const G4double  Pddbar,
const G4double  Pssbar 
)
inline

Definition at line 598 of file G4FTFParameters.hh.

600 {
601 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar );
602 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar );
603 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar + Pssbar );
604}

Referenced by G4FTFParameters().

◆ SetR2ofNuclearDestruction()

void G4FTFParameters::SetR2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 632 of file G4FTFParameters.hh.

632 {
633 R2ofNuclearDestruction = aValue;
634}

Referenced by InitForInteraction().

◆ SetRadiusOfHNinteractions2()

void G4FTFParameters::SetRadiusOfHNinteractions2 ( const G4double  Radius2)
inline

Definition at line 527 of file G4FTFParameters.hh.

527 {
528 RadiusOfHNinteractions2 = Radius2;
529}

Referenced by InitForInteraction().

◆ SetSlope()

void G4FTFParameters::SetSlope ( const G4double  Slope)
inline

Definition at line 531 of file G4FTFParameters.hh.

531 {
532 FTFSlope = 12.84 / Slope; // Slope is in GeV^-2, FTFSlope in fm^-2
533}

Referenced by InitForInteraction().

◆ SetTarMinDiffMass()

void G4FTFParameters::SetTarMinDiffMass ( const G4double  aValue)
inline

Definition at line 572 of file G4FTFParameters.hh.

572 {
573 TarMinDiffMass = aValue*CLHEP::GeV;
574}

Referenced by InitForInteraction().

◆ SetTarMinNonDiffMass()

void G4FTFParameters::SetTarMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 576 of file G4FTFParameters.hh.

576 {
577 TarMinNonDiffMass = aValue*CLHEP::GeV;
578}

Referenced by InitForInteraction().

◆ SetTotalCrossSection()

void G4FTFParameters::SetTotalCrossSection ( const G4double  Xtotal)
inline

Definition at line 498 of file G4FTFParameters.hh.

498 {
499 FTFXtotal = Xtotal;
500}

Referenced by InitForInteraction().

Member Data Documentation

◆ AvaragePt2ofElasticScattering

G4double G4FTFParameters::AvaragePt2ofElasticScattering

◆ AveragePt2

G4double G4FTFParameters::AveragePt2

Definition at line 448 of file G4FTFParameters.hh.

Referenced by GetAveragePt2(), and SetAveragePt2().

◆ CofNuclearDestruction

G4double G4FTFParameters::CofNuclearDestruction

Definition at line 460 of file G4FTFParameters.hh.

Referenced by GetCofNuclearDestruction(), and SetCofNuclearDestruction().

◆ CofNuclearDestructionPr

G4double G4FTFParameters::CofNuclearDestructionPr

Definition at line 459 of file G4FTFParameters.hh.

Referenced by GetCofNuclearDestructionPr(), and SetCofNuclearDestructionPr().

◆ DeltaProbAtQuarkExchange

G4double G4FTFParameters::DeltaProbAtQuarkExchange

Definition at line 439 of file G4FTFParameters.hh.

Referenced by GetDeltaProbAtQuarkExchange(), and SetDeltaProbAtQuarkExchange().

◆ DofNuclearDestruction

G4double G4FTFParameters::DofNuclearDestruction

Definition at line 465 of file G4FTFParameters.hh.

Referenced by GetDofNuclearDestruction(), and SetDofNuclearDestruction().

◆ ExcitationEnergyPerWoundedNucleon

G4double G4FTFParameters::ExcitationEnergyPerWoundedNucleon

◆ FTFGamma0

G4double G4FTFParameters::FTFGamma0

Definition at line 434 of file G4FTFParameters.hh.

Referenced by GammaElastic(), and SetGamma0().

◆ FTFhNcmsEnergy

G4double G4FTFParameters::FTFhNcmsEnergy

Definition at line 422 of file G4FTFParameters.hh.

Referenced by SethNcmsEnergy().

◆ FTFSlope

G4double G4FTFParameters::FTFSlope

Definition at line 432 of file G4FTFParameters.hh.

Referenced by GammaElastic(), GetSlope(), and SetSlope().

◆ FTFXannihilation

G4double G4FTFParameters::FTFXannihilation

Definition at line 428 of file G4FTFParameters.hh.

◆ FTFXelastic

G4double G4FTFParameters::FTFXelastic

Definition at line 426 of file G4FTFParameters.hh.

Referenced by GetElasticCrossSection(), and SetElastisCrossSection().

◆ FTFXinelastic

G4double G4FTFParameters::FTFXinelastic

Definition at line 427 of file G4FTFParameters.hh.

Referenced by GetInelasticCrossSection(), and SetInelasticCrossSection().

◆ FTFXtotal

G4double G4FTFParameters::FTFXtotal

Definition at line 425 of file G4FTFParameters.hh.

Referenced by GetTotalCrossSection(), and SetTotalCrossSection().

◆ MaxNumberOfCollisions

G4double G4FTFParameters::MaxNumberOfCollisions

Definition at line 456 of file G4FTFParameters.hh.

Referenced by GetMaxNumberOfCollisions(), and SetMaxNumberOfCollisions().

◆ MaxPt2ofNuclearDestruction

G4double G4FTFParameters::MaxPt2ofNuclearDestruction

◆ ProbabilityOfAnnihilation

G4double G4FTFParameters::ProbabilityOfAnnihilation

◆ ProbabilityOfElasticScatt

G4double G4FTFParameters::ProbabilityOfElasticScatt

◆ ProbLogDistr

G4double G4FTFParameters::ProbLogDistr

Definition at line 449 of file G4FTFParameters.hh.

Referenced by GetProbLogDistr(), and SetProbLogDistr().

◆ ProbLogDistrPrD

G4double G4FTFParameters::ProbLogDistrPrD

Definition at line 444 of file G4FTFParameters.hh.

Referenced by GetProbLogDistrPrD(), and SetProbLogDistrPrD().

◆ ProbOfInelInteraction

G4double G4FTFParameters::ProbOfInelInteraction

Definition at line 457 of file G4FTFParameters.hh.

Referenced by GetProbOfInteraction(), and SetProbOfInteraction().

◆ ProbOfSameQuarkExchange

G4double G4FTFParameters::ProbOfSameQuarkExchange

Definition at line 440 of file G4FTFParameters.hh.

Referenced by GetProbOfSameQuarkExchange(), and SetProbOfSameQuarkExchange().

◆ ProcParams

G4double G4FTFParameters::ProcParams[5][7]

Definition at line 437 of file G4FTFParameters.hh.

Referenced by GetProcProb(), and SetParams().

◆ ProjMinDiffMass

G4double G4FTFParameters::ProjMinDiffMass

Definition at line 442 of file G4FTFParameters.hh.

Referenced by GetProjMinDiffMass(), and SetProjMinDiffMass().

◆ ProjMinNonDiffMass

G4double G4FTFParameters::ProjMinNonDiffMass

Definition at line 443 of file G4FTFParameters.hh.

Referenced by GetProjMinNonDiffMass(), and SetProjMinNonDiffMass().

◆ Pt2kink

G4double G4FTFParameters::Pt2kink

Definition at line 452 of file G4FTFParameters.hh.

Referenced by GetPt2Kink(), and SetPt2Kink().

◆ Pt2ofNuclearDestruction

G4double G4FTFParameters::Pt2ofNuclearDestruction

Definition at line 466 of file G4FTFParameters.hh.

Referenced by GetPt2ofNuclearDestruction(), and SetPt2ofNuclearDestruction().

◆ QuarkProbabilitiesAtGluonSplitUp

std::vector< G4double > G4FTFParameters::QuarkProbabilitiesAtGluonSplitUp

◆ R2ofNuclearDestruction

G4double G4FTFParameters::R2ofNuclearDestruction

Definition at line 461 of file G4FTFParameters.hh.

Referenced by GetR2ofNuclearDestruction(), and SetR2ofNuclearDestruction().

◆ RadiusOfHNinteractions2

G4double G4FTFParameters::RadiusOfHNinteractions2

Definition at line 431 of file G4FTFParameters.hh.

Referenced by GetProbabilityOfInteraction(), and SetRadiusOfHNinteractions2().

◆ TarMinDiffMass

G4double G4FTFParameters::TarMinDiffMass

Definition at line 445 of file G4FTFParameters.hh.

Referenced by GetTarMinDiffMass(), and SetTarMinDiffMass().

◆ TarMinNonDiffMass

G4double G4FTFParameters::TarMinNonDiffMass

Definition at line 446 of file G4FTFParameters.hh.

Referenced by GetTarMinNonDiffMass(), and SetTarMinNonDiffMass().


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