IMSRG

Files in src/IMSRG.jl for IM-SRG calculations:

  • imsrg_util.jl: contains main and util functions
  • commutator.jl: functions to calculate commutators and BCH transform to carry out IMSRG flow with Magnus expansion
  • valencespace.jl: functions for Valence-space IM-SRG (VS-IMSRG) calculations to derive shell-model effective interactions/operators

Since we use the so-called Magnus formulation of IMSRG flow, $\Omega$

\[H(s) = e^{\Omega(s)}H(0)e^{-\Omega(s)} \\\\ e^{\Omega(s+ds)} = e^{\eta(s)ds}e^{\Omega(s)}\]

is needed to explicitly calculate the flow generator $\eta$. See e.g. Annual Review of Nuclear and Particle Science 69, 307-362 for more details.

In the package, temporary binary files for $\Omega$ are genereted in flowOmega directory (which will be made if not available). Although users do not need to read or write these binary files in normal use, the contents of the binary files are described for specific uses and development, e.g. re-start IMSRG flow from the temporary files.

A temporary file is generated like flowOmega/Omega_1234He4_1.bin. This contains

  • number of single particle basis (Int64), which is assumbed to be common between proton and neutron.
  • One-body matrix element for protons and neutrons (Matrix{Float64})
  • number of two-body channels (Int64)
  • dimensions of each channel (Vector{Int64})
  • Two-body matrix element for protons and neutrons (Vector{Matrix{Float64}})

Note that the all matrix elements are written in row-major order (I know Julia is column-major language).


NuclearToolkit.BCH_ProductMethod
BCH_Product(X::Op,Y::Op,Z::Op,tmpOp::Op,Nested::Op,ncomm::Vector{Int64},norms::Vector{Float64},Chan1b::chan1b,Chan2bD::chan2bD,HFobj::HamiltonianNormalOrdered,dictMono,d6j_lj,PandyaObj::PandyaObject,to;tol=1.e-4) where Op <:Operator

returns $Z$ to satisfy: $e^Z = e^X e^Y$.

$Z$ is calculated with Baker–Campbell–Hausdorff (BCH) formula: $Z = X + Y + \sum^{\infty}_{k=1} ad^{(k)}_\Omega(\eta)$ $Z = X + Y + 1/2[X, Y] + 1/12 [X,[X,Y]] + 1/12 [Y,[Y,X]] -1/24 [Y,[X,[X,Y]]] -1/720 [Y,[Y,[Y,[Y,X]]]] -1/720 [X,[X,[X,[X,Y]]]] +...$

For IMSRG flow of $H(s)$, $X=\eta(s) ds$, $Y=\Omega(s)$, and $Z=\Omega(s+ds)$

source
NuclearToolkit.BCH_TransformMethod
BCH_Transform(Omega,O0,Os,tOp,Nested,ncomm,norms,Chan1b,Chan2bD,HFobj,dictMono,d6j_lj,PandyaObj,to;tol=1.e-9,maxit=50,verbose=false)

Update $Os$ via $e^ABe^{-A} =B+[A,B]+1/2![A,[A,B]]+1/3![A,[A,[A,B]]]+...$ Note that $Os,tOp,Nested$ are set zero.

  • Omega: $\Omega(s+ds)$
  • O0: $O(s=0)$
  • Os: $O(s+ds)$
source
NuclearToolkit.Bernoulli_numberMethod
Bernoulli_number(k::Int64)::Float64

Return the k-th Bernoulli number. Some special cases are listed below

  • B(0) = 1
  • B(1) = -1/2. In some literature, B(1) = 1/2.
  • B(2) = 1/6
  • B(2n+1) = 0 (for n >=1)

For others, see the references A000367 and A002445.

source
NuclearToolkit.OpCommutator!Method
OpCommutator!(X::Op,Y::Op,ret::Op,HFobj,Chan1b,Chan2bD,dictMono,d6j_lj,PandyaObj,to) where{Op <: Operator}

overwrite ret operator by the commutator $[X,Y]$

