Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FTFParameters.hh
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//
28#ifndef G4FTFParameters_h
29#define G4FTFParameters_h 1
30
32#include <vector>
33#include "G4Types.hh"
34#include "G4Exp.hh"
35
39
40
41 // NOTE: the settings are different for:
42 // * baryons projectile
43 // * anti-baryons projectile
44 // * pions (chg or pi0) projectile
45 // * kaons projectile (pdg = +/-321, 311, 130, or 310)
46 // * "undefined" projectile - nucleon assumed
47
49 public:
50
51 //dtor
53
54 // parameters of excitation
55 //
56 // Proc=0 --> Qexchg w/o excitation
57 //
58 double GetProc0A1() const { return fProc0A1; }
59 double GetProc0B1() const { return fProc0B1; }
60 double GetProc0A2() const { return fProc0A2; }
61 double GetProc0B2() const { return fProc0B2; }
62 double GetProc0A3() const { return fProc0A3; }
63 double GetProc0Atop() const { return fProc0Atop; }
64 double GetProc0Ymin() const { return fProc0Ymin; }
65 //
66 // Proc=1 --> Qexchg w/excitation
67 //
68 double GetProc1A1() const { return fProc1A1; }
69 double GetProc1B1() const { return fProc1B1; }
70 double GetProc1A2() const { return fProc1A2; }
71 double GetProc1B2() const { return fProc1B2; }
72 double GetProc1A3() const { return fProc1A3; }
73 double GetProc1Atop() const { return fProc1Atop; }
74 double GetProc1Ymin() const { return fProc1Ymin; }
75 //
76 // Proc=2 & Proc=3 in case ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )
77 // Update: Proc=2 & Proc=3 in case ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 )
78 // (diffraction dissociation)
79 //
80 // Other parameters have a complex form for baryon projectile
81 // although they're just numbers for e.g. pions projectile
82 //
83 // Proc=2 --> Projectile diffraction
84 //
85 double GetProc2A1() const { return fProc2A1; }
86 double GetProc2B1() const { return fProc2B1; }
87 double GetProc2A2() const { return fProc2A2; }
88 double GetProc2B2() const { return fProc2B2; }
89 double GetProc2A3() const { return fProc2A3; }
90 double GetProc2Atop() const { return fProc2Atop; }
91 double GetProc2Ymin() const { return fProc2Ymin; }
92 //
93 // Proc=3 --> Target diffraction
94 //
95 double GetProc3A1() const { return fProc3A1; }
96 double GetProc3B1() const { return fProc3B1; }
97 double GetProc3A2() const { return fProc3A2; }
98 double GetProc3B2() const { return fProc3B2; }
99 double GetProc3A3() const { return fProc3A3; }
100 double GetProc3Atop() const { return fProc3Atop; }
101 double GetProc3Ymin() const { return fProc3Ymin; }
102 //
105 //
106 // Proc=4 --> Qexchg "w/additional multiplier" in excitation
107 //
108 double GetProc4A1() const { return fProc4A1; }
109 double GetProc4B1() const { return fProc4B1; }
110 double GetProc4A2() const { return fProc4A2; }
111 double GetProc4B2() const { return fProc4B2; }
112 double GetProc4A3() const { return fProc4A3; }
113 double GetProc4Atop() const { return fProc4Atop; }
114 double GetProc4Ymin() const { return fProc4Ymin; }
115 //
116 //
119 double GetProjMinDiffMass() const { return fProjMinDiffMass; }
121 double GetTgtMinDiffMass() const { return fTgtMinDiffMass; }
122 double GetTgtMinNonDiffMass() const { return fTgtMinNonDiffMass; }
123 double GetAveragePt2() const { return fAveragePt2; }
124 double GetProbLogDistrPrD() const { return fProbLogDistrPrD; }
125 double GetProbLogDistr() const { return fProbLogDistr; }
126
127 // NOTE (JVY): There is also the Pt2Kind parameter but for now it's set to 0., so we'll leave it aside
128
129 // --> FIXME !!! --> void Get/SetBaryonMaxNumberOfCollisions( const double, const double ); // 1st is Plab, 2nd - D=2.
130 //
143 //
144 // separately for baryons, mesons, etc.
145 //
148 double GetDofNuclearDestruct() const { return fDofNuclearDestruct; }
150
151 protected:
152
153 // ctor
155
156 // parameters of excitation
157 //
158 //
159 // these are for Inelastic interactions, i.e. Xinelastic=(Xtotal-Xelastix)>0.
160 // for elastic, all the A's & B's, Atop & Ymin are zeros
161 // general formula: Pp = A1*exp(B1*Y) + A2*exp(B2*Y) + A3
162 // but if Y<Ymin, then Pp=max(0.,Atop)
163 // for details, see also G4FTFParameters::GetProcProb( ProcN, y )
164 //
165 // Proc=0 --> Qexchg w/o excitation
166 double fProc0A1; // D=13.71
167 double fProc0B1; // D=1.75
168 double fProc0A2; // D=-30.69 (or -214.5 as in Doc ?)
169 double fProc0B2; // D=3. ( or 4. as in Doc ?)
170 double fProc0A3; // D=0.
171 double fProc0Atop; // D=1. ( or 0.5 as in Doc ?)
172 double fProc0Ymin; // D=0.93 (or 1.1 as in Doc ?)
173 // Proc=1 --> Qexchg w/excitation
174 double fProc1A1; // D=25.
175 double fProc1B1; // D=1.
176 double fProc1A2; // D=-50.34
177 double fProc1B2; // D=1.5
178 double fProc1A3; // D=0.
179 double fProc1Atop; // D=0.
180 double fProc1Ymin; // D=1.4
181 //
182 // NOTE: Proc #2 & 3 are projectile & target diffraction
183 // they have more complex definition of A1 & A2
184 // for *baryons* although they're just numbers for pions
185 // (example for baryons below)
186 // SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Projectile diffraction
187 // SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Target diffraction
188 //
189 // Also, for ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )
190 // projectile and/or target diffraction (dissociation) may be switched ON/OFF
193 // Proc=2 --> Projectile diffraction
194 double fProc2A1;
195 double fProc2B1;
196 double fProc2A2;
197 double fProc2B2;
198 double fProc2A3;
199 double fProc2Atop;
200 double fProc2Ymin;
201 // Proc=3 --> Target diffraction
202 double fProc3A1;
203 double fProc3B1;
204 double fProc3A2;
205 double fProc3B2;
206 double fProc3A3;
207 double fProc3Atop;
208 double fProc3Ymin;
209 // Proc=4 --> Qexchg w/additional multiplier in excitation
210 double fProc4A1; // D=0.6 (or 1. as in Doc ?)
211 double fProc4B1; // D=0.
212 double fProc4A2; // D=-1.2 (or -2.01 as in Doc ?)
213 double fProc4B2; // D=0.5
214 double fProc4A3; // D=0.
215 double fProc4Atop; // D=0.
216 double fProc4Ymin; // D=1.4
217 //
218 // parameters of participating baryon excitation
219 // NOTE: baryon ot HADRON ???
220 // NOTE: this parameters (as C++ class data members) are used for all types of hadrons
221 // but the values for a specific group of particles can be are different from
222 // another group of particles
223 // the defaults listed under coments are for baryons,
224 // and they may be different or the same for other hadrons (e.g. mesons)
225 //
227 double fProbOfSameQuarkExchange; // D=0. if A<=26, otherwise D=1.
228 double fProjMinDiffMass; // projectile, D=1.16GeV
229 double fProjMinNonDiffMass; // projectile, D=1.16GeV
230 double fTgtMinDiffMass; // target, D=1.16GeV
231 double fTgtMinNonDiffMass; // target, D=1.16GeV
232 double fAveragePt2; // D=0.3GeV**2 ( or 0.15 as in the Doc ???)
233 double fProbLogDistrPrD; // D=0.55 (or 0.6 ??? or 0.3 ???)
234 double fProbLogDistr; // D=0.55 (or 0.6 ??? or 0.3 ???)
235
236 // parameters of nuclear distruction
237 //
238 // NOTE (JVY): there're 3 cases here:
239 // * baryon projectile
240 // * anti-baryon projectile
241 // * meson projectile
242 //
243 // double fBaryonMaxNumberOfCollisions; // D=2.
244 // void SetBaryonProbOfInteraction( const double ); // ??? this is prob. of inelastic interaction
245 // that is set internally based on certain conditions...
246 // general (i.e. for used for baryons,anti-baryons, and mesons)
247 // NOTE: these parameters have stayed THE SAME for quite a while
248 double fNuclearProjDestructP1; // D=0.00481 in 10.3.ref04 !!!
249 // BUT !!! In 10.3.ref04 as well as in 10.2-seriesit's multiplied of AbsProjectileBaryonNumber
250 // which somehow is 0 for the proton projectile (see in 10.3.ref04 around lines 130-140 In G4FTFParameters.cc).
251 // For the target destr. it's multipled by the number of target nucleons (12 for Carbon).
252 // In 10.3.p01 it's set to 1. FLAT OUT for both projectile & target, no multiplications, etc.
253 // Now, make default at 1.
255 double fNuclearTgtDestructP1; // Make D=1. as in 10.3.p01
257 double fNuclearProjDestructP2; // D=4.0
258 double fNuclearProjDestructP3; // D=2.1
259 double fNuclearTgtDestructP2; // D=4.0
260 double fNuclearTgtDestructP3; // D=2.1
261 //
262 double fPt2NuclearDestructP1; // D=0.035
263 double fPt2NuclearDestructP2; // D=0.04
264 double fPt2NuclearDestructP3; // D=4.0
265 double fPt2NuclearDestructP4; // D=2.5
266 // baryons... well, in fact also mesons...
267 double fR2ofNuclearDestruct; // D=1.5*fermi*fermi
269 double fDofNuclearDestruct; // D=0.3
270 // NOTE: this parameter has changed from 1. to 9. between 10.2 and 10.4.ref04 !!!
271 // ... but that's for baryons !
272 // ... while for mesons it's 1GeV**2
273 double fMaxPt2ofNuclearDestruct; // D=9GeV**2
274
275 private:
276
277 void Reset();
278};
279
280
282 public:
283 // ctor
285};
286
288
289 public:
290
291 // ctor
293
294};
295
297
298 public:
299
300 // ctor
302
303};
304
306 public:
309
310 void InitForInteraction( const G4ParticleDefinition* , G4int theA, G4int theZ, G4double s );
311
312 // Set geometrical parameteres
313 void SethNcmsEnergy( const G4double s );
314 void SetTotalCrossSection( const G4double Xtotal );
315 void SetElastisCrossSection( const G4double Xelastic );
316 void SetInelasticCrossSection( const G4double Xinelastic );
317 void SetProbabilityOfElasticScatt( const G4double Xtotal, const G4double Xelastic );
318 void SetProbabilityOfElasticScatt( const G4double aValue );
319 void SetProbabilityOfAnnihilation( const G4double aValue );
320 void SetRadiusOfHNinteractions2( const G4double Radius2 );
321
322 void SetSlope( const G4double Slope );
323 void SetGamma0( const G4double Gamma0 );
324 G4double GammaElastic( const G4double impactsquare );
325
326 // Set parameters of elastic scattering
328
329 // Set parameters of excitations
330 void SetParams( const G4int ProcN,
331 const G4double A1, const G4double B1, const G4double A2, const G4double B2,
332 const G4double A3, const G4double Atop, const G4double Ymin );
333
334 void SetDeltaProbAtQuarkExchange( const G4double aValue );
335 void SetProbOfSameQuarkExchange( const G4double aValue );
336
337 void SetProjMinDiffMass( const G4double aValue );
338 void SetProjMinNonDiffMass( const G4double aValue );
339 //void SetProbabilityOfProjDiff( const G4double aValue );
340 void SetProbLogDistrPrD( const G4double aValue );
341
342 void SetTarMinDiffMass( const G4double aValue );
343 void SetTarMinNonDiffMass( const G4double aValue );
344 //void SetProbabilityOfTarDiff( const G4double aValue );
345
346 void SetAveragePt2( const G4double aValue );
347 void SetProbLogDistr( const G4double aValue );
348
349 // Set parameters of a string kink
350 void SetPt2Kink( const G4double aValue );
351 void SetQuarkProbabilitiesAtGluonSplitUp( const G4double Puubar, const G4double Pddbar,
352 const G4double Pssbar );
353
354 // Set parameters of nuclear destruction
355 void SetMaxNumberOfCollisions( const G4double aValue, const G4double bValue );
356 void SetProbOfInteraction( const G4double aValue );
357
358 void SetCofNuclearDestructionPr( const G4double aValue );
359 void SetCofNuclearDestruction( const G4double aValue );
360 void SetR2ofNuclearDestruction( const G4double aValue );
361
363
364 void SetDofNuclearDestruction( const G4double aValue );
365 void SetPt2ofNuclearDestruction( const G4double aValue );
366 void SetMaxPt2ofNuclearDestruction( const G4double aValue );
367
368 // Get geometrical parameteres
372
373 G4double GetProbabilityOfInteraction( const G4double impactsquare );
374 G4double GetInelasticProbability( const G4double impactsquare );
378
379 // Get parameters of elastic scattering
381
382 // Get parameters of excitations
383 G4double GetProcProb( const G4int ProcN, const G4double y );
384
387
391
394
397
398 // Get parameters of a string kink
400 std::vector< G4double > GetQuarkProbabilitiesAtGluonSplitUp();
401
402 // Get parameters of nuclear destruction
405
409
411
415
416 // JVY, July 31, 2017: Is there any reason for NOT making
417 // all the members data private ???
418 //
419 // private:
420
421 // Initial energy of hN interactions
422 G4double FTFhNcmsEnergy; // Initial hN CMS energy
423
424 // Geometrical parameteres
425 G4double FTFXtotal; // Total X in mb
426 G4double FTFXelastic; // Elastic X in mb
427 G4double FTFXinelastic; // Inelastic X in mb
428 G4double FTFXannihilation; // Annihilation X in mb
432 G4double FTFSlope; // in fm^-1
435
436 // Parameters of excitations
438
441
447
450
451 // Parameters of kink
453 std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp;
454
455 // Parameters of nuclear destruction
458
459 G4double CofNuclearDestructionPr; // Cnd of nuclear destruction of projectile nucleus
460 G4double CofNuclearDestruction; // Cnd of nuclear destruction
462
464
465 G4double DofNuclearDestruction; // Dispersion for momentum sampling
468
469 private:
470 G4LundStringFragmentation* StringMass;
471 G4double GetMinMass( const G4ParticleDefinition* aParticle );
472
473 void Reset();
474
475 // JVY, July 31, 2017: encapsulates (current set of) parameters for the baryon projectile
476 //
477 G4FTFParamCollBaryonProj fParCollBaryonProj;
478
479 // JVY, Feb 14, 2019: encapsulates (current set of) parameters for meson/pion (+/-/0) projectile
480 G4FTFParamCollMesonProj fParCollMesonProj;
481 G4FTFParamCollPionProj fParCollPionProj;
482
483 // Glauber-Gribov hN x-section
484 G4VComponentCrossSection* csGGinstance;
485};
486
487
488inline G4double G4FTFParameters::GammaElastic( const G4double impactsquare ) {
489 return ( FTFGamma0 * G4Exp( -FTFSlope * impactsquare ) );
490}
491
494}
495
496// Set geometrical parameteres
497
499 FTFXtotal = Xtotal;
500}
501
503 FTFXelastic = Xelastic;
504}
505
507 FTFXinelastic = Xinelastic;
508}
509
511 const G4double Xelastic ) {
512 if ( Xtotal == 0.0 ) {
514 } else {
515 ProbabilityOfElasticScatt = Xelastic / Xtotal;
516 }
517}
518
521}
522
525}
526
528 RadiusOfHNinteractions2 = Radius2;
529}
530
531inline void G4FTFParameters::SetSlope( const G4double Slope ) {
532 FTFSlope = 12.84 / Slope; // Slope is in GeV^-2, FTFSlope in fm^-2
533}
534
535inline void G4FTFParameters::SetGamma0( const G4double Gamma0 ) {
536 FTFGamma0 = Gamma0;
537}
538
539// Set parameters of elastic scattering
542}
543
544// Set parameters of excitations
545
546inline void G4FTFParameters::SetParams( const G4int ProcN,
547 const G4double A1, const G4double B1, const G4double A2,
548 const G4double B2, const G4double A3, const G4double Atop,
549 const G4double Ymin ) {
550 ProcParams[ProcN][0] = A1; ProcParams[ProcN][1] = B1;
551 ProcParams[ProcN][2] = A2; ProcParams[ProcN][3] = B2;
552 ProcParams[ProcN][4] = A3;
553 ProcParams[ProcN][5] = Atop; ProcParams[ProcN][6] = Ymin;
554}
555
558}
559
562}
563
565 ProjMinDiffMass = aValue*CLHEP::GeV;
566}
567
569 ProjMinNonDiffMass = aValue*CLHEP::GeV;
570}
571
572inline void G4FTFParameters::SetTarMinDiffMass( const G4double aValue ) {
573 TarMinDiffMass = aValue*CLHEP::GeV;
574}
575
577 TarMinNonDiffMass = aValue*CLHEP::GeV;
578}
579
580inline void G4FTFParameters::SetAveragePt2( const G4double aValue ) {
581 AveragePt2 = aValue*CLHEP::GeV*CLHEP::GeV;
582}
583
585 ProbLogDistrPrD = aValue;
586}
587
588inline void G4FTFParameters::SetProbLogDistr( const G4double aValue ) {
589 ProbLogDistr = aValue;
590}
591
592// Set parameters of a string kink
593
594inline void G4FTFParameters::SetPt2Kink( const G4double aValue ) {
595 Pt2kink = aValue;
596}
597
599 const G4double Pddbar,
600 const G4double Pssbar ) {
601 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar );
602 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar );
603 QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar + Pssbar );
604}
605
606// Set parameters of nuclear destruction
608 const G4double Pbound ) {
609 if ( Plab > Pbound ) {
610 MaxNumberOfCollisions = Plab/Pbound;
611 SetProbOfInteraction( -1.0 );
612 } else {
613 //MaxNumberOfCollisions = -1.0;
614 //SetProbOfInteraction( G4Exp( 0.25*(Plab-Pbound) ) );
616 SetProbOfInteraction( -1.0 );
617 }
618}
619
621 ProbOfInelInteraction = aValue;
622}
623
626}
627
629 CofNuclearDestruction = aValue;
630}
631
633 R2ofNuclearDestruction = aValue;
634}
635
638}
639
641 DofNuclearDestruction = aValue;
642}
643
646}
647
650}
651
652// Get geometrical parameteres
654 return FTFXtotal;
655}
656
658 return FTFXelastic;
659}
660
662 return FTFXinelastic;
663}
664
666 return FTFSlope;
667}
668
670 if ( RadiusOfHNinteractions2 > impactsquare ) {
671 return 1.0;
672 } else {
673 return 0.0;
674 }
675}
676
679}
680
682 G4double Gamma = GammaElastic( impactsquare );
683 return 2*Gamma - Gamma*Gamma;
684}
685
688}
689
690// Get parameters of elastic scattering
693}
694
695// Get parameters of excitations
696
699}
700
703}
704
706 return ProjMinDiffMass;
707}
708
710 return ProjMinNonDiffMass;
711}
712
714 return TarMinDiffMass;
715}
716
718 return TarMinNonDiffMass;
719}
720
722 return AveragePt2;
723}
724
726 return ProbLogDistrPrD;
727}
728
730 return ProbLogDistr;
731}
732
733// Get parameters of a string kink
734
736 return Pt2kink;
737}
738
741}
742
743// Get parameters of nuclear destruction
744
747}
748
751}
752
755}
756
759}
760
763}
764
767}
768
771}
772
775}
776
779}
780
781#endif
782
double S(double temp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double GetProc3B2() const
double GetProc4A3() const
double GetNuclearTgtDestructP3() const
double GetProc3B1() const
bool IsNuclearProjDestructP1_NBRNDEP() const
double GetProc0A2() const
double GetExciEnergyPerWoundedNucleon() const
double GetProc1Atop() const
double GetProc2A3() const
double GetNuclearProjDestructP2() const
double GetProbLogDistrPrD() const
double GetPt2NuclearDestructP4() const
double GetPt2NuclearDestructP3() const
double GetNuclearProjDestructP3() const
double GetProc4A2() const
double GetPt2NuclearDestructP1() const
double GetProc2Ymin() const
double GetProjMinNonDiffMass() const
double GetProc3Ymin() const
double GetProc2A2() const
double GetProc4B1() const
double GetProc1A1() const
double GetProc0B2() const
double GetProbLogDistr() const
double GetProc1B2() const
double GetProc3Atop() const
double GetProc1A3() const
double GetProc1B1() const
double GetProc1Ymin() const
double GetProc2B1() const
double GetNuclearTgtDestructP2() const
double GetPt2NuclearDestructP2() const
double GetR2ofNuclearDestruct() const
double GetDeltaProbAtQuarkExchange() const
double GetDofNuclearDestruct() const
double GetProc0A3() const
bool IsProjDiffDissociation() const
double GetProc0Atop() const
double GetProc2B2() const
double GetProc2A1() const
double GetProc2Atop() const
double GetProc1A2() const
double GetProc4Atop() const
double GetProc3A1() const
double GetAveragePt2() const
double GetProc4Ymin() const
double GetProc0B1() const
double GetTgtMinNonDiffMass() const
virtual ~G4FTFParamCollection()
bool IsTgtDiffDissociation() const
double GetProbOfSameQuarkExchange() const
double GetProc4B2() const
double GetProc0Ymin() const
double GetProc3A2() const
double GetProc3A3() const
bool IsNuclearTgtDestructP1_ADEP() const
double GetMaxPt2ofNuclearDestruct() const
double GetProc4A1() const
double GetTgtMinDiffMass() const
double GetNuclearTgtDestructP1() const
double GetProjMinDiffMass() const
double GetProc0A1() const
double GetNuclearProjDestructP1() const
void SetTarMinNonDiffMass(const G4double aValue)
G4double R2ofNuclearDestruction
void SetProjMinDiffMass(const G4double aValue)
G4double GetProbLogDistrPrD()
G4double GetCofNuclearDestructionPr()
void SetTotalCrossSection(const G4double Xtotal)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
G4double GetProbabilityOfAnnihilation()
G4double ExcitationEnergyPerWoundedNucleon
G4double GetAveragePt2()
G4double GetMaxNumberOfCollisions()
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void SetElastisCrossSection(const G4double Xelastic)
void SetQuarkProbabilitiesAtGluonSplitUp(const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
G4double GetPt2ofNuclearDestruction()
G4double GetProbLogDistr()
G4double GetMaxPt2ofNuclearDestruction()
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
G4double GetProjMinNonDiffMass()
G4double ProcParams[5][7]
G4double GetProbabilityOfElasticScatt()
G4double ProbOfSameQuarkExchange
G4double GetAvaragePt2ofElasticScattering()
G4double TarMinNonDiffMass
G4double GetElasticCrossSection()
G4double ProbOfInelInteraction
void SetPt2Kink(const G4double aValue)
void SetSlope(const G4double Slope)
G4double GetTarMinDiffMass()
G4double RadiusOfHNinteractions2
void SetDeltaProbAtQuarkExchange(const G4double aValue)
G4double GetInelasticProbability(const G4double impactsquare)
G4double MaxPt2ofNuclearDestruction
G4double GetTarMinNonDiffMass()
G4double GetTotalCrossSection()
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
G4double DeltaProbAtQuarkExchange
G4double GetProcProb(const G4int ProcN, const G4double y)
G4double CofNuclearDestruction
void SetTarMinDiffMass(const G4double aValue)
void SetProjMinNonDiffMass(const G4double aValue)
void SetCofNuclearDestructionPr(const G4double aValue)
void SethNcmsEnergy(const G4double s)
G4double GetPt2Kink()
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp
void SetR2ofNuclearDestruction(const G4double aValue)
G4double GammaElastic(const G4double impactsquare)
void SetProbLogDistrPrD(const G4double aValue)
G4double CofNuclearDestructionPr
void SetGamma0(const G4double Gamma0)
void SetAveragePt2(const G4double aValue)
G4double GetDeltaProbAtQuarkExchange()
G4double DofNuclearDestruction
G4double GetInelasticCrossSection()
G4double GetExcitationEnergyPerWoundedNucleon()
void SetProbabilityOfAnnihilation(const G4double aValue)
G4double GetProbOfInteraction()
void SetProbOfInteraction(const G4double aValue)
void InitForInteraction(const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
void SetDofNuclearDestruction(const G4double aValue)
G4double ProbabilityOfAnnihilation
G4double GetProjMinDiffMass()
G4double GetProbOfSameQuarkExchange()
G4double ProbabilityOfElasticScatt
G4double GetDofNuclearDestruction()
void SetRadiusOfHNinteractions2(const G4double Radius2)
void SetParams(const G4int ProcN, const G4double A1, const G4double B1, const G4double A2, const G4double B2, const G4double A3, const G4double Atop, const G4double Ymin)
G4double Pt2ofNuclearDestruction
G4double ProbLogDistrPrD
G4double GetProbabilityOfInteraction(const G4double impactsquare)
void SetPt2ofNuclearDestruction(const G4double aValue)
G4double GetR2ofNuclearDestruction()
G4double GetCofNuclearDestruction()
void SetProbLogDistr(const G4double aValue)
void SetCofNuclearDestruction(const G4double aValue)
G4double ProjMinNonDiffMass
void SetMaxNumberOfCollisions(const G4double aValue, const G4double bValue)
std::vector< G4double > GetQuarkProbabilitiesAtGluonSplitUp()
G4double ProjMinDiffMass
G4double AvaragePt2ofElasticScattering
void SetProbOfSameQuarkExchange(const G4double aValue)
G4double MaxNumberOfCollisions
G4double FTFXannihilation
void SetInelasticCrossSection(const G4double Xinelastic)