Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Data Fields
G4FTFParameters Class Reference

#include <G4FTFParameters.hh>

Public Member Functions

 G4FTFParameters (const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
 
 ~G4FTFParameters ()
 
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 SetProbabilityOfProjDiff (const G4double aValue)
 
void SetTarMinDiffMass (const G4double aValue)
 
void SetTarMinNonDiffMass (const G4double aValue)
 
void SetProbabilityOfTarDiff (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 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 GetTarMinDiffMass ()
 
G4double GetTarMinNonDiffMass ()
 
G4double GetAveragePt2 ()
 
G4double GetProbLogDistr ()
 
G4double GetPt2Kink ()
 
std::vector< G4doubleGetQuarkProbabilitiesAtGluonSplitUp ()
 
G4double GetMaxNumberOfCollisions ()
 
G4double GetProbOfInteraction ()
 
G4double GetCofNuclearDestruction ()
 
G4double GetR2ofNuclearDestruction ()
 
G4double GetExcitationEnergyPerWoundedNucleon ()
 
G4double GetDofNuclearDestruction ()
 
G4double GetPt2ofNuclearDestruction ()
 
G4double GetMaxPt2ofNuclearDestruction ()
 
 G4FTFParameters ()
 

Data Fields

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

Detailed Description

Definition at line 40 of file G4FTFParameters.hh.

Constructor & Destructor Documentation

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

Definition at line 94 of file G4FTFParameters.cc.

References python.hepunit::fermi, FTFhNcmsEnergy, FTFXannihilation, FTFxsManager, G4cout, G4endl, G4ParticleDefinition::GetBaryonNumber(), G4ChipsComponentXS::GetElasticElementCrossSection(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGEncoding(), G4ParticleDefinition::GetPDGMass(), GetSlope(), G4ChipsComponentXS::GetTotalElementCrossSection(), python.hepunit::GeV, G4KaonMinus::KaonMinus(), G4KaonPlus::KaonPlus(), python.hepunit::MeV, python.hepunit::millibarn, G4INCL::Neutron, G4Neutron::Neutron(), python.hepunit::pi, G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), ProbOfSameQuarkExchange, G4INCL::Proton, G4Proton::Proton(), SetAvaragePt2ofElasticScattering(), SetAveragePt2(), SetCofNuclearDestruction(), SetDeltaProbAtQuarkExchange(), SetDofNuclearDestruction(), SetElastisCrossSection(), SetExcitationEnergyPerWoundedNucleon(), SetGamma0(), SetInelasticCrossSection(), SetMaxNumberOfCollisions(), SetMaxPt2ofNuclearDestruction(), SetParams(), SetProbabilityOfAnnihilation(), SetProbabilityOfElasticScatt(), SetProbLogDistr(), SetProbOfSameQuarkExchange(), SetProjMinDiffMass(), SetProjMinNonDiffMass(), SetPt2Kink(), SetPt2ofNuclearDestruction(), SetQuarkProbabilitiesAtGluonSplitUp(), SetR2ofNuclearDestruction(), SetRadiusOfHNinteractions2(), SetSlope(), SetTarMinDiffMass(), SetTarMinNonDiffMass(), SetTotalCrossSection(), and sqr().

95  {
96 
97  FTFXannihilation = 0.0;
98  FTFhNcmsEnergy = 0.0;
100 
101  G4int ProjectilePDGcode = particle->GetPDGEncoding();
102  G4int ProjectileabsPDGcode = std::abs( ProjectilePDGcode );
103  G4double ProjectileMass = particle->GetPDGMass();
104  G4double ProjectileMass2 = ProjectileMass * ProjectileMass;
105 
106  G4int ProjectileBaryonNumber( 0 ), AbsProjectileBaryonNumber( 0 ), AbsProjectileCharge( 0 );
107  G4bool ProjectileIsNucleus = false;
108 
109  if ( std::abs( particle->GetBaryonNumber() ) > 1 ) { // The projectile is a nucleus
110  ProjectileIsNucleus = true;
111  ProjectileBaryonNumber = particle->GetBaryonNumber();
112  AbsProjectileBaryonNumber = std::abs( ProjectileBaryonNumber );
113  AbsProjectileCharge = G4int( particle->GetPDGCharge() );
114  if ( ProjectileBaryonNumber > 1 ) {
115  ProjectilePDGcode = 2212; ProjectileabsPDGcode = 2212; // Proton
116  } else {
117  ProjectilePDGcode = -2212; ProjectileabsPDGcode = 2212; // Anti-Proton
118  }
119  ProjectileMass = G4Proton::Proton()->GetPDGMass();
120  ProjectileMass2 = sqr( ProjectileMass );
121  }
122 
123  G4double TargetMass = G4Proton::Proton()->GetPDGMass();
124  G4double TargetMass2 = TargetMass * TargetMass;
125 
126  G4double Plab = PlabPerParticle;
127  G4double Elab = std::sqrt( Plab*Plab + ProjectileMass2 );
128  G4double KineticEnergy = Elab - ProjectileMass;
129 
130  G4double S = ProjectileMass2 + TargetMass2 + 2.0*TargetMass*Elab;
131 
132  #ifdef debugFTFparams
133  G4cout << "--------- FTF Parameters --------------" << G4endl << "Proj Plab "
134  << ProjectilePDGcode << " " << Plab << G4endl << "Mass KinE " << ProjectileMass
135  << " " << KineticEnergy << G4endl << " A Z " << theA << " " << theZ << G4endl;
136  #endif
137 
138  G4double Ylab, Xtotal, Xelastic, Xannihilation;
139  G4int NumberOfTargetNucleons;
140 
141  Ylab = 0.5 * std::log( (Elab + Plab)/(Elab - Plab) );
142 
143  G4double ECMSsqr = S/GeV/GeV;
144  G4double SqrtS = std::sqrt( S )/GeV;
145 
146  #ifdef debugFTFparams
147  G4cout << "Sqrt(s) " << SqrtS << G4endl;
148  #endif
149 
150  TargetMass /= GeV; TargetMass2 /= (GeV*GeV);
151  ProjectileMass /= GeV; ProjectileMass2 /= (GeV*GeV);
152 
153  // Andrea Dotti (13Jan2013):
154  // The following lines are changed for G4MT. Originally the code was:
155  // static G4ChipsComponentXS* _instance = new G4ChipsComponentXS(); // Witek Pokorski
156  // Note the code could go back at original if _instance could be shared among threads
157  if ( ! chipsComponentXSisInitialized ) {
158  chipsComponentXSisInitialized = true;
159  chipsComponentXSinstance = new G4ChipsComponentXS();
160  }
161  G4ChipsComponentXS* _instance = chipsComponentXSinstance;
162  FTFxsManager = _instance;
163 
164  Plab /= GeV;
165  G4double Xftf = 0.0;
166  //G4double LogPlab = std::log( Plab );
167  //G4double sqrLogPlab = LogPlab * LogPlab;
168 
169  G4int NumberOfTargetProtons = theZ;
170  G4int NumberOfTargetNeutrons = theA - theZ;
171  NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
172 
173  if ( ProjectilePDGcode == 2212 || ProjectilePDGcode == 2112 ) { // Projectile is nucleon
174 
175  G4double XtotPP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
177  G4double XtotPN = FTFxsManager->GetTotalElementCrossSection( Neutron, KineticEnergy, 1, 0 );
178  G4double XelPP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
179  G4double XelPN = FTFxsManager->GetElasticElementCrossSection( Neutron, KineticEnergy, 1, 0 );
180 
181  #ifdef debugFTFparams
182  G4cout << "XsPP " << XtotPP/millibarn << " " << XelPP/millibarn << G4endl
183  << "XsPN " << XtotPN/millibarn << " " << XelPN/millibarn << G4endl;
184  #endif
185 
186  if ( ! ProjectileIsNucleus ) { // Projectile is hadron
187  Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN ) /
188  NumberOfTargetNucleons;
189  Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN ) /
190  NumberOfTargetNucleons;
191  } else { // Projectile is a nucleus
192  Xtotal = (
193  AbsProjectileCharge * NumberOfTargetProtons * XtotPP +
194  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
195  NumberOfTargetNeutrons * XtotPP
196  +
197  ( AbsProjectileCharge * NumberOfTargetNeutrons +
198  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
199  NumberOfTargetProtons ) * XtotPN
200  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
201  Xelastic= (
202  AbsProjectileCharge * NumberOfTargetProtons * XelPP +
203  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
204  NumberOfTargetNeutrons * XelPP
205  +
206  ( AbsProjectileCharge * NumberOfTargetNeutrons +
207  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
208  NumberOfTargetProtons ) * XelPN
209  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
210  }
211 
212  Xannihilation = 0.0;
213  Xtotal /= millibarn;
214  Xelastic /= millibarn;
215 
216  } else if ( ProjectilePDGcode < -1000 ) { // Projectile is anti_baryon
217 
218  G4double X_a( 0.0 ), X_b( 0.0 ), X_c( 0.0 ), X_d( 0.0 );
219  G4double MesonProdThreshold = ProjectileMass + TargetMass +
220  ( 2.0 * 0.14 + 0.016 ); // 2 Mpi + DeltaE;
221 
222  if ( PlabPerParticle < 40.0*MeV ) { // Low energy limits. Projectile at rest.
223  Xtotal = 1512.9; // mb
224  Xelastic = 473.2; // mb
225  X_a = 625.1; // mb
226  X_b = 9.780; // mb
227  X_c = 49.989; // mb
228  X_d = 6.614; // mb
229  } else { // Total and elastic cross section of PbarP interactions a'la Arkhipov
230  G4double LogS = std::log( ECMSsqr / 33.0625 );
231  G4double Xasmpt = 36.04 + 0.304*LogS*LogS; // mb
232  LogS = std::log( SqrtS / 20.74 );
233  G4double Basmpt = 11.92 + 0.3036*LogS*LogS; // GeV^(-2)
234  G4double R0 = std::sqrt( 0.40874044*Xasmpt - Basmpt ); // GeV^(-1)
235 
236  G4double FlowF = SqrtS / std::sqrt( ECMSsqr*ECMSsqr + ProjectileMass2*ProjectileMass2 +
237  TargetMass2*TargetMass2 - 2.0*ECMSsqr*ProjectileMass2
238  - 2.0*ECMSsqr*TargetMass2
239  - 2.0*ProjectileMass2*TargetMass2 );
240 
241  Xtotal = Xasmpt * ( 1.0 + 13.55*FlowF/R0/R0/R0*
242  (1.0 - 4.47/SqrtS + 12.38/ECMSsqr - 12.43/SqrtS/ECMSsqr) ); // mb
243 
244  Xasmpt = 4.4 + 0.101*LogS*LogS; // mb
245  Xelastic = Xasmpt * ( 1.0 + 59.27*FlowF/R0/R0/R0*
246  (1.0 - 6.95/SqrtS + 23.54/ECMSsqr - 25.34/SqrtS/ECMSsqr ) ); // mb
247 
248  //G4cout << "Param Xtotal Xelastic " << Xtotal << " " << Xelastic << G4endl
249  // << "FlowF " << FlowF << " SqrtS " << SqrtS << G4endl
250  // << "Param Xelastic-NaN " << Xelastic << " "
251  // << 1.5*16.654/pow(ECMSsqr/2.176/2.176,2.2) << " " << ECMSsqr << G4endl;
252 
253  X_a = 25.0*FlowF; // mb, 3-shirts diagram
254 
255  if ( SqrtS < MesonProdThreshold ) {
256  X_b = 3.13 + 140.0*std::pow( MesonProdThreshold - SqrtS, 2.5 ); // mb anti-quark-quark annihilation
257  Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
258  } else {
259  X_b = 6.8/SqrtS; // mb anti-quark-quark annihilation
260  Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
261  }
262 
263  X_c = 2.0*FlowF*sqr( ProjectileMass + TargetMass )/ECMSsqr; // mb rearrangement
264 
265  //G4cout << "Old new Xa " << 35.*FlowF << " " << 25.*FlowF << G4endl;
266 
267  X_d = 23.3/ECMSsqr; // mb anti-quark-quark string creation
268  }
269 
270  //G4cout << "Param Xtotal Xelastic " << Xtotal << " " << Xelastic << G4endl
271  // << "Para a b c d " << X_a << " " << X_b << " " << X_c << " " << X_d << G4endl;
272  // << "Para a b c d " << X_a << " " << 5.*X_b << " " << 5.*X_c << " " << 6.*X_d
273  // << G4endl;
274 
275  G4double Xann_on_P( 0.0), Xann_on_N( 0.0 );
276 
277  if ( ProjectilePDGcode == -2212 ) { // Pbar+P/N
278  Xann_on_P = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
279  Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
280  } else if ( ProjectilePDGcode == -2112 ) { // NeutrBar+P/N
281  Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
282  Xann_on_N = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
283  } else if ( ProjectilePDGcode == -3122 ) { // LambdaBar+P/N
284  Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
285  Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
286  } else if ( ProjectilePDGcode == -3112 ) { // Sigma-Bar+P/N
287  Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
288  Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
289  } else if ( ProjectilePDGcode == -3212 ) { // Sigma0Bar+P/N
290  Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
291  Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
292  } else if ( ProjectilePDGcode == -3222 ) { // Sigma+Bar+P/N
293  Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
294  Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
295  } else if ( ProjectilePDGcode == -3312 ) { // Xi-Bar+P/N
296  Xann_on_P = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
297  Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
298  } else if ( ProjectilePDGcode == -3322 ) { // Xi0Bar+P/N
299  Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
300  Xann_on_N = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
301  } else if ( ProjectilePDGcode == -3334 ) { // Omega-Bar+P/N
302  Xann_on_P = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
303  Xann_on_N = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
304  } else {
305  G4cout << "Unknown anti-baryon for FTF annihilation" << G4endl;
306  }
307 
308  //G4cout << "Sum " << Xann_on_P << G4endl;
309 
310  if ( ! ProjectileIsNucleus ) { // Projectile is anti-baryon
311  Xannihilation = ( NumberOfTargetProtons * Xann_on_P + NumberOfTargetNeutrons * Xann_on_N )
312  / NumberOfTargetNucleons;
313  } else { // Projectile is a nucleus
314  Xannihilation = (
315  ( AbsProjectileCharge * NumberOfTargetProtons +
316  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
317  NumberOfTargetNeutrons ) * Xann_on_P
318  +
319  ( AbsProjectileCharge * NumberOfTargetNeutrons +
320  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
321  NumberOfTargetProtons ) * Xann_on_N
322  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
323  }
324 
325  //G4double Xftf = 0.0;
326  MesonProdThreshold = ProjectileMass + TargetMass + (0.14 + 0.08); // Mpi + DeltaE
327  if ( SqrtS > MesonProdThreshold ) {
328  Xftf = 36.0 * ( 1.0 - MesonProdThreshold/SqrtS );
329  }
330 
331  Xtotal = Xelastic + Xannihilation + Xftf;
332 
333  #ifdef debugFTFparams
334  G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
335  << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << G4endl
336  << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal << " "
337  << Xannihilation/(Xtotal - Xelastic) << G4endl;
338  #endif
339 
340  } else if ( ProjectilePDGcode == 211 ) { // Projectile is PionPlus
341 
342  G4double XtotPiP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
344  G4double XtotPiN = FTFxsManager->GetTotalElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
345  G4double XelPiP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
346  G4double XelPiN = FTFxsManager->GetElasticElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
347  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
348  / NumberOfTargetNucleons;
349  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
350  / NumberOfTargetNucleons;
351  Xannihilation = 0.0;
352  Xtotal /= millibarn;
353  Xelastic /= millibarn;
354 
355  } else if ( ProjectilePDGcode == -211 ) { // Projectile is PionMinus
356 
357  G4double XtotPiP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
359  G4double XtotPiN = FTFxsManager->GetTotalElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
360  G4double XelPiP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
361  G4double XelPiN = FTFxsManager->GetElasticElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
362  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
363  / NumberOfTargetNucleons;
364  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
365  / NumberOfTargetNucleons;
366  Xannihilation = 0.0;
367  Xtotal /= millibarn;
368  Xelastic /= millibarn;
369 
370  } else if ( ProjectilePDGcode == 111 ) { // Projectile is PionZero
371 
373  G4double XtotPipP = FTFxsManager->GetTotalElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
375  G4double XtotPimP = FTFxsManager->GetTotalElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
376  G4double XelPipP = FTFxsManager->GetElasticElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
377  G4double XelPimP = FTFxsManager->GetElasticElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
378  G4double XtotPiP = ( XtotPipP + XtotPimP ) / 2.0;
379  G4double XtotPiN = XtotPiP;
380  G4double XelPiP = ( XelPipP + XelPimP ) / 2.0;
381  G4double XelPiN = XelPiP;
382  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
383  / NumberOfTargetNucleons;
384  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
385  / NumberOfTargetNucleons;
386  Xannihilation = 0.0;
387  Xtotal /= millibarn;
388  Xelastic /= millibarn;
389 
390  } else if ( ProjectilePDGcode == 321 ) { // Projectile is KaonPlus
391 
392  G4double XtotKP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
394  G4double XtotKN = FTFxsManager->GetTotalElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
395  G4double XelKP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
396  G4double XelKN = FTFxsManager->GetElasticElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
397  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
398  / NumberOfTargetNucleons;
399  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
400  / NumberOfTargetNucleons;
401  Xannihilation = 0.0;
402  Xtotal /= millibarn;
403  Xelastic /= millibarn;
404 
405  } else if ( ProjectilePDGcode == -321 ) { // Projectile is KaonMinus
406 
407  G4double XtotKP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
409  G4double XtotKN = FTFxsManager->GetTotalElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
410  G4double XelKP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
411  G4double XelKN = FTFxsManager->GetElasticElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
412  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
413  / NumberOfTargetNucleons;
414  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
415  / NumberOfTargetNucleons;
416  Xannihilation = 0.0;
417  Xtotal /= millibarn;
418  Xelastic /= millibarn;
419 
420  } else if ( ProjectilePDGcode == 311 || ProjectilePDGcode == 130 ||
421  ProjectilePDGcode == 310 ) { // Projectile is KaonZero
422 
424  G4double XtotKpP = FTFxsManager->GetTotalElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
426  G4double XtotKmP = FTFxsManager->GetTotalElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
427  G4double XelKpP = FTFxsManager->GetElasticElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
428  G4double XelKmP = FTFxsManager->GetElasticElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
429  G4double XtotKP = ( XtotKpP + XtotKmP ) / 2.0;
430  G4double XtotKN = XtotKP;
431  G4double XelKP = ( XelKpP + XelKmP ) / 2.0;
432  G4double XelKN = XelKP;
433  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
434  / NumberOfTargetNucleons;
435  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
436  / NumberOfTargetNucleons;
437  Xannihilation = 0.0;
438  Xtotal /= millibarn;
439  Xelastic /= millibarn;
440 
441  } else { // Projectile is undefined, Nucleon assumed
442 
444  G4double XtotPP = FTFxsManager->GetTotalElementCrossSection( Proton, KineticEnergy, 1, 0 );
446  G4double XtotPN = FTFxsManager->GetTotalElementCrossSection( Neutron, KineticEnergy, 1, 0 );
447  G4double XelPP = FTFxsManager->GetElasticElementCrossSection( Proton, KineticEnergy, 1, 0 );
448  G4double XelPN = FTFxsManager->GetElasticElementCrossSection( Neutron, KineticEnergy, 1, 0 );
449  Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN )
450  / NumberOfTargetNucleons;
451  Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN )
452  / NumberOfTargetNucleons;
453  Xannihilation = 0.0;
454  Xtotal /= millibarn;
455  Xelastic /= millibarn;
456 
457  };
458 
459  // Geometrical parameters
460  SetTotalCrossSection( Xtotal );
461  SetElastisCrossSection( Xelastic );
462  SetInelasticCrossSection( Xtotal - Xelastic );
463 
464  //G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
465  // << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << G4endl;
466  //if (Xtotal - Xelastic != 0.0 ) {
467  // G4cout << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal
468  // << " " << Xannihilation / (Xtotal - Xelastic) << G4endl;
469  //} else {
470  // G4cout << "Plab Xelastic/Xtotal, Xann " << Plab << " " << Xelastic/Xtotal
471  // << " " << Xannihilation << G4endl;
472  //}
473  //G4int Uzhi; G4cin >> Uzhi;
474 
475  // Interactions with elastic and inelastic collisions
476  SetProbabilityOfElasticScatt( Xtotal, Xelastic );
477  SetRadiusOfHNinteractions2( Xtotal/pi/10.0 );
478  if ( Xtotal - Xelastic == 0.0 ) {
480  } else {
481  SetProbabilityOfAnnihilation( Xannihilation / (Xtotal - Xelastic) );
482  }
483 
484  // No elastic scattering
485  //SetProbabilityOfElasticScatt( Xtotal, 0.0 );
486  //SetRadiusOfHNinteractions2( (Xtotal - Xelastic)/pi/10.0 );
487  //SetProbabilityOfAnnihilation( 1.0 );
488  //SetProbabilityOfAnnihilation( 0.0 );
489 
490  SetSlope( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 ); // Slope parameter of elastic scattering
491  // (GeV/c)^(-2))
492  //G4cout << "Slope " << GetSlope() << G4endl;
493  SetGamma0( GetSlope()*Xtotal/10.0/2.0/pi );
494 
495  // Parameters of elastic scattering
496  // Gaussian parametrization of elastic scattering amplitude assumed
497  SetAvaragePt2ofElasticScattering( 1.0/( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 )*GeV*GeV );
498  //G4cout << "AvaragePt2ofElasticScattering " << GetAvaragePt2ofElasticScattering() << G4endl;
499 
500  // Parameters of excitations
501 
502  G4double Xinel = Xtotal - Xelastic; // Uzhi 25.04.2012
503  //G4cout << "Param ProjectilePDGcode " << ProjectilePDGcode << G4endl;
504 
505  if ( ProjectilePDGcode > 1000 ) { // Projectile is baryon
506 
507  // Proc# A1 B1 A2 B2 A3 Atop Ymin
508  SetParams( 0, 13.71, 1.75, -214.4 , 4.25, 0.0, 0.632, 1.45 ); // Qexchange without Exc.
509  SetParams( 1, 0.2, 0.0 , - 3.289, 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
510  SetParams( 2, 6.0/Xinel, 0.0 , -6.0/Xinel*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Projectile diffraction
511  SetParams( 3, 6.0/Xinel, 0.0 , -6.0/Xinel*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Target diffraction
512 
513  SetParams( 1, 0.3, 0.5 , - 3.289, 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
514  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
515  SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
516  SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
517  }
519  if ( NumberOfTargetNucleons > 26 ) {
521  } else {
523  }
524  SetProjMinDiffMass( 1.16 ); // GeV
525  SetProjMinNonDiffMass( 1.16 ); // GeV
526  SetTarMinDiffMass( 1.16 ); // GeV
527  SetTarMinNonDiffMass( 1.16 ); // GeV
528  SetAveragePt2( 0.3 ); // GeV^2
529  SetProbLogDistr(0.5 );
530 
531  } else if( ProjectilePDGcode < -1000 ) { // Projectile is anti_baryon
532 
533  // Proc# A1 B1 A2 B2 A3 Atop Ymin
534  SetParams( 0, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange without Exc.
535  SetParams( 1, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange with Exc.
536  if ( Xftf > 0.0 ) {
537  SetParams( 2, 6.0/Xftf, 0.0 , -6.0/Xftf*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Projectile diffraction
538  SetParams( 3, 6.0/Xftf, 0.0 , -6.0/Xftf*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Target diffraction
539  } else {
540  SetParams( 2, 0.5 , 0.0 , 0.0 , 0.0 , 0.0, 0.5 , 1000.0 ); // Projectile diffraction
541  SetParams( 3, 0.5 , 0.0 , 0.0 , 0.0 , 0.0, 0.5 , 1000.0 ); // Target diffraction
542  }
543  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
544  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
545  SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
546  }
549  SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
550  SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
551  SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
552  SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
553  SetAveragePt2( 0.3 ); // 0.15 GeV^2 // Uzhi 21.05.2012
554  SetProbLogDistr( 0.5 ); // Uzhi 21.05.2012
555 
556  } else if ( ProjectileabsPDGcode == 211 || ProjectilePDGcode == 111 ) { // Projectile is Pion
557 
558  // Proc# A1 B1 A2 B2 A3 Atop Ymin
559  SetParams( 0, 568.0 , 2.1 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
560  SetParams( 1, 6.0 , 0.6 , -26.9 , 1.1 , 0.0, 0.0 , 3.0 ); // Qexchange with Exc.
561  G4double Wprd = 0.0;
562  if ( Xinel > 0.0 ) Wprd = 0.64 *( 6.2 - 3.7*std::exp( - sqr( SqrtS - 7.0 ) / 16.0 ) ) / Xinel;
563  SetParams( 2, Wprd, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
564  G4double Wtrd = 0.0;
565  if ( Xinel > 0.0 ) Wtrd = ( 2.0 + 22.0/ECMSsqr ) / Xinel;
566  SetParams( 3, Wtrd, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
567  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
568  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
569  SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
570  }
571  SetDeltaProbAtQuarkExchange( 0.56 ); // (0.35)
572  SetProjMinDiffMass( 0.5 ); // (0.5) // GeV
573  SetProjMinNonDiffMass( 0.5 ); // (0.5) // GeV
574  SetTarMinDiffMass( 1.16 ); // GeV
575  SetTarMinNonDiffMass( 1.16 ); // GeV
576  SetAveragePt2( 0.3 ); // GeV^2
577  SetProbLogDistr( 1.0 ); // (0.0) // Uzhi 21.05.2012
578 
579  } else if ( ProjectileabsPDGcode == 321 || ProjectileabsPDGcode == 311 ||
580  ProjectilePDGcode == 130 || ProjectilePDGcode == 310 ) { // Projectile is Kaon
581 
582  // Proc# A1 B1 A2 B2 A3 Atop Ymin
583  SetParams( 0, 70.0 , 2.75, 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
584  SetParams( 1, 30.0 , 1.5 , -50.57 , 1.83, 0.0, 0.0 , 1.70 ); // Qexchange with Exc.
585  SetParams( 2, 0.6 , 0.75, -12.05 , 3.25, 0.0, 0.0 , 1.20 ); // Projectile diffraction
586  SetParams( 3, 6.0 , 1.0 , -12.08 , 1.5 , 0.1, 0.0 , 1.20 ); // Target diffraction
587  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
588  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
589  SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
590  }
592  SetProjMinDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
593  SetProjMinNonDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
594  SetTarMinDiffMass( 1.16 ); // GeV
595  SetTarMinNonDiffMass( 1.16 ); // GeV
596  SetAveragePt2( 0.3 ); // GeV^2 7 June 2011
597  SetProbLogDistr( 1.0 ); // Uzhi 5.06.2012
598 
599  } else { // Projectile is undefined, Nucleon assumed
600 
601  // Proc# A1 B1 A2 B2 A3 Atop Ymin
602  SetParams( 0, 13.71, 1.75, -214.4 , 4.25, 0.0, 0.632, 1.45 ); // Qexchange without Exc.
603  SetParams( 1, 0.2 , 0.0 , -16.445, 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
604  SetParams( 2, 6.0/Xinel, 0.0 , -6.0/Xinel*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Projectile diffraction
605  SetParams( 3, 6.0/Xinel, 0.0 , -6.0/Xinel*8.48, 2.25, 0.0, 0.0 , 0.95 ); // Target diffraction
606  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
607  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
608  SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
609  }
610  SetDeltaProbAtQuarkExchange( 0.0 ); // 7 June 2011
612  SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
613  SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
614  SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
615  SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
616  SetAveragePt2( 0.3 ); // (0.15) GeV^2 Uzhi 21.05.2012
617  SetProbLogDistr( 0.5 ); // Uzhi 21.05.2012
618 
619  }
620 
621  //if ( theA > 4 ) SetProbabilityOfProjDiff( 0.0 ); // Uzhi 6.07.2012 Closed
622  //G4cout << "Param Get Min Dif " << GetProjMinNonDiffMass() << G4endl;
623 
624  // Set parameters of a string kink
625  SetPt2Kink( 6.0*GeV*GeV );
626  G4double Puubar( 1.0/3.0 ), Pddbar( 1.0/3.0 ), Pssbar( 1.0/3.0 ); // SU(3) symmetry
627  //G4double Puubar( 0.41 ), Pddbar( 0.41 ), Pssbar( 0.18 ); // Broken SU(3) symmetry
628  SetQuarkProbabilitiesAtGluonSplitUp( Puubar, Pddbar, Pssbar );
629 
630  // Set parameters of nuclear destruction
631  if ( ProjectileabsPDGcode < 1000 ) { // Meson projectile
632  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
633  SetCofNuclearDestruction( 1.0*std::exp( 4.0*(Ylab - 2.1) )/
634  ( 1.0 + std::exp( 4.0*(Ylab - 2.1) ) ) ); // 0.62 1.0
637  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*std::exp( 4.0*(Ylab - 2.5) )/
638  ( 1.0 + std::exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV ); // 0.09
641  } else if ( ProjectilePDGcode < -1000 ) { // for anti-baryon projectile
642  //G4cout << "Nucl destruct Anti Bar" << G4endl;
643  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
644  SetCofNuclearDestruction( 1.0*std::exp( 4.0*(Ylab - 2.1) )/
645  ( 1.0 + std::exp( 4.0*(Ylab - 2.1) ) ) ); // 0.62 1.0
648  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*std::exp( 4.0*(Ylab - 2.5) )/
649  ( 1.0 + std::exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV ); // 0.09
652  if ( Plab < 2.0 ) { // 2 GeV/c
653  // For slow anti-baryon we have to garanty putting on mass-shell
656  SetDofNuclearDestruction( 0.01 );
659  //SetExcitationEnergyPerWoundedNucleon( 0.0 ); // ?????
660  }
661  } else { // Projectile baryon assumed
662  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
663  SetCofNuclearDestruction( 1.0*std::exp( 4.0*(Ylab - 2.1) )/
664  ( 1.0 + std::exp( 4.0*(Ylab - 2.1) ) ) ); // 0.62 1.0
667  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*std::exp( 4.0*(Ylab - 2.5) )/
668  ( 1.0 + std::exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV ); // 0.09
671  }
672 
673  //SetCofNuclearDestruction( 0.47*std::exp( 2.0*(Ylab - 2.5) )/( 1.0 + std::exp( 2.0*(Ylab - 2.5) ) ) );
674  //SetPt2ofNuclearDestruction( ( 0.035 + 0.1*std::exp( 4.0*(Ylab - 3.0) )/( 1.0 + std::exp( 4.0*(Ylab - 3.0) ) ) )*GeV*GeV );
675 
676  //SetMagQuarkExchange( 120.0 ); // 210.0 PipP
677  //SetSlopeQuarkExchange( 2.0 );
678  //SetDeltaProbAtQuarkExchange( 0.6 );
679  //SetProjMinDiffMass( 0.7 ); // GeV 1.1
680  //SetProjMinNonDiffMass( 0.7 ); // GeV
681  //SetProbabilityOfProjDiff( 0.0); // 0.85*std::pow( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
682  //SetTarMinDiffMass( 1.1 ); // GeV
683  //SetTarMinNonDiffMass( 1.1 ); // GeV
684  //SetProbabilityOfTarDiff( 0.0 ); // 0.85*std::pow( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
685 
686  //SetAveragePt2( 0.3 ); // GeV^2
687  //------------------------------------
688  //SetProbabilityOfElasticScatt( 1.0, 1.0); //(Xtotal, Xelastic);
689  //SetProbabilityOfProjDiff( 1.0*0.62*std::pow( s/GeV/GeV, -0.51 ) ); // 0->1
690  //SetProbabilityOfTarDiff( 4.0*0.62*std::pow( s/GeV/GeV, -0.51 ) ); // 2->4
691  //SetAveragePt2( 0.3 ); // (0.15)
692  //SetAvaragePt2ofElasticScattering( 0.0 );
693 
694  //SetMaxNumberOfCollisions( Plab, 6.0 ); //(4.0*(Plab + 0.01), Plab); // 6.0 );
695  //SetAveragePt2( 0.15 );
696  //G4cout << "Cnd " << GetCofNuclearDestruction() << G4endl;
697  //SetCofNuclearDestruction( 0.4 ); // (0.2) // (0.4)
698  //SetExcitationEnergyPerWoundedNucleon( 0.0*MeV ); // (75.0*MeV)
699  //SetDofNuclearDestruction( 0.0 );
700  //SetPt2ofNuclearDestruction( 0.0*GeV*GeV ); // (0.168*GeV*GeV)
701  //G4cout << "Pt2 " << GetPt2ofNuclearDestruction()/GeV/GeV << G4endl;
702  //G4int Uzhi; G4cin >> Uzhi;
703 
704 }
void SetProbLogDistr(const G4double aValue)
void SetTarMinNonDiffMass(const G4double aValue)
void SetDofNuclearDestruction(const G4double aValue)
void SetSlope(const G4double Slope)
G4ChipsComponentXS * FTFxsManager
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
void SetGamma0(const G4double Gamma0)
G4GLOB_DLL std::ostream G4cout
void SetAveragePt2(const G4double aValue)
void SetElastisCrossSection(const G4double Xelastic)
bool G4bool
Definition: G4Types.hh:79
void SetPt2Kink(const G4double aValue)
void SetRadiusOfHNinteractions2(const G4double Radius2)
void SetMaxNumberOfCollisions(const G4double aValue, const G4double bValue)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
static G4Proton * Proton()
Definition: G4Proton.cc:93
void SetProjMinNonDiffMass(const G4double aValue)
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
void SetQuarkProbabilitiesAtGluonSplitUp(const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double N)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
void SetInelasticCrossSection(const G4double Xinelastic)
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void SetDeltaProbAtQuarkExchange(const G4double aValue)
void SetTotalCrossSection(const G4double Xtotal)
void SetProjMinDiffMass(const G4double aValue)
G4double ProbOfSameQuarkExchange
G4double GetPDGMass() const
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
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)
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
void SetCofNuclearDestruction(const G4double aValue)
void SetProbOfSameQuarkExchange(const G4double aValue)
G4double FTFXannihilation
#define G4endl
Definition: G4ios.hh:61
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double N)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
void SetProbabilityOfAnnihilation(const G4double aValue)
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double GetPDGCharge() const
void SetPt2ofNuclearDestruction(const G4double aValue)
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
void SetTarMinDiffMass(const G4double aValue)
void SetR2ofNuclearDestruction(const G4double aValue)
G4FTFParameters::~G4FTFParameters ( )

Definition at line 83 of file G4FTFParameters.cc.

83 {}
G4FTFParameters::G4FTFParameters ( )

Definition at line 57 of file G4FTFParameters.cc.

References ProcParams.

57  :
58  FTFhNcmsEnergy( 0.0 ),
59  FTFxsManager( 0 ),
60  FTFXtotal( 0.0 ), FTFXelastic( 0.0 ), FTFXinelastic( 0.0 ), FTFXannihilation( 0.0 ),
62  RadiusOfHNinteractions2( 0.0 ), FTFSlope( 0.0 ),
65  ProjMinDiffMass( 0.0 ), ProjMinNonDiffMass( 0.0 ),
66  TarMinDiffMass( 0.0 ), TarMinNonDiffMass( 0.0 ),
67  AveragePt2( 0.0 ), ProbLogDistr( 0.0 ),
68  Pt2kink( 0.0 ),
72 {
73  for ( G4int i = 0; i < 4; i++ ) {
74  for ( G4int j = 0; j < 7; j++ ) {
75  ProcParams[i][j] = 0.0;
76  }
77  }
78 }
G4double ProjMinDiffMass
G4double ProbabilityOfAnnihilation
G4double Pt2ofNuclearDestruction
G4ChipsComponentXS * FTFxsManager
int G4int
Definition: G4Types.hh:78
G4double CofNuclearDestruction
G4double ExcitationEnergyPerWoundedNucleon
G4double DeltaProbAtQuarkExchange
G4double DofNuclearDestruction
G4double ProcParams[4][7]
G4double AvaragePt2ofElasticScattering
G4double ProbabilityOfElasticScatt
G4double RadiusOfHNinteractions2
G4double MaxPt2ofNuclearDestruction
G4double ProbOfSameQuarkExchange
G4double R2ofNuclearDestruction
G4double ProbOfInelInteraction
G4double FTFXannihilation
G4double ProjMinNonDiffMass
G4double MaxNumberOfCollisions
G4double TarMinNonDiffMass

Member Function Documentation

G4double G4FTFParameters::GammaElastic ( const G4double  impactsquare)
inline

Definition at line 207 of file G4FTFParameters.hh.

References FTFGamma0, and FTFSlope.

Referenced by GetInelasticProbability().

207  {
208  return ( FTFGamma0 * std::exp( -FTFSlope * impactsquare ) );
209 }
G4double G4FTFParameters::GetAvaragePt2ofElasticScattering ( )
inline

Definition at line 402 of file G4FTFParameters.hh.

References AvaragePt2ofElasticScattering.

Referenced by G4ElasticHNScattering::ElasticScattering().

402  {
404 }
G4double AvaragePt2ofElasticScattering
G4double G4FTFParameters::GetAveragePt2 ( )
inline

Definition at line 432 of file G4FTFParameters.hh.

References AveragePt2.

Referenced by G4FTFAnnihilation::Annihilate(), and G4DiffractiveExcitation::ExciteParticipants().

432  {
433  return AveragePt2;
434 }
G4double G4FTFParameters::GetCofNuclearDestruction ( )
inline

Definition at line 460 of file G4FTFParameters.hh.

References CofNuclearDestruction.

460  {
461  return CofNuclearDestruction;
462 }
G4double CofNuclearDestruction
G4double G4FTFParameters::GetDeltaProbAtQuarkExchange ( )
inline

Definition at line 408 of file G4FTFParameters.hh.

References DeltaProbAtQuarkExchange.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

408  {
410 }
G4double DeltaProbAtQuarkExchange
G4double G4FTFParameters::GetDofNuclearDestruction ( )
inline

Definition at line 472 of file G4FTFParameters.hh.

References DofNuclearDestruction.

472  {
473  return DofNuclearDestruction;
474 }
G4double DofNuclearDestruction
G4double G4FTFParameters::GetElasticCrossSection ( )
inline

Definition at line 368 of file G4FTFParameters.hh.

References FTFXelastic.

368  {
369  return FTFXelastic;
370 }
G4double G4FTFParameters::GetExcitationEnergyPerWoundedNucleon ( )
inline

Definition at line 468 of file G4FTFParameters.hh.

References ExcitationEnergyPerWoundedNucleon.

468  {
470 }
G4double ExcitationEnergyPerWoundedNucleon
G4double G4FTFParameters::GetInelasticCrossSection ( )
inline

Definition at line 372 of file G4FTFParameters.hh.

References FTFXinelastic.

372  {
373  return FTFXinelastic;
374 }
G4double G4FTFParameters::GetInelasticProbability ( const G4double  impactsquare)
inline

Definition at line 392 of file G4FTFParameters.hh.

References GammaElastic().

392  {
393  G4double Gamma = GammaElastic( impactsquare );
394  return 2*Gamma - Gamma*Gamma;
395 }
double G4double
Definition: G4Types.hh:76
G4double GammaElastic(const G4double impactsquare)
G4double G4FTFParameters::GetMaxNumberOfCollisions ( )
inline

Definition at line 452 of file G4FTFParameters.hh.

References MaxNumberOfCollisions.

452  {
453  return MaxNumberOfCollisions;
454 }
G4double MaxNumberOfCollisions
G4double G4FTFParameters::GetMaxPt2ofNuclearDestruction ( )
inline

Definition at line 480 of file G4FTFParameters.hh.

References MaxPt2ofNuclearDestruction.

480  {
482 }
G4double MaxPt2ofNuclearDestruction
G4double G4FTFParameters::GetProbabilityOfAnnihilation ( )
inline

Definition at line 397 of file G4FTFParameters.hh.

References ProbabilityOfAnnihilation.

397  {
399 }
G4double ProbabilityOfAnnihilation
G4double G4FTFParameters::GetProbabilityOfElasticScatt ( )
inline

Definition at line 388 of file G4FTFParameters.hh.

References ProbabilityOfElasticScatt.

388  {
390 }
G4double ProbabilityOfElasticScatt
G4double G4FTFParameters::GetProbabilityOfInteraction ( const G4double  impactsquare)
inline

Definition at line 380 of file G4FTFParameters.hh.

References RadiusOfHNinteractions2.

Referenced by G4FTFParticipants::GetList().

380  {
381  if ( RadiusOfHNinteractions2 > impactsquare ) {
382  return 1.0;
383  } else {
384  return 0.0;
385  }
386 }
G4double RadiusOfHNinteractions2
G4double G4FTFParameters::GetProbLogDistr ( )
inline

Definition at line 436 of file G4FTFParameters.hh.

References ProbLogDistr.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

436  {
437  return ProbLogDistr;
438 }
G4double G4FTFParameters::GetProbOfInteraction ( )
inline

Definition at line 456 of file G4FTFParameters.hh.

References ProbOfInelInteraction.

456  {
457  return ProbOfInelInteraction;
458 }
G4double ProbOfInelInteraction
G4double G4FTFParameters::GetProbOfSameQuarkExchange ( )
inline

Definition at line 412 of file G4FTFParameters.hh.

References ProbOfSameQuarkExchange.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

412  {
414 }
G4double ProbOfSameQuarkExchange
G4double G4FTFParameters::GetProcProb ( const G4int  ProcN,
const G4double  y 
)

Definition at line 709 of file G4FTFParameters.cc.

References ProcParams.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

709  {
710  G4double Prob( 0.0 );
711  if ( y < ProcParams[ProcN][6] ) {
712  Prob = ProcParams[ProcN][5];
713  return Prob;
714  }
715  Prob = ProcParams[ProcN][0] * std::exp( -ProcParams[ProcN][1]*y ) +
716  ProcParams[ProcN][2] * std::exp( -ProcParams[ProcN][3]*y ) +
717  ProcParams[ProcN][4];
718  return Prob;
719 }
G4double ProcParams[4][7]
double G4double
Definition: G4Types.hh:76
G4double G4FTFParameters::GetProjMinDiffMass ( )
inline

Definition at line 416 of file G4FTFParameters.hh.

References ProjMinDiffMass.

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

416  {
417  return ProjMinDiffMass;
418 }
G4double ProjMinDiffMass
G4double G4FTFParameters::GetProjMinNonDiffMass ( )
inline

Definition at line 420 of file G4FTFParameters.hh.

References ProjMinNonDiffMass.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

420  {
421  return ProjMinNonDiffMass;
422 }
G4double ProjMinNonDiffMass
G4double G4FTFParameters::GetPt2Kink ( )
inline

Definition at line 442 of file G4FTFParameters.hh.

References Pt2kink.

Referenced by G4DiffractiveExcitation::CreateStrings().

442  {
443  return Pt2kink;
444 }
G4double G4FTFParameters::GetPt2ofNuclearDestruction ( )
inline

Definition at line 476 of file G4FTFParameters.hh.

References Pt2ofNuclearDestruction.

476  {
478 }
G4double Pt2ofNuclearDestruction
std::vector< G4double > G4FTFParameters::GetQuarkProbabilitiesAtGluonSplitUp ( )
inline

Definition at line 446 of file G4FTFParameters.hh.

References QuarkProbabilitiesAtGluonSplitUp.

Referenced by G4DiffractiveExcitation::CreateStrings().

446  {
448 }
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp
G4double G4FTFParameters::GetR2ofNuclearDestruction ( )
inline

Definition at line 464 of file G4FTFParameters.hh.

References R2ofNuclearDestruction.

464  {
465  return R2ofNuclearDestruction;
466 }
G4double R2ofNuclearDestruction
G4double G4FTFParameters::GetSlope ( )
inline

Definition at line 376 of file G4FTFParameters.hh.

References FTFSlope.

Referenced by G4FTFParameters().

376  {
377  return FTFSlope;
378 }
G4double G4FTFParameters::GetTarMinDiffMass ( )
inline

Definition at line 424 of file G4FTFParameters.hh.

References TarMinDiffMass.

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

424  {
425  return TarMinDiffMass;
426 }
G4double G4FTFParameters::GetTarMinNonDiffMass ( )
inline

Definition at line 428 of file G4FTFParameters.hh.

References TarMinNonDiffMass.

Referenced by G4DiffractiveExcitation::ExciteParticipants().

428  {
429  return TarMinNonDiffMass;
430 }
G4double TarMinNonDiffMass
G4double G4FTFParameters::GetTotalCrossSection ( )
inline

Definition at line 364 of file G4FTFParameters.hh.

References FTFXtotal.

364  {
365  return FTFXtotal;
366 }
void G4FTFParameters::SetAvaragePt2ofElasticScattering ( const G4double  aPt2)
inline

Definition at line 259 of file G4FTFParameters.hh.

References AvaragePt2ofElasticScattering.

Referenced by G4FTFParameters().

259  {
261 }
G4double AvaragePt2ofElasticScattering
void G4FTFParameters::SetAveragePt2 ( const G4double  aValue)
inline

Definition at line 299 of file G4FTFParameters.hh.

References AveragePt2.

Referenced by G4FTFParameters().

299  {
300  AveragePt2 = aValue*CLHEP::GeV*CLHEP::GeV;
301 }
void G4FTFParameters::SetCofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 339 of file G4FTFParameters.hh.

References CofNuclearDestruction.

Referenced by G4FTFParameters().

339  {
340  CofNuclearDestruction = aValue;
341 }
G4double CofNuclearDestruction
void G4FTFParameters::SetDeltaProbAtQuarkExchange ( const G4double  aValue)
inline

Definition at line 275 of file G4FTFParameters.hh.

References DeltaProbAtQuarkExchange.

Referenced by G4FTFParameters().

275  {
276  DeltaProbAtQuarkExchange = aValue;
277 }
G4double DeltaProbAtQuarkExchange
void G4FTFParameters::SetDofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 351 of file G4FTFParameters.hh.

References DofNuclearDestruction.

Referenced by G4FTFParameters().

351  {
352  DofNuclearDestruction = aValue;
353 }
G4double DofNuclearDestruction
void G4FTFParameters::SetElastisCrossSection ( const G4double  Xelastic)
inline

Definition at line 221 of file G4FTFParameters.hh.

References FTFXelastic.

Referenced by G4FTFParameters().

221  {
222  FTFXelastic = Xelastic;
223 }
void G4FTFParameters::SetExcitationEnergyPerWoundedNucleon ( const G4double  aValue)
inline

Definition at line 347 of file G4FTFParameters.hh.

References ExcitationEnergyPerWoundedNucleon.

Referenced by G4FTFParameters().

347  {
349 }
G4double ExcitationEnergyPerWoundedNucleon
void G4FTFParameters::SetGamma0 ( const G4double  Gamma0)
inline

Definition at line 254 of file G4FTFParameters.hh.

References FTFGamma0.

Referenced by G4FTFParameters().

254  {
255  FTFGamma0 = Gamma0;
256 }
void G4FTFParameters::SethNcmsEnergy ( const G4double  s)
inline

Definition at line 211 of file G4FTFParameters.hh.

References FTFhNcmsEnergy.

211  {
212  FTFhNcmsEnergy = S;
213 }
void G4FTFParameters::SetInelasticCrossSection ( const G4double  Xinelastic)
inline

Definition at line 225 of file G4FTFParameters.hh.

References FTFXinelastic.

Referenced by G4FTFParameters().

225  {
226  FTFXinelastic = Xinelastic;
227 }
void G4FTFParameters::SetMaxNumberOfCollisions ( const G4double  aValue,
const G4double  bValue 
)
inline

Definition at line 322 of file G4FTFParameters.hh.

References MaxNumberOfCollisions, and SetProbOfInteraction().

Referenced by G4FTFParameters().

323  {
324  if ( Plab > Pbound ) {
325  MaxNumberOfCollisions = Plab/Pbound;
326  SetProbOfInteraction( -1.0 );
327  } else {
328  //MaxNumberOfCollisions = -1.0;
329  //SetProbOfInteraction( std::exp( 0.25*(Plab-Pbound) ) );
331  SetProbOfInteraction( -1.0 );
332  }
333 }
void SetProbOfInteraction(const G4double aValue)
G4double MaxNumberOfCollisions
void G4FTFParameters::SetMaxPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 359 of file G4FTFParameters.hh.

References MaxPt2ofNuclearDestruction.

Referenced by G4FTFParameters().

359  {
361 }
G4double MaxPt2ofNuclearDestruction
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 265 of file G4FTFParameters.hh.

References ProcParams.

Referenced by G4FTFParameters().

268  {
269  ProcParams[ProcN][0] = A1; ProcParams[ProcN][1] = B1;
270  ProcParams[ProcN][2] = A2; ProcParams[ProcN][3] = B2;
271  ProcParams[ProcN][4] = A3;
272  ProcParams[ProcN][5] = Atop; ProcParams[ProcN][6] = Ymin;
273 }
G4double ProcParams[4][7]
void G4FTFParameters::SetProbabilityOfAnnihilation ( const G4double  aValue)
inline

Definition at line 242 of file G4FTFParameters.hh.

References ProbabilityOfAnnihilation.

Referenced by G4FTFParameters().

242  {
243  ProbabilityOfAnnihilation = aValue;
244 }
G4double ProbabilityOfAnnihilation
void G4FTFParameters::SetProbabilityOfElasticScatt ( const G4double  Xtotal,
const G4double  Xelastic 
)
inline

Definition at line 229 of file G4FTFParameters.hh.

References ProbabilityOfElasticScatt.

Referenced by G4FTFParameters().

230  {
231  if ( Xtotal == 0.0 ) {
233  } else {
234  ProbabilityOfElasticScatt = Xelastic / Xtotal;
235  }
236 }
G4double ProbabilityOfElasticScatt
void G4FTFParameters::SetProbabilityOfElasticScatt ( const G4double  aValue)
inline

Definition at line 238 of file G4FTFParameters.hh.

References ProbabilityOfElasticScatt.

238  {
239  ProbabilityOfElasticScatt = aValue;
240 }
G4double ProbabilityOfElasticScatt
void G4FTFParameters::SetProbabilityOfProjDiff ( const G4double  aValue)
void G4FTFParameters::SetProbabilityOfTarDiff ( const G4double  aValue)
void G4FTFParameters::SetProbLogDistr ( const G4double  aValue)
inline

Definition at line 303 of file G4FTFParameters.hh.

References ProbLogDistr.

Referenced by G4FTFParameters().

303  {
304  ProbLogDistr = aValue;
305 }
void G4FTFParameters::SetProbOfInteraction ( const G4double  aValue)
inline

Definition at line 335 of file G4FTFParameters.hh.

References ProbOfInelInteraction.

Referenced by SetMaxNumberOfCollisions().

335  {
336  ProbOfInelInteraction = aValue;
337 }
G4double ProbOfInelInteraction
void G4FTFParameters::SetProbOfSameQuarkExchange ( const G4double  aValue)
inline

Definition at line 279 of file G4FTFParameters.hh.

References ProbOfSameQuarkExchange.

Referenced by G4FTFParameters().

279  {
280  ProbOfSameQuarkExchange = aValue;
281 }
G4double ProbOfSameQuarkExchange
void G4FTFParameters::SetProjMinDiffMass ( const G4double  aValue)
inline

Definition at line 283 of file G4FTFParameters.hh.

References ProjMinDiffMass.

Referenced by G4FTFParameters().

283  {
284  ProjMinDiffMass = aValue*CLHEP::GeV;
285 }
G4double ProjMinDiffMass
void G4FTFParameters::SetProjMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 287 of file G4FTFParameters.hh.

References ProjMinNonDiffMass.

Referenced by G4FTFParameters().

287  {
288  ProjMinNonDiffMass = aValue*CLHEP::GeV;
289 }
G4double ProjMinNonDiffMass
void G4FTFParameters::SetPt2Kink ( const G4double  aValue)
inline

Definition at line 309 of file G4FTFParameters.hh.

References Pt2kink.

Referenced by G4FTFParameters().

309  {
310  Pt2kink = aValue;
311 }
void G4FTFParameters::SetPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 355 of file G4FTFParameters.hh.

References Pt2ofNuclearDestruction.

Referenced by G4FTFParameters().

355  {
356  Pt2ofNuclearDestruction = aValue;
357 }
G4double Pt2ofNuclearDestruction
void G4FTFParameters::SetQuarkProbabilitiesAtGluonSplitUp ( const G4double  Puubar,
const G4double  Pddbar,
const G4double  Pssbar 
)
inline

Definition at line 313 of file G4FTFParameters.hh.

References QuarkProbabilitiesAtGluonSplitUp.

Referenced by G4FTFParameters().

315  {
316  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar );
317  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar );
318  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar + Pssbar );
319 }
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp
void G4FTFParameters::SetR2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 343 of file G4FTFParameters.hh.

References R2ofNuclearDestruction.

Referenced by G4FTFParameters().

343  {
344  R2ofNuclearDestruction = aValue;
345 }
G4double R2ofNuclearDestruction
void G4FTFParameters::SetRadiusOfHNinteractions2 ( const G4double  Radius2)
inline

Definition at line 246 of file G4FTFParameters.hh.

References RadiusOfHNinteractions2.

Referenced by G4FTFParameters().

246  {
247  RadiusOfHNinteractions2 = Radius2;
248 }
G4double RadiusOfHNinteractions2
void G4FTFParameters::SetSlope ( const G4double  Slope)
inline

Definition at line 250 of file G4FTFParameters.hh.

References FTFSlope.

Referenced by G4FTFParameters().

250  {
251  FTFSlope = 12.84 / Slope; // Slope is in GeV^-2, FTFSlope in fm^-2
252 }
void G4FTFParameters::SetTarMinDiffMass ( const G4double  aValue)
inline

Definition at line 291 of file G4FTFParameters.hh.

References TarMinDiffMass.

Referenced by G4FTFParameters().

291  {
292  TarMinDiffMass = aValue*CLHEP::GeV;
293 }
void G4FTFParameters::SetTarMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 295 of file G4FTFParameters.hh.

References TarMinNonDiffMass.

Referenced by G4FTFParameters().

295  {
296  TarMinNonDiffMass = aValue*CLHEP::GeV;
297 }
G4double TarMinNonDiffMass
void G4FTFParameters::SetTotalCrossSection ( const G4double  Xtotal)
inline

Definition at line 217 of file G4FTFParameters.hh.

References FTFXtotal.

Referenced by G4FTFParameters().

217  {
218  FTFXtotal = Xtotal;
219 }

Field Documentation

G4double G4FTFParameters::AvaragePt2ofElasticScattering
G4double G4FTFParameters::AveragePt2

Definition at line 179 of file G4FTFParameters.hh.

Referenced by GetAveragePt2(), and SetAveragePt2().

G4double G4FTFParameters::CofNuclearDestruction

Definition at line 190 of file G4FTFParameters.hh.

Referenced by GetCofNuclearDestruction(), and SetCofNuclearDestruction().

G4double G4FTFParameters::DeltaProbAtQuarkExchange

Definition at line 170 of file G4FTFParameters.hh.

Referenced by GetDeltaProbAtQuarkExchange(), and SetDeltaProbAtQuarkExchange().

G4double G4FTFParameters::DofNuclearDestruction

Definition at line 195 of file G4FTFParameters.hh.

Referenced by GetDofNuclearDestruction(), and SetDofNuclearDestruction().

G4double G4FTFParameters::ExcitationEnergyPerWoundedNucleon
G4double G4FTFParameters::FTFGamma0

Definition at line 165 of file G4FTFParameters.hh.

Referenced by GammaElastic(), and SetGamma0().

G4double G4FTFParameters::FTFhNcmsEnergy

Definition at line 150 of file G4FTFParameters.hh.

Referenced by G4FTFParameters(), and SethNcmsEnergy().

G4double G4FTFParameters::FTFSlope

Definition at line 163 of file G4FTFParameters.hh.

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

G4double G4FTFParameters::FTFXannihilation

Definition at line 159 of file G4FTFParameters.hh.

Referenced by G4FTFParameters().

G4double G4FTFParameters::FTFXelastic

Definition at line 157 of file G4FTFParameters.hh.

Referenced by GetElasticCrossSection(), and SetElastisCrossSection().

G4double G4FTFParameters::FTFXinelastic

Definition at line 158 of file G4FTFParameters.hh.

Referenced by GetInelasticCrossSection(), and SetInelasticCrossSection().

G4ChipsComponentXS* G4FTFParameters::FTFxsManager

Definition at line 153 of file G4FTFParameters.hh.

Referenced by G4FTFParameters().

G4double G4FTFParameters::FTFXtotal

Definition at line 156 of file G4FTFParameters.hh.

Referenced by GetTotalCrossSection(), and SetTotalCrossSection().

G4double G4FTFParameters::MaxNumberOfCollisions

Definition at line 187 of file G4FTFParameters.hh.

Referenced by GetMaxNumberOfCollisions(), and SetMaxNumberOfCollisions().

G4double G4FTFParameters::MaxPt2ofNuclearDestruction
G4double G4FTFParameters::ProbabilityOfAnnihilation
G4double G4FTFParameters::ProbabilityOfElasticScatt
G4double G4FTFParameters::ProbLogDistr

Definition at line 180 of file G4FTFParameters.hh.

Referenced by GetProbLogDistr(), and SetProbLogDistr().

G4double G4FTFParameters::ProbOfInelInteraction

Definition at line 188 of file G4FTFParameters.hh.

Referenced by GetProbOfInteraction(), and SetProbOfInteraction().

G4double G4FTFParameters::ProbOfSameQuarkExchange
G4double G4FTFParameters::ProcParams[4][7]

Definition at line 168 of file G4FTFParameters.hh.

Referenced by G4FTFParameters(), GetProcProb(), and SetParams().

G4double G4FTFParameters::ProjMinDiffMass

Definition at line 173 of file G4FTFParameters.hh.

Referenced by GetProjMinDiffMass(), and SetProjMinDiffMass().

G4double G4FTFParameters::ProjMinNonDiffMass

Definition at line 174 of file G4FTFParameters.hh.

Referenced by GetProjMinNonDiffMass(), and SetProjMinNonDiffMass().

G4double G4FTFParameters::Pt2kink

Definition at line 183 of file G4FTFParameters.hh.

Referenced by GetPt2Kink(), and SetPt2Kink().

G4double G4FTFParameters::Pt2ofNuclearDestruction

Definition at line 196 of file G4FTFParameters.hh.

Referenced by GetPt2ofNuclearDestruction(), and SetPt2ofNuclearDestruction().

std::vector< G4double > G4FTFParameters::QuarkProbabilitiesAtGluonSplitUp
G4double G4FTFParameters::R2ofNuclearDestruction

Definition at line 191 of file G4FTFParameters.hh.

Referenced by GetR2ofNuclearDestruction(), and SetR2ofNuclearDestruction().

G4double G4FTFParameters::RadiusOfHNinteractions2

Definition at line 162 of file G4FTFParameters.hh.

Referenced by GetProbabilityOfInteraction(), and SetRadiusOfHNinteractions2().

G4double G4FTFParameters::TarMinDiffMass

Definition at line 176 of file G4FTFParameters.hh.

Referenced by GetTarMinDiffMass(), and SetTarMinDiffMass().

G4double G4FTFParameters::TarMinNonDiffMass

Definition at line 177 of file G4FTFParameters.hh.

Referenced by GetTarMinNonDiffMass(), and SetTarMinNonDiffMass().


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