Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4CompetitiveFission.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 // $Id: G4CompetitiveFission.hh 98739 2016-08-09 12:56:55Z 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 
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
73  G4int FissionAtomicNumber(G4int A);
74 
75  G4double MassDistribution(G4double x, G4int A);
76 
77  // Sample Charge of fission products
78  G4int FissionCharge(G4int A, G4int Z, G4double Af);
79 
80  // Sample Kinetic energy of fission products
81  G4double FissionKineticEnergy(G4int A, G4int Z,
82  G4int Af1, G4int Zf1,
83  G4int Af2, G4int Zf2,
84  G4double U, G4double Tmax);
85 
86  inline G4double Ratio(G4double A, G4double A11, G4double B1, G4double A00);
87 
88  inline G4double SymmetricRatio(G4int A, G4double A11);
89 
90  inline G4double AsymmetricRatio(G4int A, G4double A11);
91 
92  inline G4ThreeVector IsotropicVector(G4double Magnitude);
93 
95  const G4CompetitiveFission & operator=(const G4CompetitiveFission &right);
96  G4bool operator==(const G4CompetitiveFission &right) const;
97  G4bool operator!=(const G4CompetitiveFission &right) const;
98 
99  // Maximal Kinetic Energy that can be carried by fragment
100  G4double MaximalKineticEnergy;
101 
102  // For Fission barrier
103  G4VFissionBarrier * theFissionBarrierPtr;
104  G4double FissionBarrier;
105  G4bool MyOwnFissionBarrier;
106 
107  // For Fission probability emission
108  G4VEmissionProbability * theFissionProbabilityPtr;
109  G4double FissionProbability;
110  G4bool MyOwnFissionProbability;
111 
112  // For Level Density calculation
113  G4bool MyOwnLevelDensity;
114  G4VLevelDensityParameter * theLevelDensityPtr;
115  G4double LevelDensityParameter;
116 
117  G4PairingCorrection* pairingCorrection;
118 
119  G4FissionParameters theParam;
120 
121 };
122 
124 {
125  if (MyOwnFissionBarrier) delete theFissionBarrierPtr;
126  theFissionBarrierPtr = aBarrier;
127  MyOwnFissionBarrier = false;
128 }
129 
130 inline void
132 {
133  if (MyOwnFissionProbability) delete theFissionProbabilityPtr;
134  theFissionProbabilityPtr = aFissionProb;
135  MyOwnFissionProbability = false;
136 }
137 
138 inline void
140 {
141  if (MyOwnLevelDensity) delete theLevelDensityPtr;
142  theLevelDensityPtr = aLevelDensity;
143  MyOwnLevelDensity = false;
144 }
145 
147 {
148  return FissionBarrier;
149 }
150 
152 {
153  return LevelDensityParameter;
154 }
155 
157 {
158  return MaximalKineticEnergy;
159 }
160 
161 inline
162 G4double G4CompetitiveFission::Ratio(G4double A, G4double A11,
163  G4double B1, G4double A00)
164 {
165  G4double res;
166  if (A11 >= A*0.5 && A11 <= (A00+10.0)) {
167  G4double x = (A11-A00)/A;
168  res = 1.0 - B1*x*x;
169  } else {
170  G4double x = 10.0/A;
171  res = 1.0 - B1*x*x - 2.0*x*B1*(A11-A00-10.0)/A;
172  }
173  return res;
174 }
175 
176 inline
177 G4double G4CompetitiveFission::AsymmetricRatio(G4int A, G4double A11)
178 {
179  return Ratio(G4double(A),A11,23.5,134.0);
180 }
181 
182 inline
183 G4double G4CompetitiveFission::SymmetricRatio(G4int A, G4double A11)
184 {
185  G4double A0 = G4double(A);
186  return Ratio(A0,A11,5.32,A0*0.5);
187 }
188 
189 inline
190 G4ThreeVector G4CompetitiveFission::IsotropicVector(G4double Magnitude)
191 {
192  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
193  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
195  G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
196  Magnitude*std::sin(Phi)*SinTheta,
197  Magnitude*CosTheta);
198  return Vector;
199 }
200 
201 #endif
202 
203 
void SetFissionBarrier(G4VFissionBarrier *aBarrier)
#define A00
G4double GetLevelDensityParameter(void) const
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
void SetEmissionStrategy(G4VEmissionProbability *aFissionProb)
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
bool G4bool
Definition: G4Types.hh:79
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
#define A11
double G4double
Definition: G4Types.hh:76
void SetLevelDensityParameter(G4VLevelDensityParameter *aLevelDensity)
static constexpr double twopi
Definition: SystemOfUnits.h:55
G4double GetFissionBarrier(void) const
G4double GetMaximalKineticEnergy(void) const