Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4FTFParameters.hh
이 파일의 문서화 페이지로 가기
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.hh 107868 2017-12-07 14:45:15Z gcosmo $
28 // GEANT4 tag $Name: $
29 //
30 #ifndef G4FTFParameters_h
31 #define G4FTFParameters_h 1
32 
34 
35 #include "G4Proton.hh"
36 #include "G4Neutron.hh"
37 #include "G4ChipsComponentXS.hh"
38 
39 #include "G4Exp.hh"
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 
50  public:
51 
52  //dtor
53  virtual ~G4FTFParamCollection() {}
54 
55  // parameters of excitation
56  //
57  // Proc=0 --> Qexchg w/o excitation
58  //
59  double GetProc0A1() const { return fProc0A1; }
60  double GetProc0B1() const { return fProc0B1; }
61  double GetProc0A2() const { return fProc0A2; }
62  double GetProc0B2() const { return fProc0B2; }
63  double GetProc0A3() const { return fProc0A3; }
64  double GetProc0Atop() const { return fProc0Atop; }
65  double GetProc0Ymin() const { return fProc0Ymin; }
66  //
67  // Proc=1 --> Qexchg w/excitation
68  //
69  // Proc=2 & Proc=3 for the case ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )
70  // (diffraction dissociation)
71  //
74  //
75  double GetProc1A1() const { return fProc1A1; }
76  double GetProc1B1() const { return fProc1B1; }
77  double GetProc1A2() const { return fProc1A2; }
78  double GetProc1B2() const { return fProc1B2; }
79  double GetProc1A3() const { return fProc1A3; }
80  double GetProc1Atop() const { return fProc1Atop; }
81  double GetProc1Ymin() const { return fProc0Ymin; }
82  //
83  // Proc=4 --> Qexchg "w/additional multiplier" in excitation
84  //
85  double GetProc4A1() const { return fProc4A1; }
86  double GetProc4B1() const { return fProc4B1; }
87  double GetProc4A2() const { return fProc4A2; }
88  double GetProc4B2() const { return fProc4B2; }
89  double GetProc4A3() const { return fProc4A3; }
90  double GetProc4Atop() const { return fProc4Atop; }
91  double GetProc4Ymin() const { return fProc4Ymin; }
92  //
93  //
96  double GetProjMinDiffMass() const { return fProjMinDiffMass; }
97  double GetProjMinNonDiffMass() const { return fProjMinNonDiffMass; }
98  double GetTgtMinDiffMass() const { return fTgtMinDiffMass; }
99  double GetTgtMinNonDiffMass() const { return fTgtMinNonDiffMass; }
100  double GetAveragePt2() const { return fAveragePt2; }
101  double GetProbLogDistrPrD() const { return fProbLogDistrPrD; }
102  double GetProbLogDistr() const { return fProbLogDistr; }
103 
104  // NOTE (JVY): There is also the Pt2Kind parameter but for now it's set to 0., so we'll leave it aside
105 
106  // --> FIXME !!! --> void Get/SetBaryonMaxNumberOfCollisions( const double, const double ); // 1st is Plab, 2nd - D=2.
107  //
108  // NOTE (JVY): These parameters are COMMON among various projectiles !!!
109  //
121  double GetPt2NuclearDestructP4() const { return fPt2NuclearDestructP4; }
122  //
123  // separately for baryons, mesons, etc.
124  //
125  double GetR2ofNuclearDestruct() const { return fR2ofNuclearDestruct; }
127  double GetDofNuclearDestruct() const { return fDofNuclearDestruct; }
129 
130  protected:
131 
132  // ctor
134 
135  // parameters of excitation
136  //
137  //
138  // these are for Inelastic interactions, i.e. Xinelastic=(Xtotal-Xelastix)>0.
139  // for elastic, all the A's & B's, Atop & Ymin are zeros
140  // general formula: Pp = A1*exp(B1*Y) + A2*exp(B2*Y) + A3
141  // but if Y<Ymin, then Pp=max(0.,Atop)
142  // for details, see also G4FTFParameters::GetProcProb( ProcN, y )
143  //
144  // Proc=0 --> Qexchg w/o excitation
145  double fProc0A1; // D=13.71
146  double fProc0B1; // D=1.75
147  double fProc0A2; // D=-30.69 (or -214.5 as in Doc ?)
148  double fProc0B2; // D=3. ( or 4. as in Doc ?)
149  double fProc0A3; // D=0.
150  double fProc0Atop; // D=1. ( or 0.5 as in Doc ?)
151  double fProc0Ymin; // D=0.93 (or 1.1 as in Doc ?)
152  // Proc=1 --> Qexchg w/excitation
153  double fProc1A1; // D=25.
154  double fProc1B1; // D=1.
155  double fProc1A2; // D=-50.34
156  double fProc1B2; // D=1.5
157  double fProc1A3; // D=0.
158  double fProc1Atop; // D=0.
159  double fProc1Ymin; // D=1.4
160  //
161  // NOTE: Proc #2 & 3 are projectile & target diffraction
162  // they have more complex definition of A1 & A2
163  // (see around line 540 or so)
164  // SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Projectile diffraction
165  // SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Target diffraction
166  //
167  // Also, for ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )
168  // projectile and/or target diffraction (dissociation) may be switched ON/OFF
171  //
172  // Proc=4 --> Qexchg w/additional multiplier in excitation
173  double fProc4A1; // D=0.6 (or 1. as in Doc ?)
174  double fProc4B1; // D=0.
175  double fProc4A2; // D=-1.2 (or -2.01 as in Doc ?)
176  double fProc4B2; // D=0.5
177  double fProc4A3; // D=0.
178  double fProc4Atop; // D=0.
179  double fProc4Ymin; // D=1.4
180  //
181  // parameters of participating baryon excitation
182  //
183  double fDeltaProbAtQuarkExchange; // D=0.
184  double fProbOfSameQuarkExchange; // D=0. if A<=26, otherwise D=1.
185  double fProjMinDiffMass; // projectile, D=1.16GeV
186  double fProjMinNonDiffMass; // projectile, D=1.16GeV
187  double fTgtMinDiffMass; // target, D=1.16GeV
188  double fTgtMinNonDiffMass; // target, D=1.16GeV
189  double fAveragePt2; // D=0.3GeV**2 ( or 0.15 as in the Doc ???)
190  double fProbLogDistrPrD; // D=0.6 (or 0.3 ???)
191  double fProbLogDistr; // D=0.6 (or 0.3 ???)
192 
193  // parameters of nuclear distruction
194  //
195  // NOTE (JVY): there're 3 cases here:
196  // * baryon projectile
197  // * anti-baryon projectile
198  // * meson projectile
199  //
200  // double fBaryonMaxNumberOfCollisions; // D=2.
201  // void SetBaryonProbOfInteraction( const double ); // ??? this is prob. of inelastic interaction
202  // that is set internally based on certain conditions...
203  // general (i.e. for used for baryons,anti-baryons, and mesons)
204  // NOTE: these parameters have stayed THE SAME for quite a while
205  double fNuclearProjDestructP1; // D=0.00481 in 10.3.ref04 !!!
206  // BUT !!! In 10.3.ref04 as well as in 10.2-seriesit's multiplied of AbsProjectileBaryonNumber
207  // which somehow is 0 for the proton projectile (see in 10.3.ref04 around lines 130-140 In G4FTFParameters.cc).
208  // For the target destr. it's multipled by the number of target nucleons (12 for Carbon).
209  // In 10.3.p01 it's set to 1. FLAT OUT for both projectile & target, no multiplications, etc.
210  // Now, make default at 1.
212  double fNuclearTgtDestructP1; // Make D=1. as in 10.3.p01
214  double fNuclearProjDestructP2; // D=4.0
215  double fNuclearProjDestructP3; // D=2.1
216  double fNuclearTgtDestructP2; // D=4.0
217  double fNuclearTgtDestructP3; // D=2.1
218  //
219  double fPt2NuclearDestructP1; // D=0.035
220  double fPt2NuclearDestructP2; // D=0.04
221  double fPt2NuclearDestructP3; // D=4.0
222  double fPt2NuclearDestructP4; // D=2.5
223  // baryons
224  double fR2ofNuclearDestruct; // D=1.5*fermi*fermi
225  double fExciEnergyPerWoundedNucleon; // D=40MeV
226  double fDofNuclearDestruct; // D=0.3
227  // NOTE: this parameter has changed from 1. to 9. between 10.2 and 10.4.ref04 !!!
228  double fMaxPt2ofNuclearDestruct; // D=9GeV**2
229 
230  private:
231 
232  void Reset();
233 
234 };
235 
237 
238  public:
239 
240  // ctor
242 
243 };
244 
245 
247  public:
248  G4FTFParameters();
250 
251  void InitForInteraction( const G4ParticleDefinition* , G4int theA, G4int theZ, G4double s );
252 
253  // Set geometrical parameteres
254  void SethNcmsEnergy( const G4double s );
255  void SetTotalCrossSection( const G4double Xtotal );
256  void SetElastisCrossSection( const G4double Xelastic );
257  void SetInelasticCrossSection( const G4double Xinelastic );
258  void SetProbabilityOfElasticScatt( const G4double Xtotal, const G4double Xelastic );
259  void SetProbabilityOfElasticScatt( const G4double aValue );
260  void SetProbabilityOfAnnihilation( const G4double aValue );
261  void SetRadiusOfHNinteractions2( const G4double Radius2 );
262 
263  void SetSlope( const G4double Slope );
264  void SetGamma0( const G4double Gamma0 );
265  G4double GammaElastic( const G4double impactsquare );
266 
267  // Set parameters of elastic scattering
268  void SetAvaragePt2ofElasticScattering( const G4double aPt2 );
269 
270  // Set parameters of excitations
271  void SetParams( const G4int ProcN,
272  const G4double A1, const G4double B1, const G4double A2, const G4double B2,
273  const G4double A3, const G4double Atop, const G4double Ymin );
274 
275  void SetDeltaProbAtQuarkExchange( const G4double aValue );
276  void SetProbOfSameQuarkExchange( const G4double aValue );
277 
278  void SetProjMinDiffMass( const G4double aValue );
279  void SetProjMinNonDiffMass( const G4double aValue );
280  //void SetProbabilityOfProjDiff( const G4double aValue );
281  void SetProbLogDistrPrD( const G4double aValue );
282 
283  void SetTarMinDiffMass( const G4double aValue );
284  void SetTarMinNonDiffMass( const G4double aValue );
285  //void SetProbabilityOfTarDiff( const G4double aValue );
286 
287  void SetAveragePt2( const G4double aValue );
288  void SetProbLogDistr( const G4double aValue );
289 
290  // Set parameters of a string kink
291  void SetPt2Kink( const G4double aValue );
292  void SetQuarkProbabilitiesAtGluonSplitUp( const G4double Puubar, const G4double Pddbar,
293  const G4double Pssbar );
294 
295  // Set parameters of nuclear destruction
296  void SetMaxNumberOfCollisions( const G4double aValue, const G4double bValue );
297  void SetProbOfInteraction( const G4double aValue );
298 
299  void SetCofNuclearDestructionPr( const G4double aValue );
300  void SetCofNuclearDestruction( const G4double aValue );
301  void SetR2ofNuclearDestruction( const G4double aValue );
302 
303  void SetExcitationEnergyPerWoundedNucleon( const G4double aValue );
304 
305  void SetDofNuclearDestruction( const G4double aValue );
306  void SetPt2ofNuclearDestruction( const G4double aValue );
307  void SetMaxPt2ofNuclearDestruction( const G4double aValue );
308 
309  // Get geometrical parameteres
313 
314  G4double GetProbabilityOfInteraction( const G4double impactsquare );
315  G4double GetInelasticProbability( const G4double impactsquare );
317  G4double GetSlope();
319 
320  // Get parameters of elastic scattering
322 
323  // Get parameters of excitations
324  G4double GetProcProb( const G4int ProcN, const G4double y );
325 
328 
332 
335 
338 
339  // Get parameters of a string kink
341  std::vector< G4double > GetQuarkProbabilitiesAtGluonSplitUp();
342 
343  // Get parameters of nuclear destruction
346 
350 
352 
356 
357  // JVY, July 31, 2017: Is there any reason for NOT making
358  // all the members data private ???
359  //
360  // private:
361 
362  // Initial energy of hN interactions
363  G4double FTFhNcmsEnergy; // Initial hN CMS energy
364 
365  // hN cross section manager
367 
368  // Geometrical parameteres
369  G4double FTFXtotal; // Total X in mb
370  G4double FTFXelastic; // Elastic X in mb
371  G4double FTFXinelastic; // Inelastic X in mb
372  G4double FTFXannihilation; // Annihilation X in mb
373  G4double ProbabilityOfAnnihilation; // Xannih/Xinelast
375  G4double RadiusOfHNinteractions2; // Xtot/pi, in fm^2
376  G4double FTFSlope; // in fm^-1
379 
380  // Parameters of excitations
382 
385 
391 
394 
395  // Parameters of kink
397  std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp;
398 
399  // Parameters of nuclear destruction
402 
403  G4double CofNuclearDestructionPr; // Cnd of nuclear destruction of projectile nucleus
404  G4double CofNuclearDestruction; // Cnd of nuclear destruction
406 
408 
409  G4double DofNuclearDestruction; // D for momentum sampling
412 
413  private:
414 
415  void Reset();
416 
417  // JVY, July 31, 2017: encapsulates (current set of) parameters for the baryon projectile
418  //
420 
421  // G4-MT changes
422  private:
425 };
426 
427 
428 inline G4double G4FTFParameters::GammaElastic( const G4double impactsquare ) {
429  return ( FTFGamma0 * G4Exp( -FTFSlope * impactsquare ) );
430 }
431 
433  FTFhNcmsEnergy = S;
434 }
435 
436 // Set geometrical parameteres
437 
438 inline void G4FTFParameters::SetTotalCrossSection( const G4double Xtotal ) {
439  FTFXtotal = Xtotal;
440 }
441 
442 inline void G4FTFParameters::SetElastisCrossSection( const G4double Xelastic ) {
443  FTFXelastic = Xelastic;
444 }
445 
446 inline void G4FTFParameters::SetInelasticCrossSection( const G4double Xinelastic ) {
447  FTFXinelastic = Xinelastic;
448 }
449 
451  const G4double Xelastic ) {
452  if ( Xtotal == 0.0 ) {
454  } else {
455  ProbabilityOfElasticScatt = Xelastic / Xtotal;
456  }
457 }
458 
460  ProbabilityOfElasticScatt = aValue;
461 }
462 
464  ProbabilityOfAnnihilation = aValue;
465 }
466 
468  RadiusOfHNinteractions2 = Radius2;
469 }
470 
471 inline void G4FTFParameters::SetSlope( const G4double Slope ) {
472  FTFSlope = 12.84 / Slope; // Slope is in GeV^-2, FTFSlope in fm^-2
473 }
474 
475 inline void G4FTFParameters::SetGamma0( const G4double Gamma0 ) {
476  FTFGamma0 = Gamma0;
477 }
478 
479 // Set parameters of elastic scattering
482 }
483 
484 // Set parameters of excitations
485 
486 inline void G4FTFParameters::SetParams( const G4int ProcN,
487  const G4double A1, const G4double B1, const G4double A2,
488  const G4double B2, const G4double A3, const G4double Atop,
489  const G4double Ymin ) {
490  ProcParams[ProcN][0] = A1; ProcParams[ProcN][1] = B1;
491  ProcParams[ProcN][2] = A2; ProcParams[ProcN][3] = B2;
492  ProcParams[ProcN][4] = A3;
493  ProcParams[ProcN][5] = Atop; ProcParams[ProcN][6] = Ymin;
494 }
495 
497  DeltaProbAtQuarkExchange = aValue;
498 }
499 
501  ProbOfSameQuarkExchange = aValue;
502 }
503 
504 inline void G4FTFParameters::SetProjMinDiffMass( const G4double aValue ) {
505  ProjMinDiffMass = aValue*CLHEP::GeV;
506 }
507 
510 }
511 
512 inline void G4FTFParameters::SetTarMinDiffMass( const G4double aValue ) {
513  TarMinDiffMass = aValue*CLHEP::GeV;
514 }
515 
516 inline void G4FTFParameters::SetTarMinNonDiffMass( const G4double aValue ) {
517  TarMinNonDiffMass = aValue*CLHEP::GeV;
518 }
519 
520 inline void G4FTFParameters::SetAveragePt2( const G4double aValue ) {
522 }
523 
524 inline void G4FTFParameters::SetProbLogDistrPrD( const G4double aValue ) {
525  ProbLogDistrPrD = aValue;
526 }
527 
528 inline void G4FTFParameters::SetProbLogDistr( const G4double aValue ) {
529  ProbLogDistr = aValue;
530 }
531 
532 // Set parameters of a string kink
533 
534 inline void G4FTFParameters::SetPt2Kink( const G4double aValue ) {
535  Pt2kink = aValue;
536 }
537 
539  const G4double Pddbar,
540  const G4double Pssbar ) {
541  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar );
542  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar );
543  QuarkProbabilitiesAtGluonSplitUp.push_back( Puubar + Pddbar + Pssbar );
544 }
545 
546 // Set parameters of nuclear destruction
548  const G4double Pbound ) {
549  if ( Plab > Pbound ) {
550  MaxNumberOfCollisions = Plab/Pbound;
551  SetProbOfInteraction( -1.0 );
552  } else {
553  //MaxNumberOfCollisions = -1.0;
554  //SetProbOfInteraction( G4Exp( 0.25*(Plab-Pbound) ) );
556  SetProbOfInteraction( -1.0 );
557  }
558 }
559 
560 inline void G4FTFParameters::SetProbOfInteraction( const G4double aValue ) {
561  ProbOfInelInteraction = aValue;
562 }
563 
565  CofNuclearDestructionPr = aValue;
566 }
567 
569  CofNuclearDestruction = aValue;
570 }
571 
573  R2ofNuclearDestruction = aValue;
574 }
575 
578 }
579 
581  DofNuclearDestruction = aValue;
582 }
583 
585  Pt2ofNuclearDestruction = aValue;
586 }
587 
590 }
591 
592 // Get geometrical parameteres
594  return FTFXtotal;
595 }
596 
598  return FTFXelastic;
599 }
600 
602  return FTFXinelastic;
603 }
604 
606  return FTFSlope;
607 }
608 
610  if ( RadiusOfHNinteractions2 > impactsquare ) {
611  return 1.0;
612  } else {
613  return 0.0;
614  }
615 }
616 
619 }
620 
622  G4double Gamma = GammaElastic( impactsquare );
623  return 2*Gamma - Gamma*Gamma;
624 }
625 
628 }
629 
630 // Get parameters of elastic scattering
633 }
634 
635 // Get parameters of excitations
636 
639 }
640 
643 }
644 
646  return ProjMinDiffMass;
647 }
648 
650  return ProjMinNonDiffMass;
651 }
652 
654  return TarMinDiffMass;
655 }
656 
658  return TarMinNonDiffMass;
659 }
660 
662  return AveragePt2;
663 }
664 
666  return ProbLogDistrPrD;
667 }
668 
670  return ProbLogDistr;
671 }
672 
673 // Get parameters of a string kink
674 
676  return Pt2kink;
677 }
678 
681 }
682 
683 // Get parameters of nuclear destruction
684 
686  return MaxNumberOfCollisions;
687 }
688 
690  return ProbOfInelInteraction;
691 }
692 
695 }
696 
698  return CofNuclearDestruction;
699 }
700 
702  return R2ofNuclearDestruction;
703 }
704 
707 }
708 
710  return DofNuclearDestruction;
711 }
712 
715 }
716 
719 }
720 
721 #endif
722 
G4double GetProcProb(const G4int ProcN, const G4double y)
G4double GetMaxPt2ofNuclearDestruction()
void SetMaxNumberOfCollisions(const G4double aValue, const G4double bValue)
void SetSlope(const G4double Slope)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
double GetProc0A2() const
G4double GetProbOfSameQuarkExchange()
void SetCofNuclearDestructionPr(const G4double aValue)
G4double R2ofNuclearDestruction
bool IsProjDiffDissociation() const
double GetProc1B1() const
G4double GetAvaragePt2ofElasticScattering()
void SetElastisCrossSection(const G4double Xelastic)
G4double GetProjMinDiffMass()
double GetProc1Atop() const
void SetPt2Kink(const G4double aValue)
G4double GetInelasticCrossSection()
G4double GetCofNuclearDestruction()
G4double ProbOfInelInteraction
G4double GetTotalCrossSection()
double GetProbOfSameQuarkExchange() const
G4double GetCofNuclearDestructionPr()
void SetProjMinNonDiffMass(const G4double aValue)
double GetProc0B2() const
double GetProc0A1() const
G4double GetProbabilityOfAnnihilation()
bool IsTgtDiffDissociation() const
G4double DofNuclearDestruction
void SethNcmsEnergy(const G4double s)
Float_t y
Definition: compare.C:6
double GetPt2NuclearDestructP4() const
G4double GetTarMinNonDiffMass()
G4double GetPt2ofNuclearDestruction()
void SetCofNuclearDestruction(const G4double aValue)
double GetNuclearProjDestructP1() const
double GetProc4A1() const
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
double GetPt2NuclearDestructP3() const
G4double GetProbLogDistrPrD()
G4double ProbLogDistrPrD
double GetProc1A1() const
G4double GetProbLogDistr()
static G4ThreadLocal G4ChipsComponentXS * chipsComponentXSinstance
G4double AvaragePt2ofElasticScattering
double GetProc4A3() const
G4double MaxNumberOfCollisions
double GetNuclearTgtDestructP1() const
std::vector< G4double > QuarkProbabilitiesAtGluonSplitUp
void InitForInteraction(const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
G4double GetDeltaProbAtQuarkExchange()
double GetNuclearTgtDestructP2() const
double GetProbLogDistrPrD() const
#define G4ThreadLocal
Definition: tls.hh:69
double GetPt2NuclearDestructP1() const
G4double GetProbabilityOfElasticScatt()
double S(double temp)
void SetTarMinDiffMass(const G4double aValue)
static G4ThreadLocal bool chipsComponentXSisInitialized
G4double ProbabilityOfAnnihilation
int Ymin
void SetR2ofNuclearDestruction(const G4double aValue)
void SetProbLogDistrPrD(const G4double aValue)
double GetProc1Ymin() const
const XML_Char * s
Definition: expat.h:262
double GetProc0A3() const
G4double GetInelasticProbability(const G4double impactsquare)
double G4double
Definition: G4Types.hh:76
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
G4double GetAveragePt2()
G4double GetProjMinNonDiffMass()
void SetProbOfSameQuarkExchange(const G4double aValue)
double GetProc4Ymin() const
double GetProc4Atop() const
double GetProc0Atop() const
void SetProbOfInteraction(const G4double aValue)
double GetProc4A2() const
bool IsNuclearProjDestructP1_NBRNDEP() const
void SetTotalCrossSection(const G4double Xtotal)
double GetProjMinNonDiffMass() const
double GetProbLogDistr() const
double GetTgtMinNonDiffMass() const
void SetRadiusOfHNinteractions2(const G4double Radius2)
G4double CofNuclearDestruction
double GetTgtMinDiffMass() const
double GetProc4B1() const
G4double GetMaxNumberOfCollisions()
G4double ExcitationEnergyPerWoundedNucleon
double GetProc0Ymin() const
double GetProc1A2() const
double GetDeltaProbAtQuarkExchange() const
G4double TarMinNonDiffMass
G4double ProcParams[5][7]
G4double GetProbabilityOfInteraction(const G4double impactsquare)
G4double MaxPt2ofNuclearDestruction
void SetProbabilityOfAnnihilation(const G4double aValue)
void SetProbLogDistr(const G4double aValue)
G4double ProbOfSameQuarkExchange
double GetNuclearProjDestructP3() const
G4double ProbabilityOfElasticScatt
G4double ProjMinNonDiffMass
G4double RadiusOfHNinteractions2
bool IsNuclearTgtDestructP1_ADEP() const
double GetProc1A3() const
int G4int
Definition: G4Types.hh:78
double GetDofNuclearDestruct() const
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
void SetQuarkProbabilitiesAtGluonSplitUp(const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
G4double Pt2ofNuclearDestruction
G4double GetProbOfInteraction()
double GetProc1B2() const
virtual ~G4FTFParamCollection()
double GetR2ofNuclearDestruct() const
G4double ProjMinDiffMass
G4double CofNuclearDestructionPr
void SetProjMinDiffMass(const G4double aValue)
void SetDofNuclearDestruction(const G4double aValue)
double GetNuclearProjDestructP2() const
G4ChipsComponentXS * FTFxsManager
void SetAveragePt2(const G4double aValue)
G4double GetElasticCrossSection()
void SetPt2ofNuclearDestruction(const G4double aValue)
double GetMaxPt2ofNuclearDestruct() const
G4double GetPt2Kink()
double GetNuclearTgtDestructP3() const
G4double DeltaProbAtQuarkExchange
double GetProc0B1() const
G4double GetExcitationEnergyPerWoundedNucleon()
double GetExciEnergyPerWoundedNucleon() const
double GetPt2NuclearDestructP2() const
void SetTarMinNonDiffMass(const G4double aValue)
G4double GetR2ofNuclearDestruction()
void SetInelasticCrossSection(const G4double Xinelastic)
G4FTFParamCollBaryonProj fParCollBaryonProj
G4double FTFXannihilation
G4double GammaElastic(const G4double impactsquare)
void SetGamma0(const G4double Gamma0)
std::vector< G4double > GetQuarkProbabilitiesAtGluonSplitUp()
static constexpr double GeV
void SetDeltaProbAtQuarkExchange(const G4double aValue)
double GetProjMinDiffMass() const
double GetProc4B2() const
G4double GetTarMinDiffMass()
double GetAveragePt2() const
G4double GetDofNuclearDestruction()
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)