ShellModel
Documentation for ShellModel.
ShellModel.HbitT1
ShellModel.Hbitpn
ShellModel.all_perm!
ShellModel.def_mstates
ShellModel.embedding_crsp
ShellModel.main
ShellModel.occ
ShellModel.occ
ShellModel.readsnt
ShellModel.samplerun
ShellModel.solveEC
ShellModel.tf_truncation_occ!
ShellModel.transit_main
ShellModel.HbitT1
— Methodfunction 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}
ShellModel.Hbitpn
— FunctionHbitpn(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}
ShellModel.all_perm!
— Methodall_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
ShellModel.def_mstates
— Methoddefmstates(psps,n_sps)
to define the single particle states specified by m_z
ShellModel.embedding_crsp
— MethodTo store correspondance between index for vectors w/ and w/o truncation
ShellModel.main
— Methodmain(sntf,targetnuc,targetJNs; savewav=false,q=1,isblock=false,isshow=false,numhistory=3,lm=100,ls=20,tol=1.e-6, inwf="",mdimmode=false,calcmoment = 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 nucleusnum_ev
: number of eigenstates to be evaluatedtarget_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 methodsis_block=false
whether or not to use Block algorithmsave_wav=false
whether or not to save wavefunction fileis_show = true
to show elapsed time & allocationslm = 100
number of Lanczos vectors to storels = 15
number of vectors to be used for Thick-Restarttol= 1.e-6
tolerance for convergence check in the Lanczos methodin_wf=""
path to initial w.f. (for preprocessing)mdimmode=false
true
=> calculate only the M-scheme dimensioncalc_moment=false
true
=> calculate mu&Q momentsgfactors=[1.0,0.0,5.586,-3.826]
angular momentum and spin g-factorseffcgarge=[1.5,0.5]
effective charges
ShellModel.occ
— Methodocc function to embedding truncated vector to full space one
ShellModel.occ
— Methodocc(modelspace,Mtot,truncation_params)
prepare bit representations of proton/neutron Slater determinants => pbits/nbits
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),...]
ShellModel.readsnt
— Methodreadsnt(sntf,Anum)
To read interaction file in ".snt" format.
sntf
: path to the interaction fileAnum
: mass number (used for "scaling" of TBMEs)
The current version supports only "properly ordered" interaction file in .snt format.
A .snt file can be ordered to be a<=b,c<=d,a<=c for V(abcd;J) by the Python script "ShellModel.jl/src/makeorderedsnt.py"(, which will be replaced by Julia implementation...).
ShellModel.samplerun
— Methodsamplerun()
Sample scripts to calculate
- (a) 10 lowest states of 28Si in sd shell
- (b) 10 lowest states of 28Si with J=0
- (c) EC estimates of 10 lowest J=0 states of 28Si
- (d) (b) with the preprocessing
To run samplerun()
, you need a copy of ShellModel.jl in your environment. Try e.g.,
git clone https://github.com/SotaYoshida/ShellModel.jl
and then, execute samplerun()
in the repository:
>julia using ShellModel
>julia samplerun()
or
julia -t 12 sample_run.jl
ShellModel.solveEC
— MethodsolveEC(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 nucleustJNs
: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 wavefunctionverbose=false
:to print (stdout) approx. energies for each interactioncalc_moment=true
: to calculate mu&Q momentswpath="./"
: path to sample eigenvectors to construct approx. w.f.is_show=false
: to show TimerOutputgfactors=[1.0,0.0,5.586,-3.826]
: g-factors to evaluate momentseffcharge=[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
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.
ShellModel.tf_truncation_occ!
— Methodtruncation scheme by specifying occupation numbers for a given [n,l,j,tz]
ShellModel.transit_main
— Methodtransit_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 filetarget_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 TimerOutputgfactors=[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 sizeis_block=false
:to use Block algorithm