source
NuclearToolkit.calcZbar!Method
calcZbar!(Xbar,Ybar,PhaseMat,PhaseMatY,tmpMat,hy,nph_kets,nKets_cc,Zlefthalf,Zrighthalf,hz)
  • Xbar: (nKets_cc, 2*nph_kets) matrix or SubArray
  • Ybar: (2*nph_kets,nKets_cc) matrix or SubArray
  • PhaseMatY: (nph_kets,nKets_cc) matrix or SubArray
  • Zbar: (nKets_cc, 2*nKets_cc) matrix or SubArray
source
NuclearToolkit.comm121ss!Method
comm121ss!(X,Y,ret,HFobj,Chan1b,Chan2b,dictMono,PandyaObj)

returns $[X_1,Y_2] - [Y_1,X_2]$, whose elements are given as

$[X_1,Y_2]_{ij} = \frac{1}{2j_i+1}\sum_{ab} (n_a \bar{n}_b) \sum_{J} (2J+1) (X_{ab} Y^J_{biaj} - X_{ba} Y^J_{aibj})$

source
NuclearToolkit.comm220ss!Method
comm220ss!(X::Op,Y::Op,Z::Op,HFobj::HamiltonianNormalOrdered,Chan2b::Vector{chan2b}) where Op<:Operator

$[X_2,Y_2]_0 = 2 \sum_{J}(2J+1) \mathrm{Tr}(X_{hh'pp'} Y_{pp'hh'})$

source
NuclearToolkit.comm221ss!Method
comm221ss!(X::Op,Y::Op,ret::Op,HFobj::HamiltonianNormalOrdered,Chan1b::chan1b,Chan2bD::chan2bD,PandyaObj::PandyaObject) where Op<:Operator

returns $[X_2,Y_2]_1 - [Y_2,X_2]_1$, whose elements are given as

