Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4FTFParameters Class Reference

#include <G4FTFParameters.hh>

Collaboration diagram for G4FTFParameters:

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 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 ()
 
 G4FTFParameters ()
 

Public Attributes

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 [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 42 of file G4FTFParameters.hh.

Constructor & Destructor Documentation

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

Definition at line 98 of file G4FTFParameters.cc.

99  :
100  FTFhNcmsEnergy( 0.0 ),
101  FTFxsManager( 0 ),
102  FTFXtotal( 0.0 ), FTFXelastic( 0.0 ), FTFXinelastic( 0.0 ), FTFXannihilation( 0.0 ),
104  RadiusOfHNinteractions2( 0.0 ), FTFSlope( 0.0 ),
108  TarMinDiffMass( 0.0 ), TarMinNonDiffMass( 0.0 ),
109  AveragePt2( 0.0 ), ProbLogDistr( 0.0 ),
110  Pt2kink( 0.0 ),
115 {
116  for ( G4int i = 0; i < 4; i++ ) {
117  for ( G4int j = 0; j < 7; j++ ) {
118  ProcParams[i][j] = 0.0;
119  }
120  }
121 
122  G4int ProjectilePDGcode = particle->GetPDGEncoding();
123  G4int ProjectileabsPDGcode = std::abs( ProjectilePDGcode );
124  G4double ProjectileMass = particle->GetPDGMass();
125  G4double ProjectileMass2 = ProjectileMass * ProjectileMass;
126 
127  G4int ProjectileBaryonNumber( 0 ), AbsProjectileBaryonNumber( 0 ), AbsProjectileCharge( 0 );
128  G4bool ProjectileIsNucleus = false;
129 
130  if ( std::abs( particle->GetBaryonNumber() ) > 1 ) { // The projectile is a nucleus
131  ProjectileIsNucleus = true;
132  ProjectileBaryonNumber = particle->GetBaryonNumber();
133  AbsProjectileBaryonNumber = std::abs( ProjectileBaryonNumber );
134  AbsProjectileCharge = G4int( particle->GetPDGCharge() );
135  if ( ProjectileBaryonNumber > 1 ) {
136  ProjectilePDGcode = 2212; ProjectileabsPDGcode = 2212; // Proton
137  } else {
138  ProjectilePDGcode = -2212; ProjectileabsPDGcode = 2212; // Anti-Proton
139  }
140  ProjectileMass = G4Proton::Proton()->GetPDGMass();
141  ProjectileMass2 = sqr( ProjectileMass );
142  }
143 
144  G4double TargetMass = G4Proton::Proton()->GetPDGMass();
145  G4double TargetMass2 = TargetMass * TargetMass;
146 
147  G4double Plab = PlabPerParticle;
148  G4double Elab = std::sqrt( Plab*Plab + ProjectileMass2 );
149  G4double KineticEnergy = Elab - ProjectileMass;
150 
151  G4double S = ProjectileMass2 + TargetMass2 + 2.0*TargetMass*Elab;
152 
153  #ifdef debugFTFparams
154  G4cout << "--------- FTF Parameters --------------" << G4endl << "Proj Plab "
155  << ProjectilePDGcode << " " << Plab << G4endl << "Mass KinE " << ProjectileMass
156  << " " << KineticEnergy << G4endl << " A Z " << theA << " " << theZ << G4endl;
157  #endif
158 
159  G4double Ylab, Xtotal, Xelastic, Xannihilation;
160  G4int NumberOfTargetNucleons;
161 
162  Ylab = 0.5 * G4Log( (Elab + Plab)/(Elab - Plab) );
163 
164  G4double ECMSsqr = S/GeV/GeV;
165  G4double SqrtS = std::sqrt( S )/GeV;
166 
167  #ifdef debugFTFparams
168  G4cout << "Sqrt(s) " << SqrtS << G4endl;
169  #endif
170 
171  TargetMass /= GeV; TargetMass2 /= (GeV*GeV);
172  ProjectileMass /= GeV; ProjectileMass2 /= (GeV*GeV);
173 
174  // Andrea Dotti (13Jan2013):
175  // The following lines are changed for G4MT. Originally the code was:
176  // static G4ChipsComponentXS* _instance = new G4ChipsComponentXS(); // Witek Pokorski
177  // Note the code could go back at original if _instance could be shared among threads
178  if ( ! chipsComponentXSisInitialized ) {
179  chipsComponentXSisInitialized = true;
180  chipsComponentXSinstance = new G4ChipsComponentXS();
181  }
182  G4ChipsComponentXS* _instance = chipsComponentXSinstance;
183  FTFxsManager = _instance;
184 
185  Plab /= GeV;
186  G4double Xftf = 0.0;
187 
188  G4int NumberOfTargetProtons = theZ;
189  G4int NumberOfTargetNeutrons = theA - theZ;
190  NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
191 
192  if ( ProjectilePDGcode == 2212 || ProjectilePDGcode == 2112 ) { // Projectile is nucleon
194  G4double XtotPP = FTFxsManager->GetTotalElementCrossSection( Proton, KineticEnergy, 1, 0 ); //ALB
195 
197  G4double XtotPN = FTFxsManager->GetTotalElementCrossSection( Neutron, KineticEnergy, 1, 0 ); //ALB
198  G4double XelPP = FTFxsManager->GetElasticElementCrossSection( Proton, KineticEnergy, 1, 0 );
199  G4double XelPN = FTFxsManager->GetElasticElementCrossSection( Neutron, KineticEnergy, 1, 0 );
200 
201  #ifdef debugFTFparams
202  G4cout << "XsPP " << XtotPP/millibarn << " " << XelPP/millibarn << G4endl
203  << "XsPN " << XtotPN/millibarn << " " << XelPN/millibarn << G4endl;
204  #endif
205 
206  if ( ! ProjectileIsNucleus ) { // Projectile is hadron
207  Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN ) /
208  NumberOfTargetNucleons;
209  Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN ) /
210  NumberOfTargetNucleons;
211  } else { // Projectile is a nucleus
212  Xtotal = (
213  AbsProjectileCharge * NumberOfTargetProtons * XtotPP +
214  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
215  NumberOfTargetNeutrons * XtotPP
216  +
217  ( AbsProjectileCharge * NumberOfTargetNeutrons +
218  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
219  NumberOfTargetProtons ) * XtotPN
220  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
221  Xelastic= (
222  AbsProjectileCharge * NumberOfTargetProtons * XelPP +
223  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
224  NumberOfTargetNeutrons * XelPP
225  +
226  ( AbsProjectileCharge * NumberOfTargetNeutrons +
227  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
228  NumberOfTargetProtons ) * XelPN
229  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
230  }
231 
232  Xannihilation = 0.0;
233  Xtotal /= millibarn;
234  Xelastic /= millibarn;
235 
236  } else if ( ProjectilePDGcode < -1000 ) { // Projectile is anti_baryon
237 
238  G4double X_a( 0.0 ), X_b( 0.0 ), X_c( 0.0 ), X_d( 0.0 );
239  G4double MesonProdThreshold = ProjectileMass + TargetMass +
240  ( 2.0 * 0.14 + 0.016 ); // 2 Mpi + DeltaE;
241 
242  if ( PlabPerParticle < 40.0*MeV ) { // Low energy limits. Projectile at rest.
243  Xtotal = 1512.9; // mb
244  Xelastic = 473.2; // mb
245  X_a = 625.1; // mb
246  X_b = 9.780; // mb
247  X_c = 49.989; // mb
248  X_d = 6.614; // mb
249  } else { // Total and elastic cross section of PbarP interactions a'la Arkhipov
250  G4double LogS = G4Log( ECMSsqr / 33.0625 );
251  G4double Xasmpt = 36.04 + 0.304*LogS*LogS; // mb
252  LogS = G4Log( SqrtS / 20.74 );
253  G4double Basmpt = 11.92 + 0.3036*LogS*LogS; // GeV^(-2)
254  G4double R0 = std::sqrt( 0.40874044*Xasmpt - Basmpt ); // GeV^(-1)
255 
256  G4double FlowF = SqrtS / std::sqrt( ECMSsqr*ECMSsqr + ProjectileMass2*ProjectileMass2 +
257  TargetMass2*TargetMass2 - 2.0*ECMSsqr*ProjectileMass2
258  - 2.0*ECMSsqr*TargetMass2
259  - 2.0*ProjectileMass2*TargetMass2 );
260 
261  Xtotal = Xasmpt * ( 1.0 + 13.55*FlowF/R0/R0/R0*
262  (1.0 - 4.47/SqrtS + 12.38/ECMSsqr - 12.43/SqrtS/ECMSsqr) ); // mb
263 
264  Xasmpt = 4.4 + 0.101*LogS*LogS; // mb
265  Xelastic = Xasmpt * ( 1.0 + 59.27*FlowF/R0/R0/R0*
266  (1.0 - 6.95/SqrtS + 23.54/ECMSsqr - 25.34/SqrtS/ECMSsqr ) ); // mb
267 
268  //G4cout << "Param Xtotal Xelastic " << Xtotal << " " << Xelastic << G4endl
269  // << "FlowF " << FlowF << " SqrtS " << SqrtS << G4endl
270  // << "Param Xelastic-NaN " << Xelastic << " "
271  // << 1.5*16.654/pow(ECMSsqr/2.176/2.176,2.2) << " " << ECMSsqr << G4endl;
272 
273  X_a = 25.0*FlowF; // mb, 3-shirts diagram
274 
275  if ( SqrtS < MesonProdThreshold ) {
276  X_b = 3.13 + 140.0*G4Pow::GetInstance()->powA( MesonProdThreshold - SqrtS, 2.5 ); // mb anti-quark-quark annihilation
277  Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
278  } else {
279  X_b = 6.8/SqrtS; // mb anti-quark-quark annihilation
280  Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
281  }
282 
283  X_c = 2.0*FlowF*sqr( ProjectileMass + TargetMass )/ECMSsqr; // mb rearrangement
284 
285  X_d = 23.3/ECMSsqr; // mb anti-quark-quark string creation
286  }
287 
288  //G4cout << "Param Xtotal Xelastic " << Xtotal << " " << Xelastic << G4endl
289  // << "Para a b c d " << X_a << " " << X_b << " " << X_c << " " << X_d << G4endl;
290  // << "Para a b c d " << X_a << " " << 5.*X_b << " " << 5.*X_c << " " << 6.*X_d
291  // << G4endl;
292 
293  G4double Xann_on_P( 0.0), Xann_on_N( 0.0 );
294 
295  if ( ProjectilePDGcode == -2212 ) { // Pbar+P/N
296  Xann_on_P = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
297  Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
298  } else if ( ProjectilePDGcode == -2112 ) { // NeutrBar+P/N
299  Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
300  Xann_on_N = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
301  } else if ( ProjectilePDGcode == -3122 ) { // LambdaBar+P/N
302  Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
303  Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
304  } else if ( ProjectilePDGcode == -3112 ) { // Sigma-Bar+P/N
305  Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
306  Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
307  } else if ( ProjectilePDGcode == -3212 ) { // Sigma0Bar+P/N
308  Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
309  Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
310  } else if ( ProjectilePDGcode == -3222 ) { // Sigma+Bar+P/N
311  Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
312  Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
313  } else if ( ProjectilePDGcode == -3312 ) { // Xi-Bar+P/N
314  Xann_on_P = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
315  Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
316  } else if ( ProjectilePDGcode == -3322 ) { // Xi0Bar+P/N
317  Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
318  Xann_on_N = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
319  } else if ( ProjectilePDGcode == -3334 ) { // Omega-Bar+P/N
320  Xann_on_P = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
321  Xann_on_N = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
322  } else {
323  G4cout << "Unknown anti-baryon for FTF annihilation" << G4endl;
324  }
325 
326  //G4cout << "Sum " << Xann_on_P << G4endl;
327 
328  if ( ! ProjectileIsNucleus ) { // Projectile is anti-baryon
329  Xannihilation = ( NumberOfTargetProtons * Xann_on_P + NumberOfTargetNeutrons * Xann_on_N )
330  / NumberOfTargetNucleons;
331  } else { // Projectile is a nucleus
332  Xannihilation = (
333  ( AbsProjectileCharge * NumberOfTargetProtons +
334  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
335  NumberOfTargetNeutrons ) * Xann_on_P
336  +
337  ( AbsProjectileCharge * NumberOfTargetNeutrons +
338  ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
339  NumberOfTargetProtons ) * Xann_on_N
340  ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
341  }
342 
343  //G4double Xftf = 0.0;
344  MesonProdThreshold = ProjectileMass + TargetMass + (0.14 + 0.08); // Mpi + DeltaE
345  if ( SqrtS > MesonProdThreshold ) {
346  Xftf = 36.0 * ( 1.0 - MesonProdThreshold/SqrtS );
347  }
348 
349  Xtotal = Xelastic + Xannihilation + Xftf;
350 
351  #ifdef debugFTFparams
352  G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
353  << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << G4endl
354  << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal << " "
355  << Xannihilation/(Xtotal - Xelastic) << G4endl;
356  #endif
357 
358  } else if ( ProjectilePDGcode == 211 ) { // Projectile is PionPlus
359 
360  G4double XtotPiP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
362  G4double XtotPiN = FTFxsManager->GetTotalElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
363  G4double XelPiP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
364  G4double XelPiN = FTFxsManager->GetElasticElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
365  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
366  / NumberOfTargetNucleons;
367  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
368  / NumberOfTargetNucleons;
369  Xannihilation = 0.0;
370  Xtotal /= millibarn;
371  Xelastic /= millibarn;
372 
373  } else if ( ProjectilePDGcode == -211 ) { // Projectile is PionMinus
374 
375  G4double XtotPiP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
377  G4double XtotPiN = FTFxsManager->GetTotalElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
378  G4double XelPiP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
379  G4double XelPiN = FTFxsManager->GetElasticElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
380  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
381  / NumberOfTargetNucleons;
382  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
383  / NumberOfTargetNucleons;
384  Xannihilation = 0.0;
385  Xtotal /= millibarn;
386  Xelastic /= millibarn;
387 
388  } else if ( ProjectilePDGcode == 111 ) { // Projectile is PionZero
389 
391  G4double XtotPipP = FTFxsManager->GetTotalElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
393  G4double XtotPimP = FTFxsManager->GetTotalElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
394  G4double XelPipP = FTFxsManager->GetElasticElementCrossSection( PionPlus, KineticEnergy, 1, 0 );
395  G4double XelPimP = FTFxsManager->GetElasticElementCrossSection( PionMinus, KineticEnergy, 1, 0 );
396  G4double XtotPiP = ( XtotPipP + XtotPimP ) / 2.0;
397  G4double XtotPiN = XtotPiP;
398  G4double XelPiP = ( XelPipP + XelPimP ) / 2.0;
399  G4double XelPiN = XelPiP;
400  Xtotal = ( NumberOfTargetProtons * XtotPiP + NumberOfTargetNeutrons * XtotPiN )
401  / NumberOfTargetNucleons;
402  Xelastic = ( NumberOfTargetProtons * XelPiP + NumberOfTargetNeutrons * XelPiN )
403  / NumberOfTargetNucleons;
404  Xannihilation = 0.0;
405  Xtotal /= millibarn;
406  Xelastic /= millibarn;
407 
408  } else if ( ProjectilePDGcode == 321 ) { // Projectile is KaonPlus
409 
410  G4double XtotKP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
412  G4double XtotKN = FTFxsManager->GetTotalElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
413  G4double XelKP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
414  G4double XelKN = FTFxsManager->GetElasticElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
415  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
416  / NumberOfTargetNucleons;
417  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
418  / NumberOfTargetNucleons;
419  Xannihilation = 0.0;
420  Xtotal /= millibarn;
421  Xelastic /= millibarn;
422 
423  } else if ( ProjectilePDGcode == -321 ) { // Projectile is KaonMinus
424 
425  G4double XtotKP = FTFxsManager->GetTotalElementCrossSection( particle, KineticEnergy, 1, 0 );
427  G4double XtotKN = FTFxsManager->GetTotalElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
428  G4double XelKP = FTFxsManager->GetElasticElementCrossSection( particle, KineticEnergy, 1, 0 );
429  G4double XelKN = FTFxsManager->GetElasticElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
430  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
431  / NumberOfTargetNucleons;
432  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
433  / NumberOfTargetNucleons;
434  Xannihilation = 0.0;
435  Xtotal /= millibarn;
436  Xelastic /= millibarn;
437 
438  } else if ( ProjectilePDGcode == 311 || ProjectilePDGcode == 130 ||
439  ProjectilePDGcode == 310 ) { // Projectile is KaonZero
440 
442  G4double XtotKpP = FTFxsManager->GetTotalElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
444  G4double XtotKmP = FTFxsManager->GetTotalElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
445  G4double XelKpP = FTFxsManager->GetElasticElementCrossSection( KaonPlus, KineticEnergy, 1, 0 );
446  G4double XelKmP = FTFxsManager->GetElasticElementCrossSection( KaonMinus, KineticEnergy, 1, 0 );
447  G4double XtotKP = ( XtotKpP + XtotKmP ) / 2.0;
448  G4double XtotKN = XtotKP;
449  G4double XelKP = ( XelKpP + XelKmP ) / 2.0;
450  G4double XelKN = XelKP;
451  Xtotal = ( NumberOfTargetProtons * XtotKP + NumberOfTargetNeutrons * XtotKN )
452  / NumberOfTargetNucleons;
453  Xelastic = ( NumberOfTargetProtons * XelKP + NumberOfTargetNeutrons * XelKN )
454  / NumberOfTargetNucleons;
455  Xannihilation = 0.0;
456  Xtotal /= millibarn;
457  Xelastic /= millibarn;
458 
459  } else { // Projectile is undefined, Nucleon assumed
460 
462  G4double XtotPP = FTFxsManager->GetTotalElementCrossSection( Proton, KineticEnergy, 1, 0 );
464  G4double XtotPN = FTFxsManager->GetTotalElementCrossSection( Neutron, KineticEnergy, 1, 0 );
465  G4double XelPP = FTFxsManager->GetElasticElementCrossSection( Proton, KineticEnergy, 1, 0 );
466  G4double XelPN = FTFxsManager->GetElasticElementCrossSection( Neutron, KineticEnergy, 1, 0 );
467  Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN )
468  / NumberOfTargetNucleons;
469  Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN )
470  / NumberOfTargetNucleons;
471  Xannihilation = 0.0;
472  Xtotal /= millibarn;
473  Xelastic /= millibarn;
474 
475  };
476 
477  // Geometrical parameters
478  SetTotalCrossSection( Xtotal );
479  SetElastisCrossSection( Xelastic );
480  SetInelasticCrossSection( Xtotal - Xelastic );
481 
482  //G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
483  // << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << G4endl;
484  //if (Xtotal - Xelastic != 0.0 ) {
485  // G4cout << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal
486  // << " " << Xannihilation / (Xtotal - Xelastic) << G4endl;
487  //} else {
488  // G4cout << "Plab Xelastic/Xtotal, Xann " << Plab << " " << Xelastic/Xtotal
489  // << " " << Xannihilation << G4endl;
490  //}
491  //G4int Uzhi; G4cin >> Uzhi;
492 
493  // Interactions with elastic and inelastic collisions
494  SetProbabilityOfElasticScatt( Xtotal, Xelastic );
495  SetRadiusOfHNinteractions2( Xtotal/pi/10.0 );
496  if ( Xtotal - Xelastic == 0.0 ) {
498  } else {
499  SetProbabilityOfAnnihilation( Xannihilation / (Xtotal - Xelastic) );
500  }
501 
502  // No elastic scattering
503  //SetProbabilityOfElasticScatt( Xtotal, 0.0 );
504  //SetRadiusOfHNinteractions2( (Xtotal - Xelastic)/pi/10.0 );
505  //SetProbabilityOfAnnihilation( 1.0 );
506  //SetProbabilityOfAnnihilation( 0.0 );
507 
508  SetSlope( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 ); // Slope parameter of elastic scattering
509  // (GeV/c)^(-2))
510  //G4cout << "Slope " << GetSlope() << G4endl;
511  SetGamma0( GetSlope()*Xtotal/10.0/2.0/pi );
512 
513  // Parameters of elastic scattering
514  // Gaussian parametrization of elastic scattering amplitude assumed
515  SetAvaragePt2ofElasticScattering( 1.0/( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 )*GeV*GeV );
516  //G4cout << "AvaragePt2ofElasticScattering " << GetAvaragePt2ofElasticScattering() << G4endl;
517 
518  // Parameters of excitations
519 
520  G4double Xinel = Xtotal - Xelastic;
521  //G4cout << "Param ProjectilePDGcode " << ProjectilePDGcode << G4endl;
522 
523  if ( ProjectilePDGcode > 1000 ) { // Projectile is baryon
524  // Proc# A1 B1 A2 B2 A3 Atop Ymin
525  SetParams( 0, 13.71, 1.75, -214.5, 4.25, 0.0, 0.5 , 1.1 ); // Qexchange without Exc.
526  SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
527  if( Xinel > 0.) {
528  SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Projectile diffraction
529  SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Target diffraction
530  SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 );// Qexchange with Exc. Additional multiply
531  } else {
532  SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
533  SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
534  SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
535  }
536 
537  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
538  // It is not decided what to do with diffraction dissociation in Had-Nucl and Nucl-Nucl interactions
539  SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
540  //SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
541  }
542 
544  if ( NumberOfTargetNucleons > 26 ) {
546  } else {
548  }
549  SetProjMinDiffMass( 1.16 ); // GeV
550  SetProjMinNonDiffMass( 1.16 ); // GeV
551  SetTarMinDiffMass( 1.16 ); // GeV
552  SetTarMinNonDiffMass( 1.16 ); // GeV
553  SetAveragePt2( 0.15 ); // GeV^2
554  SetProbLogDistrPrD( 0.3 ); // Before it was: 0.5
555  SetProbLogDistr(0.3 ); // Before it was: 0.5
556 
557  } else if( ProjectilePDGcode < -1000 ) { // Projectile is anti_baryon
558 
559  // Proc# A1 B1 A2 B2 A3 Atop Ymin
560  SetParams( 0, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange without Exc.
561  SetParams( 1, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange with Exc.
562  if( Xinel > 0.) {
563  SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
564  SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
565  SetParams( 4, 1.0, 0.0 , 0.0, 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
566  } else {
567  SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
568  SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
569  SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
570  }
571 
572  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
573  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
574  //SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
575  }
576 
579  SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
580  SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
581  SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
582  SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
583  SetAveragePt2( 0.15 ); // GeV^2
584  SetProbLogDistrPrD( 0.3 );
585  SetProbLogDistr( 0.3 );
586 
587  } else if ( ProjectileabsPDGcode == 211 || ProjectilePDGcode == 111 ) { // Projectile is Pion
588 
589  // Proc# A1 B1 A2 B2 A3 Atop Ymin
590  SetParams( 0, 720.0, 2.5 , 2.3 , 1.0, 0., 1. , 2.7 );
591  SetParams( 1, 12.87, 0.5 , -44.91, 1.0, 0., 0. , 2.5 );
592  SetParams( 2, 0.086, 0. , -0.3 , 0.5, 0., 0. , 2.5 );
593  SetParams( 3, 32.8, 1.0 , -114.5 , 1.5, 0.084, 0. , 2.5 );
594  SetParams( 4, 1.0, 0.0 , -3.49, 0.5, 0.0, 0. , 2.5 ); // Qexchange with Exc. Additional multiply
595 
596  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
597  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
598  //SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
599  }
600 
601  SetDeltaProbAtQuarkExchange( 0.56 ); // (0.35)
602  SetProjMinDiffMass( 0.5 ); // (0.5) // GeV
603  SetProjMinNonDiffMass( 0.5 ); // (0.5) // GeV
604  SetTarMinDiffMass( 1.16 ); // GeV
605  SetTarMinNonDiffMass( 1.16 ); // GeV
606  SetAveragePt2( 0.15 ); // GeV^2
607  SetProbLogDistrPrD( 0.3 );
608  SetProbLogDistr( 0.3 );
609 
610  } else if ( ProjectileabsPDGcode == 321 || ProjectileabsPDGcode == 311 ||
611  ProjectilePDGcode == 130 || ProjectilePDGcode == 310 ) { // Projectile is Kaon
612 
613  // Proc# A1 B1 A2 B2 A3 Atop Ymin
614  SetParams( 0, 60.0 , 2.5 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
615  SetParams( 1, 6.0 , 1.0 , -24.33 , 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
616  SetParams( 2, 2.76, 1.2 , -22.5 , 2.7 ,0.04, 0.0 , 1.40 ); // Projectile diffraction
617  SetParams( 3, 1.09, 0.5 , -8.88 , 2. ,0.05, 0.0 , 1.40 ); // Target diffraction
618  SetParams( 4, 1.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
619 
620  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
621  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
622  //SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
623  }
624 
626  SetProjMinDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
627  SetProjMinNonDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
628  SetTarMinDiffMass( 1.16 ); // GeV
629  SetTarMinNonDiffMass( 1.16 ); // GeV
630  SetAveragePt2( 0.15 ); // GeV^2
631  SetProbLogDistrPrD( 0.5 );
632  SetProbLogDistr( 0.3 );
633 
634  } else { // Projectile is undefined, Nucleon assumed
635 
636  // Proc# A1 B1 A2 B2 A3 Atop Ymin
637  SetParams( 0, 13.71, 1.75, -214.5, 4.25, 0.0, 0.5 , 1.1 ); // Qexchange without Exc.
638  SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
639  if( Xinel > 0.) {
640  SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
641  SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
642  SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiply
643  } else {
644  SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
645  SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
646  SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
647  }
648  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
649  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
650  //SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
651  }
652  SetDeltaProbAtQuarkExchange( 0.0 ); // 7 June 2011
654  SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
655  SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
656  SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
657  SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
658  SetAveragePt2( 0.15 ); // GeV^2
659  SetProbLogDistrPrD( 0.3 );
660  SetProbLogDistr( 0.3 );
661 
662  }
663 
664  // Set parameters of a string kink
665  // SetPt2Kink( 6.0*GeV*GeV );
666  SetPt2Kink( 0.0*GeV*GeV ); // Uzhi Oct 2014 to switch off kinky strings
667  G4double Puubar( 1.0/3.0 ), Pddbar( 1.0/3.0 ), Pssbar( 1.0/3.0 ); // SU(3) symmetry
668  //G4double Puubar( 0.41 ), Pddbar( 0.41 ), Pssbar( 0.18 ); // Broken SU(3) symmetry
669  SetQuarkProbabilitiesAtGluonSplitUp( Puubar, Pddbar, Pssbar );
670 
671  // Set parameters of nuclear destruction
672  if ( ProjectileabsPDGcode < 1000 ) { // Meson projectile
673  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
674  //AR-18May2016 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
675  SetCofNuclearDestruction( 1.0* // AR-18May2016
676  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
679  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
680  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
683  } else if ( ProjectilePDGcode < -1000 ) { // for anti-baryon projectile
684  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
685  //AR-18May2016 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
686  SetCofNuclearDestruction( 1.0* // AR-18May2016
687  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
690  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
691  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
694  if ( Plab < 2.0 ) { // 2 GeV/c
695  // For slow anti-baryon we have to garanty putting on mass-shell
698  SetDofNuclearDestruction( 0.01 );
701  //SetExcitationEnergyPerWoundedNucleon( 0.0 ); // ?????
702  }
703  } else { // Projectile baryon assumed
704  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
705  //AR-18May2016 SetCofNuclearDestructionPr( 0.00481*G4double(AbsProjectileBaryonNumber)* // Uzhi 3.05.2015
706  SetCofNuclearDestructionPr( 1.0* // AR-18May2016
707  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
708  //AR-18May2016 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
709  SetCofNuclearDestruction( 1.0* // AR-18May2016
710  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
713  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
714  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
717  }
718 
719  //SetCofNuclearDestruction( 0.47*G4Exp( 2.0*(Ylab - 2.5) )/( 1.0 + G4Exp( 2.0*(Ylab - 2.5) ) ) );
720  //SetPt2ofNuclearDestruction( ( 0.035 + 0.1*G4Exp( 4.0*(Ylab - 3.0) )/( 1.0 + G4Exp( 4.0*(Ylab - 3.0) ) ) )*GeV*GeV );
721 
722  //SetMagQuarkExchange( 120.0 ); // 210.0 PipP
723  //SetSlopeQuarkExchange( 2.0 );
724  //SetDeltaProbAtQuarkExchange( 0.6 );
725  //SetProjMinDiffMass( 0.7 ); // GeV 1.1
726  //SetProjMinNonDiffMass( 0.7 ); // GeV
727  //SetProbabilityOfProjDiff( 0.0); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
728  //SetTarMinDiffMass( 1.1 ); // GeV
729  //SetTarMinNonDiffMass( 1.1 ); // GeV
730  //SetProbabilityOfTarDiff( 0.0 ); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
731 
732  //SetAveragePt2( 0.0 ); // GeV^2 0.3
733  //------------------------------------
734  //SetProbabilityOfElasticScatt( 1.0, 1.0); //(Xtotal, Xelastic);
735  //SetProbabilityOfProjDiff( 1.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 0->1
736  //SetProbabilityOfTarDiff( 4.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 2->4
737  //SetAveragePt2( 0.3 ); // (0.15)
738  //SetAvaragePt2ofElasticScattering( 0.0 );
739 
740  //SetMaxNumberOfCollisions( Plab, 6.0 ); //(4.0*(Plab + 0.01), Plab); // 6.0 );
741  //SetAveragePt2( 0.15 );
742  //SetCofNuclearDestruction(-1.);//( 0.75 ); // (0.25)
743  //SetExcitationEnergyPerWoundedNucleon(0.);//( 30.0*MeV ); // (75.0*MeV)
744  //SetDofNuclearDestruction(0.);//( 0.2 ); //0.4 // 0.3 0.5
745 
746  //SetPt2ofNuclearDestruction(0.);//(2.*0.075*GeV*GeV); //( 0.3*GeV*GeV ); // (0.168*GeV*GeV)
747  //SetMaxNumberOfCollisions( Plab, 78.0 ); // 3.0 )
748 
749  //G4cout << "Cnd " << GetCofNuclearDestruction() << G4endl;
750  //G4cout << "Dnd " << GetDofNuclearDestruction() << G4endl;
751  //G4cout << "Pt2 " << GetPt2ofNuclearDestruction()/GeV/GeV << G4endl;
752  //G4int Uzhi; G4cin >> Uzhi;
753 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
void SetProbLogDistr(const G4double aValue)
G4double ProjMinDiffMass
double S(double temp)
G4double ProbabilityOfAnnihilation
G4double Pt2ofNuclearDestruction
void SetCofNuclearDestructionPr(const G4double aValue)
G4double ProbLogDistrPrD
void SetProbLogDistrPrD(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
G4double CofNuclearDestruction
G4double ExcitationEnergyPerWoundedNucleon
G4double DeltaProbAtQuarkExchange
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
void SetGamma0(const G4double Gamma0)
G4GLOB_DLL std::ostream G4cout
void SetAveragePt2(const G4double aValue)
G4double DofNuclearDestruction
void SetElastisCrossSection(const G4double Xelastic)
bool G4bool
Definition: G4Types.hh:79
void SetPt2Kink(const G4double aValue)
G4double AvaragePt2ofElasticScattering
void SetRadiusOfHNinteractions2(const G4double Radius2)
G4double ProbabilityOfElasticScatt
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)
G4double RadiusOfHNinteractions2
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void SetDeltaProbAtQuarkExchange(const G4double aValue)
void SetTotalCrossSection(const G4double Xtotal)
void SetProjMinDiffMass(const G4double aValue)
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double MaxPt2ofNuclearDestruction
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)
G4double R2ofNuclearDestruction
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
G4double ProcParams[5][7]
void SetCofNuclearDestruction(const G4double aValue)
G4double ProbOfInelInteraction
static constexpr double GeV
Definition: G4SIunits.hh:217
void SetProbOfSameQuarkExchange(const G4double aValue)
G4double FTFXannihilation
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static constexpr double pi
Definition: G4SIunits.hh:75
G4double ProjMinNonDiffMass
virtual G4double GetTotalElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double N)
T sqr(const T &x)
Definition: templates.hh:145
G4double MaxNumberOfCollisions
double G4double
Definition: G4Types.hh:76
void SetProbabilityOfAnnihilation(const G4double aValue)
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
G4double CofNuclearDestructionPr
static constexpr double fermi
Definition: G4SIunits.hh:103
G4double GetPDGCharge() const
void SetPt2ofNuclearDestruction(const G4double aValue)
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
void SetTarMinDiffMass(const G4double aValue)
static constexpr double millibarn
Definition: G4SIunits.hh:106
void SetR2ofNuclearDestruction(const G4double aValue)
G4double TarMinNonDiffMass

Here is the call graph for this function:

G4FTFParameters::~G4FTFParameters ( )

Definition at line 87 of file G4FTFParameters.cc.

87 {}
G4FTFParameters::G4FTFParameters ( )

Definition at line 60 of file G4FTFParameters.cc.

60  :
61  FTFhNcmsEnergy( 0.0 ),
62  FTFxsManager( 0 ),
63  FTFXtotal( 0.0 ), FTFXelastic( 0.0 ), FTFXinelastic( 0.0 ), FTFXannihilation( 0.0 ),
65  RadiusOfHNinteractions2( 0.0 ), FTFSlope( 0.0 ),
69  TarMinDiffMass( 0.0 ), TarMinNonDiffMass( 0.0 ),
70  AveragePt2( 0.0 ), ProbLogDistr( 0.0 ),
71  Pt2kink( 0.0 ),
76 {
77  for ( G4int i = 0; i < 4; i++ ) {
78  for ( G4int j = 0; j < 7; j++ ) {
79  ProcParams[i][j] = 0.0;
80  }
81  }
82 }
G4double ProjMinDiffMass
G4double ProbabilityOfAnnihilation
G4double Pt2ofNuclearDestruction
G4double ProbLogDistrPrD
G4ChipsComponentXS * FTFxsManager
int G4int
Definition: G4Types.hh:78
G4double CofNuclearDestruction
G4double ExcitationEnergyPerWoundedNucleon
G4double DeltaProbAtQuarkExchange
G4double DofNuclearDestruction
G4double AvaragePt2ofElasticScattering
G4double ProbabilityOfElasticScatt
G4double RadiusOfHNinteractions2
G4double MaxPt2ofNuclearDestruction
G4double ProbOfSameQuarkExchange
G4double R2ofNuclearDestruction
G4double ProcParams[5][7]
G4double ProbOfInelInteraction
G4double FTFXannihilation
G4double ProjMinNonDiffMass
G4double MaxNumberOfCollisions
G4double CofNuclearDestructionPr
G4double TarMinNonDiffMass

Member Function Documentation

G4double G4FTFParameters::GammaElastic ( const G4double  impactsquare)
inline

Definition at line 212 of file G4FTFParameters.hh.

212  {
213  return ( FTFGamma0 * G4Exp( -FTFSlope * impactsquare ) );
214 }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4FTFParameters::GetAvaragePt2ofElasticScattering ( )
inline

Definition at line 415 of file G4FTFParameters.hh.

415  {
417 }
G4double AvaragePt2ofElasticScattering

Here is the caller graph for this function:

G4double G4FTFParameters::GetAveragePt2 ( )
inline

Definition at line 445 of file G4FTFParameters.hh.

445  {
446  return AveragePt2;
447 }

Here is the caller graph for this function:

G4double G4FTFParameters::GetCofNuclearDestruction ( )
inline

Definition at line 481 of file G4FTFParameters.hh.

481  {
482  return CofNuclearDestruction;
483 }
G4double CofNuclearDestruction
G4double G4FTFParameters::GetCofNuclearDestructionPr ( )
inline

Definition at line 477 of file G4FTFParameters.hh.

477  {
479 }
G4double CofNuclearDestructionPr
G4double G4FTFParameters::GetDeltaProbAtQuarkExchange ( )
inline

Definition at line 421 of file G4FTFParameters.hh.

421  {
423 }
G4double DeltaProbAtQuarkExchange

Here is the caller graph for this function:

G4double G4FTFParameters::GetDofNuclearDestruction ( )
inline

Definition at line 493 of file G4FTFParameters.hh.

493  {
494  return DofNuclearDestruction;
495 }
G4double DofNuclearDestruction
G4double G4FTFParameters::GetElasticCrossSection ( )
inline

Definition at line 381 of file G4FTFParameters.hh.

381  {
382  return FTFXelastic;
383 }
G4double G4FTFParameters::GetExcitationEnergyPerWoundedNucleon ( )
inline

Definition at line 489 of file G4FTFParameters.hh.

489  {
491 }
G4double ExcitationEnergyPerWoundedNucleon
G4double G4FTFParameters::GetInelasticCrossSection ( )
inline

Definition at line 385 of file G4FTFParameters.hh.

385  {
386  return FTFXinelastic;
387 }
G4double G4FTFParameters::GetInelasticProbability ( const G4double  impactsquare)
inline

Definition at line 405 of file G4FTFParameters.hh.

405  {
406  G4double Gamma = GammaElastic( impactsquare );
407  return 2*Gamma - Gamma*Gamma;
408 }
double G4double
Definition: G4Types.hh:76
G4double GammaElastic(const G4double impactsquare)

Here is the call graph for this function:

G4double G4FTFParameters::GetMaxNumberOfCollisions ( )
inline

Definition at line 469 of file G4FTFParameters.hh.

469  {
470  return MaxNumberOfCollisions;
471 }
G4double MaxNumberOfCollisions
G4double G4FTFParameters::GetMaxPt2ofNuclearDestruction ( )
inline

Definition at line 501 of file G4FTFParameters.hh.

501  {
503 }
G4double MaxPt2ofNuclearDestruction
G4double G4FTFParameters::GetProbabilityOfAnnihilation ( )
inline

Definition at line 410 of file G4FTFParameters.hh.

410  {
412 }
G4double ProbabilityOfAnnihilation
G4double G4FTFParameters::GetProbabilityOfElasticScatt ( )
inline

Definition at line 401 of file G4FTFParameters.hh.

401  {
403 }
G4double ProbabilityOfElasticScatt
G4double G4FTFParameters::GetProbabilityOfInteraction ( const G4double  impactsquare)
inline

Definition at line 393 of file G4FTFParameters.hh.

393  {
394  if ( RadiusOfHNinteractions2 > impactsquare ) {
395  return 1.0;
396  } else {
397  return 0.0;
398  }
399 }
G4double RadiusOfHNinteractions2

Here is the caller graph for this function:

G4double G4FTFParameters::GetProbLogDistr ( )
inline

Definition at line 453 of file G4FTFParameters.hh.

453  {
454  return ProbLogDistr;
455 }

Here is the caller graph for this function:

G4double G4FTFParameters::GetProbLogDistrPrD ( )
inline

Definition at line 449 of file G4FTFParameters.hh.

449  {
450  return ProbLogDistrPrD;
451 }
G4double ProbLogDistrPrD
G4double G4FTFParameters::GetProbOfInteraction ( )
inline

Definition at line 473 of file G4FTFParameters.hh.

473  {
474  return ProbOfInelInteraction;
475 }
G4double ProbOfInelInteraction
G4double G4FTFParameters::GetProbOfSameQuarkExchange ( )
inline

Definition at line 425 of file G4FTFParameters.hh.

425  {
427 }
G4double ProbOfSameQuarkExchange

Here is the caller graph for this function:

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

Definition at line 758 of file G4FTFParameters.cc.

758  {
759  G4double Prob( 0.0 );
760  if ( y < ProcParams[ProcN][6] ) {
761  Prob = ProcParams[ProcN][5];
762  if(Prob < 0.) Prob=0.;
763  return Prob;
764  }
765  Prob = ProcParams[ProcN][0] * G4Exp( -ProcParams[ProcN][1]*y ) +
766  ProcParams[ProcN][2] * G4Exp( -ProcParams[ProcN][3]*y ) +
767  ProcParams[ProcN][4];
768  if(Prob < 0.) Prob=0.;
769  return Prob;
770 }
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double ProcParams[5][7]
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4FTFParameters::GetProjMinDiffMass ( )
inline

Definition at line 429 of file G4FTFParameters.hh.

429  {
430  return ProjMinDiffMass;
431 }
G4double ProjMinDiffMass

Here is the caller graph for this function:

G4double G4FTFParameters::GetProjMinNonDiffMass ( )
inline

Definition at line 433 of file G4FTFParameters.hh.

433  {
434  return ProjMinNonDiffMass;
435 }
G4double ProjMinNonDiffMass

Here is the caller graph for this function:

G4double G4FTFParameters::GetPt2Kink ( )
inline

Definition at line 459 of file G4FTFParameters.hh.

459  {
460  return Pt2kink;
461 }

Here is the caller graph for this function:

G4double G4FTFParameters::GetPt2ofNuclearDestruction ( )
inline

Definition at line 497 of file G4FTFParameters.hh.

497  {
499 }
G4double Pt2ofNuclearDestruction
std::vector< G4double > G4FTFParameters::GetQuarkProbabilitiesAtGluonSplitUp ( )
inline

Definition at line 463 of file G4FTFParameters.hh.

463  {
465 }
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp

Here is the caller graph for this function:

G4double G4FTFParameters::GetR2ofNuclearDestruction ( )
inline

Definition at line 485 of file G4FTFParameters.hh.

485  {
486  return R2ofNuclearDestruction;
487 }
G4double R2ofNuclearDestruction
G4double G4FTFParameters::GetSlope ( )
inline

Definition at line 389 of file G4FTFParameters.hh.

389  {
390  return FTFSlope;
391 }

Here is the caller graph for this function:

G4double G4FTFParameters::GetTarMinDiffMass ( )
inline

Definition at line 437 of file G4FTFParameters.hh.

437  {
438  return TarMinDiffMass;
439 }

Here is the caller graph for this function:

G4double G4FTFParameters::GetTarMinNonDiffMass ( )
inline

Definition at line 441 of file G4FTFParameters.hh.

441  {
442  return TarMinNonDiffMass;
443 }
G4double TarMinNonDiffMass

Here is the caller graph for this function:

G4double G4FTFParameters::GetTotalCrossSection ( )
inline

Definition at line 377 of file G4FTFParameters.hh.

377  {
378  return FTFXtotal;
379 }
void G4FTFParameters::SetAvaragePt2ofElasticScattering ( const G4double  aPt2)
inline

Definition at line 264 of file G4FTFParameters.hh.

264  {
266 }
G4double AvaragePt2ofElasticScattering

Here is the caller graph for this function:

void G4FTFParameters::SetAveragePt2 ( const G4double  aValue)
inline

Definition at line 304 of file G4FTFParameters.hh.

304  {
306 }
static constexpr double GeV

Here is the caller graph for this function:

void G4FTFParameters::SetCofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 352 of file G4FTFParameters.hh.

352  {
353  CofNuclearDestruction = aValue;
354 }
G4double CofNuclearDestruction

Here is the caller graph for this function:

void G4FTFParameters::SetCofNuclearDestructionPr ( const G4double  aValue)
inline

Definition at line 348 of file G4FTFParameters.hh.

348  {
349  CofNuclearDestructionPr = aValue;
350 }
G4double CofNuclearDestructionPr

Here is the caller graph for this function:

void G4FTFParameters::SetDeltaProbAtQuarkExchange ( const G4double  aValue)
inline

Definition at line 280 of file G4FTFParameters.hh.

280  {
281  DeltaProbAtQuarkExchange = aValue;
282 }
G4double DeltaProbAtQuarkExchange

Here is the caller graph for this function:

void G4FTFParameters::SetDofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 364 of file G4FTFParameters.hh.

364  {
365  DofNuclearDestruction = aValue;
366 }
G4double DofNuclearDestruction

Here is the caller graph for this function:

void G4FTFParameters::SetElastisCrossSection ( const G4double  Xelastic)
inline

Definition at line 226 of file G4FTFParameters.hh.

226  {
227  FTFXelastic = Xelastic;
228 }

Here is the caller graph for this function:

void G4FTFParameters::SetExcitationEnergyPerWoundedNucleon ( const G4double  aValue)
inline

Definition at line 360 of file G4FTFParameters.hh.

360  {
362 }
G4double ExcitationEnergyPerWoundedNucleon

Here is the caller graph for this function:

void G4FTFParameters::SetGamma0 ( const G4double  Gamma0)
inline

Definition at line 259 of file G4FTFParameters.hh.

259  {
260  FTFGamma0 = Gamma0;
261 }

Here is the caller graph for this function:

void G4FTFParameters::SethNcmsEnergy ( const G4double  s)
inline

Definition at line 216 of file G4FTFParameters.hh.

216  {
217  FTFhNcmsEnergy = S;
218 }
double S(double temp)

Here is the call graph for this function:

void G4FTFParameters::SetInelasticCrossSection ( const G4double  Xinelastic)
inline

Definition at line 230 of file G4FTFParameters.hh.

230  {
231  FTFXinelastic = Xinelastic;
232 }

Here is the caller graph for this function:

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

Definition at line 331 of file G4FTFParameters.hh.

332  {
333  if ( Plab > Pbound ) {
334  MaxNumberOfCollisions = Plab/Pbound;
335  SetProbOfInteraction( -1.0 );
336  } else {
337  //MaxNumberOfCollisions = -1.0;
338  //SetProbOfInteraction( G4Exp( 0.25*(Plab-Pbound) ) );
340  SetProbOfInteraction( -1.0 );
341  }
342 }
void SetProbOfInteraction(const G4double aValue)
G4double MaxNumberOfCollisions

Here is the call graph for this function:

Here is the caller graph for this function:

void G4FTFParameters::SetMaxPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 372 of file G4FTFParameters.hh.

372  {
374 }
G4double MaxPt2ofNuclearDestruction

Here is the caller graph for this function:

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 270 of file G4FTFParameters.hh.

273  {
274  ProcParams[ProcN][0] = A1; ProcParams[ProcN][1] = B1;
275  ProcParams[ProcN][2] = A2; ProcParams[ProcN][3] = B2;
276  ProcParams[ProcN][4] = A3;
277  ProcParams[ProcN][5] = Atop; ProcParams[ProcN][6] = Ymin;
278 }
G4double ProcParams[5][7]

Here is the caller graph for this function:

void G4FTFParameters::SetProbabilityOfAnnihilation ( const G4double  aValue)
inline

Definition at line 247 of file G4FTFParameters.hh.

247  {
248  ProbabilityOfAnnihilation = aValue;
249 }
G4double ProbabilityOfAnnihilation

Here is the caller graph for this function:

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

Definition at line 234 of file G4FTFParameters.hh.

235  {
236  if ( Xtotal == 0.0 ) {
238  } else {
239  ProbabilityOfElasticScatt = Xelastic / Xtotal;
240  }
241 }
G4double ProbabilityOfElasticScatt

Here is the caller graph for this function:

void G4FTFParameters::SetProbabilityOfElasticScatt ( const G4double  aValue)
inline

Definition at line 243 of file G4FTFParameters.hh.

243  {
244  ProbabilityOfElasticScatt = aValue;
245 }
G4double ProbabilityOfElasticScatt
void G4FTFParameters::SetProbLogDistr ( const G4double  aValue)
inline

Definition at line 312 of file G4FTFParameters.hh.

312  {
313  ProbLogDistr = aValue;
314 }

Here is the caller graph for this function:

void G4FTFParameters::SetProbLogDistrPrD ( const G4double  aValue)
inline

Definition at line 308 of file G4FTFParameters.hh.

308  {
309  ProbLogDistrPrD = aValue;
310 }
G4double ProbLogDistrPrD

Here is the caller graph for this function:

void G4FTFParameters::SetProbOfInteraction ( const G4double  aValue)
inline

Definition at line 344 of file G4FTFParameters.hh.

344  {
345  ProbOfInelInteraction = aValue;
346 }
G4double ProbOfInelInteraction

Here is the caller graph for this function:

void G4FTFParameters::SetProbOfSameQuarkExchange ( const G4double  aValue)
inline

Definition at line 284 of file G4FTFParameters.hh.

284  {
285  ProbOfSameQuarkExchange = aValue;
286 }
G4double ProbOfSameQuarkExchange

Here is the caller graph for this function:

void G4FTFParameters::SetProjMinDiffMass ( const G4double  aValue)
inline

Definition at line 288 of file G4FTFParameters.hh.

288  {
289  ProjMinDiffMass = aValue*CLHEP::GeV;
290 }
G4double ProjMinDiffMass
static constexpr double GeV

Here is the caller graph for this function:

void G4FTFParameters::SetProjMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 292 of file G4FTFParameters.hh.

292  {
294 }
static constexpr double GeV
G4double ProjMinNonDiffMass

Here is the caller graph for this function:

void G4FTFParameters::SetPt2Kink ( const G4double  aValue)
inline

Definition at line 318 of file G4FTFParameters.hh.

318  {
319  Pt2kink = aValue;
320 }

Here is the caller graph for this function:

void G4FTFParameters::SetPt2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 368 of file G4FTFParameters.hh.

368  {
369  Pt2ofNuclearDestruction = aValue;
370 }
G4double Pt2ofNuclearDestruction

Here is the caller graph for this function:

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

Definition at line 322 of file G4FTFParameters.hh.

324  {
325  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar );
326  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar );
327  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar + Pssbar );
328 }
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp

Here is the caller graph for this function:

void G4FTFParameters::SetR2ofNuclearDestruction ( const G4double  aValue)
inline

Definition at line 356 of file G4FTFParameters.hh.

356  {
357  R2ofNuclearDestruction = aValue;
358 }
G4double R2ofNuclearDestruction

Here is the caller graph for this function:

void G4FTFParameters::SetRadiusOfHNinteractions2 ( const G4double  Radius2)
inline

Definition at line 251 of file G4FTFParameters.hh.

251  {
252  RadiusOfHNinteractions2 = Radius2;
253 }
G4double RadiusOfHNinteractions2

Here is the caller graph for this function:

void G4FTFParameters::SetSlope ( const G4double  Slope)
inline

Definition at line 255 of file G4FTFParameters.hh.

255  {
256  FTFSlope = 12.84 / Slope; // Slope is in GeV^-2, FTFSlope in fm^-2
257 }

Here is the caller graph for this function:

void G4FTFParameters::SetTarMinDiffMass ( const G4double  aValue)
inline

Definition at line 296 of file G4FTFParameters.hh.

296  {
297  TarMinDiffMass = aValue*CLHEP::GeV;
298 }
static constexpr double GeV

Here is the caller graph for this function:

void G4FTFParameters::SetTarMinNonDiffMass ( const G4double  aValue)
inline

Definition at line 300 of file G4FTFParameters.hh.

300  {
301  TarMinNonDiffMass = aValue*CLHEP::GeV;
302 }
static constexpr double GeV
G4double TarMinNonDiffMass

Here is the caller graph for this function:

void G4FTFParameters::SetTotalCrossSection ( const G4double  Xtotal)
inline

Definition at line 222 of file G4FTFParameters.hh.

222  {
223  FTFXtotal = Xtotal;
224 }

Here is the caller graph for this function:

Member Data Documentation

G4double G4FTFParameters::AvaragePt2ofElasticScattering

Definition at line 169 of file G4FTFParameters.hh.

G4double G4FTFParameters::AveragePt2

Definition at line 184 of file G4FTFParameters.hh.

G4double G4FTFParameters::CofNuclearDestruction

Definition at line 196 of file G4FTFParameters.hh.

G4double G4FTFParameters::CofNuclearDestructionPr

Definition at line 195 of file G4FTFParameters.hh.

G4double G4FTFParameters::DeltaProbAtQuarkExchange

Definition at line 175 of file G4FTFParameters.hh.

G4double G4FTFParameters::DofNuclearDestruction

Definition at line 201 of file G4FTFParameters.hh.

G4double G4FTFParameters::ExcitationEnergyPerWoundedNucleon

Definition at line 199 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFGamma0

Definition at line 170 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFhNcmsEnergy

Definition at line 155 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFSlope

Definition at line 168 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFXannihilation

Definition at line 164 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFXelastic

Definition at line 162 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFXinelastic

Definition at line 163 of file G4FTFParameters.hh.

G4ChipsComponentXS* G4FTFParameters::FTFxsManager

Definition at line 158 of file G4FTFParameters.hh.

G4double G4FTFParameters::FTFXtotal

Definition at line 161 of file G4FTFParameters.hh.

G4double G4FTFParameters::MaxNumberOfCollisions

Definition at line 192 of file G4FTFParameters.hh.

G4double G4FTFParameters::MaxPt2ofNuclearDestruction

Definition at line 203 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbabilityOfAnnihilation

Definition at line 165 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbabilityOfElasticScatt

Definition at line 166 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbLogDistr

Definition at line 185 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbLogDistrPrD

Definition at line 180 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbOfInelInteraction

Definition at line 193 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProbOfSameQuarkExchange

Definition at line 176 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProcParams[5][7]

Definition at line 173 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProjMinDiffMass

Definition at line 178 of file G4FTFParameters.hh.

G4double G4FTFParameters::ProjMinNonDiffMass

Definition at line 179 of file G4FTFParameters.hh.

G4double G4FTFParameters::Pt2kink

Definition at line 188 of file G4FTFParameters.hh.

G4double G4FTFParameters::Pt2ofNuclearDestruction

Definition at line 202 of file G4FTFParameters.hh.

std::vector< G4double > G4FTFParameters::QuarkProbabilitiesAtGluonSplitUp

Definition at line 189 of file G4FTFParameters.hh.

G4double G4FTFParameters::R2ofNuclearDestruction

Definition at line 197 of file G4FTFParameters.hh.

G4double G4FTFParameters::RadiusOfHNinteractions2

Definition at line 167 of file G4FTFParameters.hh.

G4double G4FTFParameters::TarMinDiffMass

Definition at line 181 of file G4FTFParameters.hh.

G4double G4FTFParameters::TarMinNonDiffMass

Definition at line 182 of file G4FTFParameters.hh.


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