Geant4 Cross Reference |
1 // 1 // 2 // ******************************************* 2 // ******************************************************************** 3 // * License and Disclaimer 3 // * License and Disclaimer * 4 // * 4 // * * 5 // * The Geant4 software is copyright of th 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/ 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. 9 // * include a list of copyright holders. * 10 // * 10 // * * 11 // * Neither the authors of this software syst 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing fin 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warran 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assum 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitatio 16 // * for the full disclaimer and the limitation of liability. * 17 // * 17 // * * 18 // * This code implementation is the result 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboratio 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distri 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you ag 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publicati 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Sof 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************* 24 // ******************************************************************** 25 // 25 // 26 // 26 // 27 // 27 // 28 // Hadronic Process: Nuclear De-excitations 28 // Hadronic Process: Nuclear De-excitations 29 // by V. Lara (Oct 1998) 29 // by V. Lara (Oct 1998) 30 30 31 #ifndef G4CompetitiveFission_h 31 #ifndef G4CompetitiveFission_h 32 #define G4CompetitiveFission_h 1 32 #define G4CompetitiveFission_h 1 33 33 34 #include "G4VEvaporationChannel.hh" 34 #include "G4VEvaporationChannel.hh" 35 #include "G4Fragment.hh" 35 #include "G4Fragment.hh" 36 #include "G4VEmissionProbability.hh" 36 #include "G4VEmissionProbability.hh" 37 #include "G4FissionParameters.hh" 37 #include "G4FissionParameters.hh" 38 #include <CLHEP/Units/SystemOfUnits.h> 38 #include <CLHEP/Units/SystemOfUnits.h> 39 #include "G4Exp.hh" 39 #include "G4Exp.hh" 40 40 41 class G4VFissionBarrier; 41 class G4VFissionBarrier; 42 class G4VEmissionProbability; 42 class G4VEmissionProbability; 43 class G4VLevelDensityParameter; 43 class G4VLevelDensityParameter; 44 class G4PairingCorrection; 44 class G4PairingCorrection; 45 45 46 class G4CompetitiveFission : public G4VEvapora 46 class G4CompetitiveFission : public G4VEvaporationChannel 47 { 47 { 48 public: 48 public: 49 49 50 G4CompetitiveFission(); << 50 explicit G4CompetitiveFission(); 51 ~G4CompetitiveFission() override; 51 ~G4CompetitiveFission() override; 52 52 53 void Initialise() override; << 54 << 55 G4Fragment* EmittedFragment(G4Fragment* theN 53 G4Fragment* EmittedFragment(G4Fragment* theNucleus) override; 56 54 57 G4double GetEmissionProbability(G4Fragment* 55 G4double GetEmissionProbability(G4Fragment* theNucleus) override; 58 56 59 void SetFissionBarrier(G4VFissionBarrier * a 57 void SetFissionBarrier(G4VFissionBarrier * aBarrier); 60 58 61 void SetEmissionStrategy(G4VEmissionProbabil 59 void SetEmissionStrategy(G4VEmissionProbability * aFissionProb); 62 60 63 void SetLevelDensityParameter(G4VLevelDensit 61 void SetLevelDensityParameter(G4VLevelDensityParameter * aLevelDensity); 64 62 65 inline G4double GetFissionBarrier(void) cons 63 inline G4double GetFissionBarrier(void) const; 66 64 67 inline G4double GetLevelDensityParameter(voi 65 inline G4double GetLevelDensityParameter(void) const; 68 66 69 inline G4double GetMaximalKineticEnergy(void 67 inline G4double GetMaximalKineticEnergy(void) const; 70 68 71 G4CompetitiveFission(const G4CompetitiveFiss << 72 const G4CompetitiveFission & operator=(const << 73 G4bool operator==(const G4CompetitiveFission << 74 G4bool operator!=(const G4CompetitiveFission << 75 << 76 private: 69 private: 77 70 78 // Sample AtomicNumber of Fission products 71 // Sample AtomicNumber of Fission products 79 G4int FissionAtomicNumber(G4int A); 72 G4int FissionAtomicNumber(G4int A); 80 73 81 G4double MassDistribution(G4double x, G4int 74 G4double MassDistribution(G4double x, G4int A); 82 75 83 // Sample Charge of fission products 76 // Sample Charge of fission products 84 G4int FissionCharge(G4int A, G4int Z, G4doub 77 G4int FissionCharge(G4int A, G4int Z, G4double Af); 85 78 86 // Sample Kinetic energy of fission products 79 // Sample Kinetic energy of fission products 87 G4double FissionKineticEnergy(G4int A, G4int 80 G4double FissionKineticEnergy(G4int A, G4int Z, 88 G4int Af1, G4int Zf1, 81 G4int Af1, G4int Zf1, 89 G4int Af2, G4int Zf2, 82 G4int Af2, G4int Zf2, 90 G4double U, G4double Tmax); 83 G4double U, G4double Tmax); 91 84 92 inline G4double Ratio(G4double A, G4double A 85 inline G4double Ratio(G4double A, G4double A11, 93 G4double B1, G4double 86 G4double B1, G4double A00) const; 94 87 95 inline G4double SymmetricRatio(G4int A, G4do 88 inline G4double SymmetricRatio(G4int A, G4double A11) const; 96 89 97 inline G4double AsymmetricRatio(G4int A, G4d 90 inline G4double AsymmetricRatio(G4int A, G4double A11) const; 98 91 99 inline G4double LocalExp(G4double x) const; 92 inline G4double LocalExp(G4double x) const; 100 93 >> 94 G4CompetitiveFission(const G4CompetitiveFission &right); >> 95 const G4CompetitiveFission & operator=(const G4CompetitiveFission &right); >> 96 G4bool operator==(const G4CompetitiveFission &right) const; >> 97 G4bool operator!=(const G4CompetitiveFission &right) const; >> 98 101 // Maximal Kinetic Energy that can be carrie 99 // Maximal Kinetic Energy that can be carried by fragment 102 G4double maxKineticEnergy{0.0}; << 100 G4double maxKineticEnergy; 103 G4double fissionBarrier{0.0}; << 101 G4double fissionBarrier; 104 G4double fissionProbability{0.0}; << 102 G4double fissionProbability; 105 G4double fFactor{1.0}; << 106 103 107 // For Fission barrier 104 // For Fission barrier 108 G4VFissionBarrier* theFissionBarrierPtr; 105 G4VFissionBarrier* theFissionBarrierPtr; 109 106 110 // For Fission probability emission 107 // For Fission probability emission 111 G4VEmissionProbability* theFissionProbabilit 108 G4VEmissionProbability* theFissionProbabilityPtr; 112 109 113 // For Level Density calculation 110 // For Level Density calculation 114 G4VLevelDensityParameter* theLevelDensityPtr 111 G4VLevelDensityParameter* theLevelDensityPtr; 115 G4PairingCorrection* pairingCorrection; 112 G4PairingCorrection* pairingCorrection; 116 113 117 G4bool myOwnFissionProbability{true}; << 114 G4bool myOwnFissionProbability; 118 G4bool myOwnFissionBarrier{true}; << 115 G4bool myOwnFissionBarrier; 119 G4bool myOwnLevelDensity{true}; << 116 G4bool myOwnLevelDensity; 120 117 121 G4FissionParameters theParam; 118 G4FissionParameters theParam; 122 << 119 123 G4int theSecID; // Creator model ID for the << 124 G4bool isInitialised{false}; << 125 }; 120 }; 126 121 127 inline G4double G4CompetitiveFission::GetFissi 122 inline G4double G4CompetitiveFission::GetFissionBarrier(void) const 128 { 123 { 129 return fissionBarrier; 124 return fissionBarrier; 130 } 125 } 131 126 132 inline G4double G4CompetitiveFission::GetMaxim 127 inline G4double G4CompetitiveFission::GetMaximalKineticEnergy(void) const 133 { 128 { 134 return maxKineticEnergy; 129 return maxKineticEnergy; 135 } 130 } 136 131 137 inline 132 inline 138 G4double G4CompetitiveFission::Ratio(G4double 133 G4double G4CompetitiveFission::Ratio(G4double A, G4double A11, 139 G4double B1, G4double A00) const 134 G4double B1, G4double A00) const 140 { 135 { 141 G4double res; 136 G4double res; 142 if (A11 >= A*0.5 && A11 <= (A00+10.0)) { 137 if (A11 >= A*0.5 && A11 <= (A00+10.0)) { 143 G4double x = (A11-A00)/A; 138 G4double x = (A11-A00)/A; 144 res = 1.0 - B1*x*x; 139 res = 1.0 - B1*x*x; 145 } else { 140 } else { 146 G4double x = 10.0/A; 141 G4double x = 10.0/A; 147 res = 1.0 - B1*x*x - 2.0*x*B1*(A11-A00-10. 142 res = 1.0 - B1*x*x - 2.0*x*B1*(A11-A00-10.0)/A; 148 } 143 } 149 return res; 144 return res; 150 } 145 } 151 146 152 inline 147 inline 153 G4double G4CompetitiveFission::AsymmetricRatio 148 G4double G4CompetitiveFission::AsymmetricRatio(G4int A, G4double A11) const 154 { 149 { 155 return Ratio(G4double(A),A11,23.5,134.0); 150 return Ratio(G4double(A),A11,23.5,134.0); 156 } 151 } 157 152 158 inline 153 inline 159 G4double G4CompetitiveFission::SymmetricRatio( 154 G4double G4CompetitiveFission::SymmetricRatio(G4int A, G4double A11) const 160 { 155 { 161 G4double A0 = G4double(A); 156 G4double A0 = G4double(A); 162 return Ratio(A0,A11,5.32,A0*0.5); 157 return Ratio(A0,A11,5.32,A0*0.5); 163 } 158 } 164 159 165 inline G4double G4CompetitiveFission::LocalExp 160 inline G4double G4CompetitiveFission::LocalExp(G4double x) const 166 { 161 { 167 return (std::abs(x) < 8.) ? G4Exp(-0.5*x*x) 162 return (std::abs(x) < 8.) ? G4Exp(-0.5*x*x) : 0.0; 168 } 163 } 169 164 170 #endif 165 #endif 171 166 172 167 173 168