Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4CompetitiveFission.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: G4CompetitiveFission.hh 107060 2017-11-01 15:00:04Z gcosmo $
28 //
29 // Hadronic Process: Nuclear De-excitations
30 // by V. Lara (Oct 1998)
31 
32 #ifndef G4CompetitiveFission_h
33 #define G4CompetitiveFission_h 1
34 
35 #include "G4VEvaporationChannel.hh"
36 #include "G4Fragment.hh"
37 #include "G4VFissionBarrier.hh"
38 #include "G4FissionBarrier.hh"
40 #include "G4FissionProbability.hh"
43 #include "G4FissionParameters.hh"
44 #include "Randomize.hh"
46 
48 {
49 public:
50 
51  explicit G4CompetitiveFission();
52  virtual ~G4CompetitiveFission();
53 
54  virtual G4Fragment* EmittedFragment(G4Fragment* theNucleus);
55 
56  virtual G4double GetEmissionProbability(G4Fragment* theNucleus);
57 
58  inline void SetFissionBarrier(G4VFissionBarrier * aBarrier);
59 
60  inline void SetEmissionStrategy(G4VEmissionProbability * aFissionProb);
61 
62  inline void SetLevelDensityParameter(G4VLevelDensityParameter * aLevelDensity);
63 
64  inline G4double GetFissionBarrier(void) const;
65 
66  inline G4double GetLevelDensityParameter(void) const;
67 
68  inline G4double GetMaximalKineticEnergy(void) const;
69 
70 private:
71 
72  // Sample AtomicNumber of Fission products
74 
76 
77  // Sample Charge of fission products
79 
80  // Sample Kinetic energy of fission products
82  G4int Af1, G4int Zf1,
83  G4int Af2, G4int Zf2,
84  G4double U, G4double Tmax);
85 
87 
88  inline G4double SymmetricRatio(G4int A, G4double A11);
89 
90  inline G4double AsymmetricRatio(G4int A, G4double A11);
91 
93  const G4CompetitiveFission & operator=(const G4CompetitiveFission &right) = delete;
94  G4bool operator==(const G4CompetitiveFission &right) const = delete;
95  G4bool operator!=(const G4CompetitiveFission &right) const = delete;
96 
97  // Maximal Kinetic Energy that can be carried by fragment
99 
100  // For Fission barrier
104 
105  // For Fission probability emission
109 
110  // For Level Density calculation
114 
116 
118 
119 };
120 
122 {
124  theFissionBarrierPtr = aBarrier;
125  MyOwnFissionBarrier = false;
126 }
127 
128 inline void
130 {
132  theFissionProbabilityPtr = aFissionProb;
133  MyOwnFissionProbability = false;
134 }
135 
136 inline void
138 {
140  theLevelDensityPtr = aLevelDensity;
141  MyOwnLevelDensity = false;
142 }
143 
145 {
146  return FissionBarrier;
147 }
148 
150 {
151  return LevelDensityParameter;
152 }
153 
155 {
156  return MaximalKineticEnergy;
157 }
158 
159 inline
161  G4double B1, G4double A00)
162 {
163  G4double res;
164  if (A11 >= A*0.5 && A11 <= (A00+10.0)) {
165  G4double x = (A11-A00)/A;
166  res = 1.0 - B1*x*x;
167  } else {
168  G4double x = 10.0/A;
169  res = 1.0 - B1*x*x - 2.0*x*B1*(A11-A00-10.0)/A;
170  }
171  return res;
172 }
173 
174 inline
176 {
177  return Ratio(G4double(A),A11,23.5,134.0);
178 }
179 
180 inline
182 {
183  G4double A0 = G4double(A);
184  return Ratio(A0,A11,5.32,A0*0.5);
185 }
186 
187 #endif
188 
189 
Float_t x
Definition: compare.C:6
G4double AsymmetricRatio(G4int A, G4double A11)
void SetEmissionStrategy(G4VEmissionProbability *aFissionProb)
G4int FissionCharge(G4int A, G4int Z, G4double Af)
void SetFissionBarrier(G4VFissionBarrier *aBarrier)
#define A11
G4double GetLevelDensityParameter(void) const
G4FissionParameters theParam
Float_t Z
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
G4VEmissionProbability * theFissionProbabilityPtr
G4double FissionKineticEnergy(G4int A, G4int Z, G4int Af1, G4int Zf1, G4int Af2, G4int Zf2, G4double U, G4double Tmax)
void SetLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
double A(double temperature)
G4PairingCorrection * pairingCorrection
G4VLevelDensityParameter * theLevelDensityPtr
G4VFissionBarrier * theFissionBarrierPtr
const G4CompetitiveFission & operator=(const G4CompetitiveFission &right)=delete
int G4int
Definition: G4Types.hh:78
G4double GetMaximalKineticEnergy(void) const
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)
G4bool operator!=(const G4CompetitiveFission &right) const =delete
#define A00
G4double Ratio(G4double A, G4double A11, G4double B1, G4double A00)
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
G4int FissionAtomicNumber(G4int A)
G4double MassDistribution(G4double x, G4int A)
G4bool operator==(const G4CompetitiveFission &right) const =delete
G4double SymmetricRatio(G4int A, G4double A11)
G4double GetFissionBarrier(void) const