Geant4-11
G4StatMFFragment.cc
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// Hadronic Process: Nuclear De-excitations
29// by V. Lara
30
31#include "G4StatMFFragment.hh"
34#include "G4Pow.hh"
35
36// Copy constructor
38{
39 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::copy_constructor meant to not be accessible");
40}
41
42// Operators
43
46{
47 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator= meant to not be accessible");
48 return *this;
49}
50
52{
53// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator== meant to not be accessible");
54 return false;
55}
56
58{
59// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator!= meant to not be accessible");
60 return true;
61}
62
64{
65 G4double res = 0.0;
66 if (theZ >= 1) {
68 }
69 return res;
70}
71
73{
74 if (theA < 1 || theZ < 0 || theZ > theA) {
75 G4cout << "G4StatMFFragment::GetEnergy: A = " << theA
76 << ", Z = " << theZ << G4endl;
77 throw G4HadronicException(__FILE__, __LINE__,
78 "G4StatMFFragment::GetEnergy: Wrong values for A and Z!");
79 }
81
82 if (theA < 4) return BulkEnergy - GetCoulombEnergy();
83
84 G4double SurfaceEnergy;
85 if (G4StatMFParameters::DBetaDT(T) == 0.0) SurfaceEnergy = 0.0;
86 else SurfaceEnergy = 2.5*G4Pow::GetInstance()->Z23(theA)*T*T*
90
91 G4double ExchangeEnergy = theA*T*T/GetInvLevelDensity();
92 if (theA != 4) ExchangeEnergy += SurfaceEnergy;
93 return BulkEnergy + ExchangeEnergy - GetCoulombEnergy();
94}
95
97{
98 G4double res = 0.0;
99 if (theA > 1) {
100 res = G4StatMFParameters::GetEpsilon0()*(1.0+3.0/(theA - 1.0));
101 }
102 return res;
103}
104
106{
109 G4LorentzVector FourMomentum(_momentum,std::sqrt(_momentum.mag2()+(M+U)*(M+U)));
110 G4Fragment * theFragment = new G4Fragment(theA, theZ, FourMomentum);
111 return theFragment;
112}
113
115{
116 if (theA <= 3) return 0.0;
117
118 G4double BulkEnergy = theA*T*T/GetInvLevelDensity();
119
120 // if it is an alpha particle: done
121 if (theA == 4) return BulkEnergy;
122
123 // Term connected with surface energy
124 G4double SurfaceEnergy = 0.0;
126 if (std::abs(q) > 1.0e-20) {
127 SurfaceEnergy = 2.5*G4Pow::GetInstance()->Z23(theA)
129 }
130 return BulkEnergy + SurfaceEnergy;
131}
#define M(row, col)
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double mag2() const
static G4double GetMassExcess(const G4int A, const G4int Z)
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double Z23(G4int Z) const
Definition: G4Pow.hh:125
G4bool operator==(const G4StatMFFragment &right) const
G4ThreeVector _momentum
G4StatMFFragment & operator=(const G4StatMFFragment &right)
G4double GetInvLevelDensity(void) const
G4double CalcExcitationEnergy(const G4double T)
G4Fragment * GetFragment(const G4double T)
G4double GetNuclearMass(void)
G4double GetEnergy(const G4double T) const
G4bool operator!=(const G4StatMFFragment &right) const
G4double GetCoulombEnergy(void) const
static G4double DBetaDT(G4double T)
static G4double GetBeta0()
static G4double Beta(G4double T)
static G4double GetCoulomb()
static G4double GetCriticalTemp()
static G4double GetEpsilon0()