ShellModel

Functions for shell-model calculations

NuclearToolkit.HbitT1Method
function HbitT1(p_sps::Array{Array{Int64,1}},n_sps::Array{Array{Int64,1}},
            mstates_p::Array{Array{Int64,1},1},mstates_n::Array{Array{Int64,1},1},
            labels::Array{Array{Array{Int64,1},1},1},TBMEs::Array{Array{Float64,1}})

make bit representation of T=1 (proton-proton&neutron-neutron) interactions for each {m_z}

NuclearToolkit.HbitpnFunction
Hbitpn(p_sps::Array{Array{Int64,1}},n_sps::Array{Array{Int64,1}},
       mstates_p::Array{Array{Int64,1},1},mstates_n::Array{Array{Int64,1},1},
       labels::Array{Array{Int64,1}},TBMEs::Array{Float64,1},zeroME=false)

make bit representation of T=0 (proton-neutron) interactions for each {m_z}

NuclearToolkit.all_perm!Method
all_perm!(ln::Int64,num_valence::Int64,occs::Array{Array{Bool,1}})

make all possible permutation of 'bits'

Example: If 2 protons and 1 neutron are in the 0p-shell space, valence orbits(0p1/2,0p3/2) => -1/2, 1/2, -3/2, -1/2, 1/2, 3/2

configurations are represented like:

proton: 000011, 000101, ..., 110000

neutron: 000001, 000010, ..., 100000

NuclearToolkit.def_mstatesMethod

defmstates(psps,n_sps)

to define the single particle states specified by [n,l,j,tz,mz,p(n)idx]. The last elements pidx and nidx to represent original index of j states, [n,l,j,tz].

NuclearToolkit.main_smMethod

mainsm(sntf,targetnuc,numev,targetJ; savewav=false,q=1,isblock=false,isshow=false,numhistory=3,lm=100,ls=20,tol=1.e-8, inwf="",mdimmode=false,calcmoment=false, visualize_occ=false, gfactors=[1.0,0.0,5.586,-3.826],effcharge=[1.5,0.5])

Digonalize the model-space Hamiltonian

Arguments:

  • sntf: path to input interaction file (.snt fmt)
  • target_nuc: target nucleus
  • num_ev: number of eigenstates to be evaluated
  • target_J: target total J (specified by e.g. [0]). Set to [] if you want lowest states with any J. Note that J should be doubled (J=0=>[0], J=1/2=>[1], J=1=>[2],...)

Optional arguments:

  • q=1 block size for Block-Lanczos methods
  • is_block=false whether or not to use Block algorithm
  • save_wav=false whether or not to save wavefunction file
  • is_show = true to show elapsed time & allocations
  • lm = 100 number of Lanczos vectors to store
  • ls = 20 number of vectors to be used for Thick-Restart
  • tol= 1.e-8 tolerance for convergence check in the Lanczos method
  • in_wf="" path to initial w.f. (for preprocessing)
  • mdimmode=false true => calculate only the M-scheme dimension
  • calc_moment=false true => calculate mu&Q moments
  • visualize_occ=false true => visualize all configurations to be considered
  • gfactors=[1.0,0.0,5.586,-3.826] angular momentum and spin g-factors
  • effcgarge=[1.5,0.5] effective charges
NuclearToolkit.occMethod
occ(p_sps::Array{Array{Int64,1}},mstates_p::Array{Array{Int64,1}},mzp::Array{Int64,1},num_vp::Int64,
    n_sps::Array{Array{Int64,1}},mstates_n::Array{Array{Int64,1}},mzn::Array{Int64,1},num_vn::Int64,Mtot::Int64)

prepare bit representations of proton/neutron Slater determinants => pbits/nbits

joccp, joccn: corresponding occupation numbers for a "j" state, which is used for one-body operation and OBTDs.

Mps/Mns: total M for proton/neutron "blocks"

For 6Li in the p shell and M=0, Mps = [-3,-1,1,3] & Mns = [3,1,-1,-3] blocks => [ (Mp,Mn)=(-3,3),(Mp,Mn)=(-1,1),...]

tdims: array of cumulative number of M-scheme dimensions for "blocks"