$[X_2,Y_2]_{ij} = 1/(2[j_i]) \sum_{abc}\sum_{J}[J](n'_an'_bn_c-n_an_bn'_c)(X_{2,ciab}Y_{2,abcj}-Y_{2,ciab}X_{2,abcj})$

source
NuclearToolkit.single_121Method
single_121(a,b,i,j,o1b,o2bs,sps,key,targetDict;verbose=false)

Function to calculate 121part $[X_1,Y_2]-[Y_1,X_2]$ for given i,j and a,b.

$\sum_{J} [J]^2 (o_{1,ab}o_{2,biaj} - o_{1,ba}o_{2,aibj})$

source
NuclearToolkit.GatherOmegaMethod
GatherOmega(Omega,nOmega,gatherer,tmpOp,Nested,H0,Hs,ncomm,norms,Chan1b,Chan2bD,HFobj,IMSRGobj,dictMono,d6j_lj,PandyaObj,maxnormOmega,to)

This may not be used now.

source
NuclearToolkit.Get2bDenominatorMethod
Get2bDenominator(ch,pnrank,a,b,i,j,na,nb,ni,nj,f,Delta,dictMono,key;verbose=false)

$f_{aa}+f_{bb}-f_{ii}-f_{jj}+G_{abij} +\Delta$

with $G_{abij} = \Gamma_{abab} + \Gamma_{ijij} - (\Gamma_{aiai} + \Gamma_{bjbj} + [a \leftrightarrow b])$

source
NuclearToolkit.calc_Eta_atan!Method
calc_Eta_atan!(HFobj::HamiltonianNormalOrdered,IMSRGobj::IMSRGobject,Chan2b::Vector{chan2b},dictMono,norms)

calc. $\eta(s)$ with atan generator

source
NuclearToolkit.flow_OperatorsMethod
flow_Operators(binfo,HFobj,IMSRGobj,PandyaObj,Chan1b,Chan2bD,d6j_lj,dictMono,Operators,to)

consistent IMSRG flow of scaler operators (Rp2) using written Omega

source
NuclearToolkit.imsrg_mainMethod
imsrg_main(binfo::basedat,Chan1b,Chan2bD,HFobj,dictMono,dWS,valencespace,Operators,to; core_generator_type="atan",valence_generator_type="shell-model-atan",denominatorDelta=0.0)

Arguments

  • binfo::basedat struct basedat(nuc::nuclei,sntf::String,hw::Float,emax::Int)
  • Chan1b::chan1b struct for one-body stuffs
  • Chan2bD::chan2bD struct for two-body stuffs (e.g., dict to get idx from JPT)
  • HFobj::HamiltonianNormalOrdered struct HNO, which includes info. of HF solution (HF energy, occupation, f,Gamma,...)
  • dictMono::Dict dictionary to get Vmonopole
  • dWS dictionary of preallocated wigner-symbols
  • valencespace to specify valence space
  • Operators::Vector{String} non-Hamiltonian operators
  • to TimerOutput object to measure runtime&memory allocations

Optional Arguments

  • delete_Ops if true, delete Operators with current pid after IMSRGflow
  • core_generator_type only the "atan" is available
  • valence_generator_type only the "shell-model-atan" is available
  • denominatorDelta::Float denominator Delta, which is needed for multi-major shell decoupling
  • debugmode=0: 0: no debug, 1: debug, 2: debug with more info
  • restart_from_files: files to be read for restart 1st one is for IMSRG and 2nd one is for VSIMSRG
source
NuclearToolkit.init_IMSRGobjectMethod
init_IMSRGobject(HFobj;smax=500.0,dsmax=0.5,maxnormOmega=0.25,tol=1.e-6,eta_criterion=1.e-6,denominatorDelta=0.0)

Constructor for IMSRGobject

  • H0::Operator for starting point of BCH product
  • H::Operator Hamiltonian $H(s)$
  • s::Vector{Float} current $s$ and $ds$
  • smax::Float maximum $s$
  • dsmax::Float maximum $ds$
  • maxnormOmega::Float maximum ||Omega||
  • eta::Operator generator of IMSRG flow (antihermite Operator)
  • Omega::Operator generator of IMSRG flow (antihermite Operator)
  • eta_criterion::Float ||eta|| to check convergence
  • denominatorDelta::Float64 parameter for multi-major shell decoupling
  • n_written_omega::Int # of written Omega by splitting to solve IMSRGflow
  • Ncomm::Vector{Int} # of commutator evaluated during IMSRG flow
source
NuclearToolkit.make_PandyaKetsMethod

make_PandyaKets(emax,HFobj)

To prepare "kets" for Pandya transformation. For ordinary two-body channels, kets like |i,j=i;J=odd> with `={n,l,j,tz} are hindered, but necessary for Pandya transformation.

source
NuclearToolkit.prep_PandyaLookupMethod
prep_PandyaLookup(binfo::basedat,HFobj::HamiltonianNormalOrdered,Chan1b::chan1b,Chan2bD::chan2bD;rank_J=0,rank_T=0,parity=0,ofst=1000)

constructor of utils for Pandya transformation and others numbersPandya:[ch,nKetcc,nhh,nph] for ich (channel index of Chan2b_Pandya)

source
NuclearToolkit.print_flowstatusFunction

print_flowstatus(istep,s,ncomm,norms,IMSRGobj)

Function to print flowstatus s,E0,1b&2b norm for Omega, 1b&2b norm for Eta, Ncomm, nwritten

source
NuclearToolkit.check_array_valencespaceMethod
check_array_valencespace(valencespace::Vector{Vector{Int64}},HFobj,v)

check valence space and overwrtie SingleParticleState.v/q

specified by array of sps (e.g., [[0,1,1,-1],[0,1,3,-1], [0,1,1,1],[0,1,3,1]])

source
NuclearToolkit.check_major_valencespaceMethod
check_major_valencespace(str::String,HFobj,v)

Function to check valence space and overwrite v and q fields of SingleParticleState The valencespace is specified by argument str (e.g. "p-shell")

source
NuclearToolkit.write_vs_sntMethod
write_vs_snt(binfo,HFobj,IMSRGobj,Operators,effOps,Chan1b,Chan2bD,vspace)

Function to write out valence space effective interaction in snt (KSHELL/ShellModel.jl) format.

source