Geant4  10.02.p01
G4FTFParameters.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id: G4FTFParameters.cc 91775 2015-08-05 14:42:39Z gcosmo $
28 // GEANT4 tag $Name: $
29 //
30 
31 #include <utility>
32 
33 #include "G4FTFParameters.hh"
34 
35 #include "G4ios.hh"
36 #include "G4PhysicalConstants.hh"
37 #include "G4SystemOfUnits.hh"
38 
39 #include "G4ParticleDefinition.hh"
40 
41 #include "G4Proton.hh"
42 #include "G4Neutron.hh"
43 
44 #include "G4PionPlus.hh"
45 #include "G4PionMinus.hh"
46 #include "G4KaonPlus.hh"
47 #include "G4KaonMinus.hh"
48 
49 #include "G4Exp.hh"
50 #include "G4Log.hh"
51 #include "G4Pow.hh"
52 
53 //============================================================================
54 
55 //#define debugFTFparams
56 
57 
58 //============================================================================
59 
61  FTFhNcmsEnergy( 0.0 ),
62  FTFxsManager( 0 ),
63  FTFXtotal( 0.0 ), FTFXelastic( 0.0 ), FTFXinelastic( 0.0 ), FTFXannihilation( 0.0 ),
64  ProbabilityOfAnnihilation( 0.0 ), ProbabilityOfElasticScatt( 0.0 ),
65  RadiusOfHNinteractions2( 0.0 ), FTFSlope( 0.0 ),
66  AvaragePt2ofElasticScattering( 0.0 ), FTFGamma0( 0.0 ),
67  DeltaProbAtQuarkExchange( 0.0 ), ProbOfSameQuarkExchange( 0.0 ),
68  ProjMinDiffMass( 0.0 ), ProjMinNonDiffMass( 0.0 ), ProbLogDistrPrD(0.0),
69  TarMinDiffMass( 0.0 ), TarMinNonDiffMass( 0.0 ),
70  AveragePt2( 0.0 ), ProbLogDistr( 0.0 ),
71  Pt2kink( 0.0 ),
72  MaxNumberOfCollisions( 0.0 ), ProbOfInelInteraction( 0.0 ),
73  CofNuclearDestructionPr( 0.0 ), CofNuclearDestruction( 0.0 ),
74  R2ofNuclearDestruction( 0.0 ), ExcitationEnergyPerWoundedNucleon( 0.0 ),
75  DofNuclearDestruction( 0.0 ), Pt2ofNuclearDestruction( 0.0 ), MaxPt2ofNuclearDestruction( 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 }
83 
84 
85 //============================================================================
86 
88 
89 
90 //============================================================================
91 
94 
95 
96 //============================================================================
97 
99  G4int theA, G4int theZ, G4double PlabPerParticle ) :
100  FTFhNcmsEnergy( 0.0 ),
101  FTFxsManager( 0 ),
102  FTFXtotal( 0.0 ), FTFXelastic( 0.0 ), FTFXinelastic( 0.0 ), FTFXannihilation( 0.0 ),
103  ProbabilityOfAnnihilation( 0.0 ), ProbabilityOfElasticScatt( 0.0 ),
104  RadiusOfHNinteractions2( 0.0 ), FTFSlope( 0.0 ),
105  AvaragePt2ofElasticScattering( 0.0 ), FTFGamma0( 0.0 ),
106  DeltaProbAtQuarkExchange( 0.0 ), ProbOfSameQuarkExchange( 0.0 ),
107  ProjMinDiffMass( 0.0 ), ProjMinNonDiffMass( 0.0 ), ProbLogDistrPrD(0.0),
108  TarMinDiffMass( 0.0 ), TarMinNonDiffMass( 0.0 ),
109  AveragePt2( 0.0 ), ProbLogDistr( 0.0 ),
110  Pt2kink( 0.0 ),
111  MaxNumberOfCollisions( 0.0 ), ProbOfInelInteraction( 0.0 ),
112  CofNuclearDestructionPr( 0.0 ), CofNuclearDestruction( 0.0 ),
113  R2ofNuclearDestruction( 0.0 ), ExcitationEnergyPerWoundedNucleon( 0.0 ),
114  DofNuclearDestruction( 0.0 ), Pt2ofNuclearDestruction( 0.0 ), MaxPt2ofNuclearDestruction( 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
181  }
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; // Uzhi 25.04.2012
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 // Uzhi Oct 2014
554  SetProbLogDistrPrD( 0.3 ); // Uzhi Oct 2014 0.5
555  SetProbLogDistr(0.3 ); // 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 // Uzhi Oct 2014
584  SetProbLogDistrPrD( 0.3 ); // Uzhi Oct 2014
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 // Uzhi Oct 2014
607  SetProbLogDistrPrD( 0.3 ); // Uzhi Oct 2014
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  }
625  SetProjMinDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
626  SetProjMinNonDiffMass( 0.7 ); // (1.4) // (0.7) // GeV
627  SetTarMinDiffMass( 1.16 ); // GeV
628  SetTarMinNonDiffMass( 1.16 ); // GeV
629  SetAveragePt2( 0.15 ); // GeV^2 // Uzhi Oct 2014
630  SetProbLogDistrPrD( 0.5 ); // Uzhi Oct 2014
631  SetProbLogDistr( 0.3 ); // Uzhi 5.06.2012
632 
633  } else { // Projectile is undefined, Nucleon assumed
634 
635  // Proc# A1 B1 A2 B2 A3 Atop Ymin
636  SetParams( 0, 13.71, 1.75, -214.5, 4.25, 0.0, 0.5 , 1.1 ); // Qexchange without Exc.
637  SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
638  if( Xinel > 0.) {
639  SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
640  SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
641  SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiply
642  } else {
643  SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
644  SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
645  SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
646  }
647  if ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 ) {
648  SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
649 // SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
650  }
651  SetDeltaProbAtQuarkExchange( 0.0 ); // 7 June 2011
653  SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
654  SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
655  SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
656  SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
657  SetAveragePt2( 0.15 ); // GeV^2 // Uzhi Oct 2014
658  SetProbLogDistrPrD( 0.3 ); // Uzhi Oct 2014
659  SetProbLogDistr( 0.3 );
660 
661  }
662 
663  // Set parameters of a string kink
664 // SetPt2Kink( 6.0*GeV*GeV ); // Uzhi Oct 2014
665  SetPt2Kink( 0.0*GeV*GeV ); // Uzhi Oct 2014 // Uzhi to switch off kinky strings
666  G4double Puubar( 1.0/3.0 ), Pddbar( 1.0/3.0 ), Pssbar( 1.0/3.0 ); // SU(3) symmetry
667 // G4double Puubar( 0.41 ), Pddbar( 0.41 ), Pssbar( 0.18 ); // Broken SU(3) symmetry
668  SetQuarkProbabilitiesAtGluonSplitUp( Puubar, Pddbar, Pssbar );
669 
670  // Set parameters of nuclear destruction
671  if ( ProjectileabsPDGcode < 1000 ) { // Meson projectile
672  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
673  SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
674  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
677  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
678  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
679  SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
680  SetExcitationEnergyPerWoundedNucleon( 40.0*MeV ); // Uzhi March 2015: 100 -> 40
681  } else if ( ProjectilePDGcode < -1000 ) { // for anti-baryon projectile
682  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
683  SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
684  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
687  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
688  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
689  SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
690  SetExcitationEnergyPerWoundedNucleon( 40.0*MeV ); // Uzhi March 2015: 100 -> 20
691  if ( Plab < 2.0 ) { // 2 GeV/c
692  // For slow anti-baryon we have to garanty putting on mass-shell
694  SetR2ofNuclearDestruction( 1.5*fermi*fermi );
695  SetDofNuclearDestruction( 0.01 );
696  SetPt2ofNuclearDestruction( 0.035*GeV*GeV );
697  SetMaxPt2ofNuclearDestruction( 0.04*GeV*GeV );
698  //SetExcitationEnergyPerWoundedNucleon( 0.0 ); // ?????
699  }
700  } else { // Projectile baryon assumed
701  SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
702  SetCofNuclearDestructionPr( 0.00481*G4double(AbsProjectileBaryonNumber)* // Uzhi 3.05.2015
703  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
704  SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)* // Uzhi 3.05.2015
705  G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
708  SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
709  ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
710  SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
711  SetExcitationEnergyPerWoundedNucleon( 40.0*MeV ); // Uzhi March 2015: 100 -> 40
712  }
713 
714  //SetCofNuclearDestruction( 0.47*G4Exp( 2.0*(Ylab - 2.5) )/( 1.0 + G4Exp( 2.0*(Ylab - 2.5) ) ) );
715  //SetPt2ofNuclearDestruction( ( 0.035 + 0.1*G4Exp( 4.0*(Ylab - 3.0) )/( 1.0 + G4Exp( 4.0*(Ylab - 3.0) ) ) )*GeV*GeV );
716 
717  //SetMagQuarkExchange( 120.0 ); // 210.0 PipP
718  //SetSlopeQuarkExchange( 2.0 );
719  //SetDeltaProbAtQuarkExchange( 0.6 );
720  //SetProjMinDiffMass( 0.7 ); // GeV 1.1
721  //SetProjMinNonDiffMass( 0.7 ); // GeV
722  //SetProbabilityOfProjDiff( 0.0); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
723  //SetTarMinDiffMass( 1.1 ); // GeV
724  //SetTarMinNonDiffMass( 1.1 ); // GeV
725  //SetProbabilityOfTarDiff( 0.0 ); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
726 
727 //SetAveragePt2( 0.0 ); // GeV^2 0.3
728  //------------------------------------
729 //SetProbabilityOfElasticScatt( 1.0, 1.0); //(Xtotal, Xelastic);
730 //SetProbabilityOfProjDiff( 1.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 0->1
731 //SetProbabilityOfTarDiff( 4.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 2->4
732 //SetAveragePt2( 0.3 ); // (0.15)
733 //SetAvaragePt2ofElasticScattering( 0.0 );
734 
735 //SetMaxNumberOfCollisions( Plab, 6.0 ); //(4.0*(Plab + 0.01), Plab); // 6.0 );
736 //SetAveragePt2( 0.15 );
737 //SetCofNuclearDestruction(-1.);//( 0.75 ); // (0.25)
738 //SetExcitationEnergyPerWoundedNucleon(0.);//( 30.0*MeV ); // (75.0*MeV)
739 //SetDofNuclearDestruction(0.);//( 0.2 ); //0.4 // 0.3 0.5
740 
741 //SetPt2ofNuclearDestruction(0.);//(2.*0.075*GeV*GeV); //( 0.3*GeV*GeV ); // (0.168*GeV*GeV)
742 //SetMaxNumberOfCollisions( Plab, 78.0 ); // 3.0 )
743 
744 //G4cout << "Cnd " << GetCofNuclearDestruction() << G4endl;
745 //G4cout << "Dnd " << GetDofNuclearDestruction() << G4endl;
746 //G4cout << "Pt2 " << GetPt2ofNuclearDestruction()/GeV/GeV << G4endl;
747 //G4int Uzhi; G4cin >> Uzhi;
748 
749 }
750 
751 
752 //============================================================================
753 
755  G4double Prob( 0.0 );
756  if ( y < ProcParams[ProcN][6] ) {
757  Prob = ProcParams[ProcN][5];
758  if(Prob < 0.) Prob=0.; // Uzhi Oct 2014
759  return Prob;
760  }
761  Prob = ProcParams[ProcN][0] * G4Exp( -ProcParams[ProcN][1]*y ) +
762  ProcParams[ProcN][2] * G4Exp( -ProcParams[ProcN][3]*y ) +
763  ProcParams[ProcN][4];
764  if(Prob < 0.) Prob=0.; // Uzhi Oct 2014
765  return Prob;
766 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
static const double MeV
Definition: G4SIunits.hh:211
static G4ThreadLocal G4ChipsComponentXS * chipsComponentXSinstance
void SetProbLogDistr(const G4double aValue)
double S(double temp)
void SetCofNuclearDestructionPr(const G4double aValue)
void SetProbLogDistrPrD(const G4double aValue)
#define G4ThreadLocal
Definition: tls.hh:89
void SetTarMinNonDiffMass(const G4double aValue)
void SetDofNuclearDestruction(const G4double aValue)
void SetSlope(const G4double Slope)
G4ChipsComponentXS * FTFxsManager
int G4int
Definition: G4Types.hh:78
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)
static const double GeV
Definition: G4SIunits.hh:214
virtual G4double GetElasticElementCrossSection(const G4ParticleDefinition *aParticle, G4double kinEnergy, G4int Z, G4double N)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
G4double GetProcProb(const G4int ProcN, const G4double y)
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 G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetPDGMass() const
static const double pi
Definition: G4SIunits.hh:74
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
G4double ProcParams[5][7]
void SetCofNuclearDestruction(const G4double aValue)
static const double millibarn
Definition: G4SIunits.hh:105
void SetProbOfSameQuarkExchange(const G4double aValue)
#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
static G4ThreadLocal bool chipsComponentXSisInitialized
void SetPt2ofNuclearDestruction(const G4double aValue)
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
void SetTarMinDiffMass(const G4double aValue)
static const double fermi
Definition: G4SIunits.hh:102
void SetR2ofNuclearDestruction(const G4double aValue)