tdims =[ # of possible configurations of (-3,3), # of possible configurations of (-1,1),...]

NuclearToolkit.readsmsntMethod
readsmsnt(sntf,Anum)

To read interaction file in ".snt" format.

  • sntf: path to the interaction file
  • Anum: mass number (used for "scaling" of TBMEs)
Note

The current version supports only "properly ordered",like $a \leq b,c \leq d,a \leq c$ for $V(abcd;J)$, snt file.

NuclearToolkit.transit_mainMethod
transit_main(sntf,target_nuc,jl2,jr2,in_wfs;
             num_ev_l=100,num_ev_r=100,q=1,is_block=false,is_show=true,
             calc_EM=true,gfactors=[1.0,0.0,5.586,-3.826],eff_charge=[1.5,0.5])

Calculate the M1&E2 transitions for two wavefunctions

Arguments

  • sntf:path to the interaction file
  • target_nuc:target nucleus in string e.g., "Si28"
  • jl2:J*2 for the left w.f.
  • jr2:J*2 for the right w.f.
  • in_wfs:["path to left wf","path to right wf"]

Optional arguments

  • num_ev_l=100:upper limit of the number of eigenvectors for the left w.f.
  • num_ev_r=100:upper limit of the number of eigenvectors for the right w.f.
  • is_show=true:to display the TimerOutput
  • gfactors=[1.0,0.0,5.586,-3.826]:g factors [glp,gln,gsp,gsn]
  • eff_charge=[1.5,0.5]:effective charges [ep,en]

Optional arguments (not used, but for future developments)

  • q=1:block size
  • is_block=false:to use Block algorithm
NuclearToolkit.solveECMethod
solveEC(Hs,target_nuc,tJNs;
        write_appwav=false,verbose=false,calc_moment=true,wpath="./",is_show=false,
        gfactors = [1.0,0.0,5.586,-3.826],effcharge=[1.5,0.5],exact_logf="")

To solve EC (generalized eigenvalue problem) to approximate the eigenpairs for a given interaction.

\[H \vec{v} = \lambda N \vec{v}\]

Transition densities and overlap matrix for H and N are read from "tdmat/" directory (To be changed to more flexible)

Arguments:

  • Hs:array of paths to interaction files (.snt fmt)
  • target_nuc: target nucleus
  • tJNs:array of target total J (doubled) and number of eigenstates to be evaluated e.g., [ [0,5], [2,5] ], when you want five lowest J=0 and J=1 states.

Optional arguments:

  • write_appwav=false:write out the approximate wavefunction
  • verbose=false:to print (stdout) approx. energies for each interaction
  • calc_moment=true: to calculate mu&Q moments
  • wpath="./": path to sample eigenvectors to construct approx. w.f.
  • is_show=false: to show TimerOutput
  • gfactors=[1.0,0.0,5.586,-3.826]: g-factors to evaluate moments
  • effcharge=[1.5,0.5]:effective charges to evaluate moments

Optional arguments for author's own purpose

  • exact_logf="":path to logfile for E(EC) vs E(Exact) plot
Note

All the effective interactions must be in "the same order" and must be consistent with interaction file from which the transition density matrices were made.

NuclearToolkit.read_kshell_summaryMethod
read_kshell_summary(fns::Vector{String};targetJpi="",nuc="")

Reading KSHELL summary file from specified paths, fns. In some beta-decay studies, multiple candidates for parent g.s. will be considered. In that case, please specify optional argument targetJpi e.g., targetJpi="Si36", otherwise the state havbing the lowest energy in summary file(s) will be regarded as the parent ground state. Example:

Energy levels

N    J prty N_Jp    T     E(MeV)  Ex(MeV)  log-file

1     0 +     1     6   -319.906    0.000  log_Ar48_SDPFSDG_j0p.txt 
Note

J is doubled in old versions of KSHELL (kshell_ui.py).

NuclearToolkit.Ecoulomb_SMMethod
Ecoulomb_SM(Z,N)

Additional Coulomb energy (MeV)

\[E_C(Z,N) = 0.7 \frac{Z(Z-1) -0.76 [Z(Z-1)]^{2/3}}{R_c} \\ R_c = e^{1.5/A} A^{1/3} \left( 0.946 -0.573 \left( \frac{2T}{A} \right)^2 \right) \\ T = |Z-N|/2 \]

used in references.

  • S.Yoshida et al., Phys. Rev. C 97, 054321 (2018).
  • J. Duflo and A. P. Zuker, Phys. Rev. C 52, R23 (1995).
  • E. Caurier et al., Phys. Rev. C 59, 2033 (1999).
NuclearToolkit.Fermi_functionMethod
Fermi_function(Z,We,R)

\[F(Z,W) = 4 (2p_e R)^{-2(1-\gamma_1)} \exp \left( \pi y \frac{|\Gamma(\gamma_1 + iy)|^2}{ |\Gamma(2\gamma_1 + 1)|^2} \right) \\ \gamma_k = \sqrt{k^2 - \alpha^2 Z^2} \\ y = \alpha Z W / p_e\]

NuclearToolkit.Fermi_integralFunction
Fermi_integral(Qval,Ex,pZ,A,nmesh=40)

Fermi integrals are evaluated with Eqs. in Chapter 5 of "Weak Interactions and Nuclear Beta Decay" by H. F. Schopper.]

\[f_0 = \int^{W_0}_{1} F(Z,W) \sqrt{W^2-1} W(W_0-W)^2 dW, \\ W_0 = Q(\beta^-) + 1 - E_\mathrm{ex.}\]

Note that Ws are in natural unit, divided by $m_ec^2$.

NuclearToolkit.ME_FF_func!Method
ME_FF_func!(chan,qfactors,C_J,Mred,lambda_Ce,mass_n_nu,dict)

Function to calc FF matrix elements. Since M0rs,M1r,M1rs,M2rs (M0sp,M1p) are evaluated in fm in KSHELL, it is needed to multiply 1/lambdaCe (lambdaCe) to get transition matrix element in natural unit.

NuclearToolkit.eval_betadecay_from_kshell_logMethod
eval_betadecay_from_kshell_log(fns_sum_parent::Vector{String},fns_sum_daughter::Vector{String},fns_GT::Vector{String},fns_FF::Vector{String},parentJpi::String;pnuc="",dnuc="",verbose=false)

Arguments

  • fns_sum_parent::Vector{String} vector of paths to KSHELL-summary files for the parent nucleus
  • fns_sum_daughter::Vector{String} vector of paths to KSHELL-summary files for the daughter nucleus
  • fns_GT::Vector{String} vector of paths to KSHELL-log files for Gamow-Teller transitions.
  • fns_FF::Vector{String} vector of paths to KSHELL-log files for First-Forbidden transitions.

Optional Arguments

  • pnuc::String to specify the parent nuclei explicitly
  • dnuc::String to specify the daughter nuclei explicitly
  • verbose::Bool option to show GT/FF transitions for each state