Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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// $Id$
28//
29// Hadronic Process: Nuclear De-excitations
30// by V. Lara
31
32#include "G4StatMFFragment.hh"
35
36// Copy constructor
38{
39 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::copy_constructor meant to not be accessable");
40}
41
42// Operators
43
44G4StatMFFragment & G4StatMFFragment::
45operator=(const G4StatMFFragment & )
46{
47 throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator= meant to not be accessable");
48 return *this;
49}
50
51
53{
54// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator== meant to not be accessable");
55 return false;
56}
57
58
60{
61// throw G4HadronicException(__FILE__, __LINE__, "G4StatMFFragment::operator!= meant to not be accessable");
62 return true;
63}
64
65
66
68{
69 if (theZ <= 0.1) return 0.0;
70 G4double Coulomb = (3./5.)*(elm_coupling*theZ*theZ)*
71 std::pow(1.0+G4StatMFParameters::GetKappaCoulomb(),1./3.)/
72 (G4StatMFParameters::Getr0()*std::pow(theA,1./3.));
73
74 return Coulomb;
75}
76
77
79{
80 if (theA < 1 || theZ < 0 || theZ > theA) {
81 G4cerr << "G4StatMFFragment::GetEnergy: A = " << theA
82 << ", Z = " << theZ << G4endl;
83 throw G4HadronicException(__FILE__, __LINE__,
84 "G4StatMFFragment::GetEnergy: Wrong values for A and Z!");
85 }
86 G4double BulkEnergy = G4NucleiProperties::GetMassExcess(static_cast<G4int>(theA),
87 static_cast<G4int>(theZ));
88
89 if (theA < 4) return BulkEnergy - GetCoulombEnergy();
90
91 G4double SurfaceEnergy;
92 if (G4StatMFParameters::DBetaDT(T) == 0.0) SurfaceEnergy = 0.0;
93 else SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*T*T*
97
98
99 G4double ExchangeEnergy = theA*T*T/GetInvLevelDensity();
100 if (theA != 4) ExchangeEnergy += SurfaceEnergy;
101
102 return BulkEnergy + ExchangeEnergy - GetCoulombEnergy();
103
104}
105
106
108{
109 // Calculate Inverse Density Level
110 // Epsilon0*(1 + 3 /(Af - 1))
111 if (theA == 1) return 0.0;
112 else return
113 G4StatMFParameters::GetEpsilon0()*(1.0+3.0/(theA - 1.0));
114}
115
116
117
119{
120 G4double U = CalcExcitationEnergy(T);
121
123
124 G4LorentzVector FourMomentum(_momentum,std::sqrt(_momentum.mag2()+(M+U)*(M+U)));
125
126 G4Fragment * theFragment = new G4Fragment(static_cast<G4int>(theA),static_cast<G4int>(theZ),FourMomentum);
127
128 return theFragment;
129}
130
131
132G4double G4StatMFFragment::CalcExcitationEnergy(const G4double T)
133{
134 if (theA <= 3) return 0.0;
135
136 G4double BulkEnergy = theA*T*T/GetInvLevelDensity();
137
138 // if it is an alpha particle: done
139 if (theA == 4) return BulkEnergy;
140
141 // Term connected with surface energy
142 G4double SurfaceEnergy = 0.0;
143 if (std::abs(G4StatMFParameters::DBetaDT(T)) > 1.0e-20)
144// SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*T*T*G4StatMFParameters::GetBeta0()/
145// (G4StatMFParameters::GetCriticalTemp()*G4StatMFParameters::GetCriticalTemp());
146 SurfaceEnergy = (5./2.)*std::pow(theA,2.0/3.0)*(G4StatMFParameters::Beta(T) -
148
149 return BulkEnergy + SurfaceEnergy;
150}
151
152
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
double mag2() const
static G4double GetMassExcess(const G4int A, const G4int Z)
G4bool operator==(const G4StatMFFragment &right) const
G4double GetInvLevelDensity(void) const
G4Fragment * GetFragment(const G4double T)
G4StatMFFragment(G4int anA, G4int aZ)
G4double GetNuclearMass(void)
G4double GetEnergy(const G4double T) const
G4bool operator!=(const G4StatMFFragment &right) const
G4double GetCoulombEnergy(void) const
static G4double Getr0()
static G4double GetEpsilon0()
static G4double GetKappaCoulomb()
static G4double GetBeta0()
static G4double GetCriticalTemp()
static G4double DBetaDT(const G4double T)
static G4double Beta(const G4double T)