Snake is a research code implementing ground-state DMRG and adaptive real-time tDMRG for 1D lattice models and Wilson-chain representations of quantum impurity problems. Model inputs are generated by MATLAB scripts and written as binary files into a model/ directory; the C++ executable reads ./model/*.dat and writes results to results/.
This repository contains an implementation of the methods and models used in:
- Cheng Guo, Using Density Matrix Renormalization Group Method to Study the Real-Time Dynamics of 1D Quantum Systems and Quantum Impurity Systems (Master thesis in Chinese, 2008)
- C. Guo, A. Weichselbaum, S. Kehrein, T. Xiang, J. von Delft, Density matrix renormalization group study of a quantum impurity model with Landau–Zener time-dependent Hamiltonian, Phys. Rev. B 79, 115137 (2009)
- Ground states: “infinite-system” (warmup/growth) DMRG + finite-system sweeps in a fixed good-quantum-number sector.
- Real-time dynamics: adaptive tDMRG using a second-order Suzuki–Trotter decomposition implemented as a right sweep and left sweep of nearest-neighbor two-site gates.
- Impurity support: the first bond (impurity ↔ first bath site) may have a time-dependent gate
U_imp(t)provided step-by-step. - Workflow: generate
model/with MATLAB → runSnakein that folder → analyzeresults/with MATLAB.
The default executable (cppsrc/snake.cpp) runs a fixed pipeline: iDMRG → fDMRG → tDMRG. There are no command-line flags; to change the pipeline (e.g. “ground state only”), edit cppsrc/snake.cpp and rebuild.
Dependencies:
- BLAS + LAPACK
- ARPACK
- LAPACK++ (
lapackpp) headers and library - C++ compiler (CMake also enables Fortran for typical BLAS/LAPACK/ARPACK link stacks)
Build with CMake:
mkdir -p build
cd build
cmake ..
cmake --build . -jThis produces build/cppsrc/Snake.
LAPACK++ path note: cppsrc/CMakeLists.txt assumes
- headers:
$HOME/usr/include/lapackpp - libs:
$HOME/usr/lib
Adjust cppsrc/CMakeLists.txt (or provide symlinks) if your installation differs.
From the repo root:
% gen_model(Omega, omegaratio, Lambda, BathL, z)
gen_model(0.05, 1.0, 2.0, 50, 1);gen_model.m creates a folder named like ratio...Lambda...BathL...z.../ containing:
model/(binary input files)Snake(copied from./build/cppsrc/Snake)
Notes:
gen_model.musessystem('cp ...'); it assumes a Unix-like shell (Linux/macOS, or Windows via WSL). Edit that line if needed.- The shipped generator is an example; for the dissipative Landau–Zener / resonant-level model you typically create your own generator (see “Physics background” and “Input contract”).
cd ratio1Lambda2BathL50z1 # name depends on generator parameters
./SnakeSnake reads inputs from ./model/ relative to your current working directory, so you can also run ../../build/cppsrc/Snake as long as you cd into the run folder first.
From the directory containing your run folders:
plotsigmaz.m(time trace fromresults/sigmaz_t.dat)plotspectral.m(extracts a late-time average fromsigmaz_t.datacross parameter sweeps)plotvisibility.m(coherence magnitude fromresults/rdm.dat)
Standard DMRG finds an approximation to the ground state of a 1D nearest-neighbor Hamiltonian by optimizing an MPS and truncating using reduced density matrices.
In Snake:
- “iDMRG” (
cppsrc/InfiniteDMRG.*) is the classic infinite-system warmup: grow to the target finite length by repeatedly adding two sites and truncating. - “fDMRG” (
cppsrc/FiniteDMRG.*) then performs finite-system sweeps to variationally improve the state (default: 1 sweep, hard-coded).
The superblock ground state at each step is computed with ARPACK (see cppsrc/SuperChain_iDMRG.cpp).
For a nearest-neighbor Hamiltonian
In Snake:
- the MATLAB generator precomputes the half-step gates
$e^{-i(\tau/2)h_{i,i+1}}$ and writes them tomodel/rt_T0.dat; - a time-dependent half-step gate for the first bond (impurity + first bath site) is provided step-by-step in
model/rt_H1_T0.dat; -
cppsrc/SuperChain_tDMRG.cppapplies these gates in a left→right sweep and a right→left sweep, updating the basis using truncation matrices saved during fDMRG.
The code writes a local observable labelled sigma_z(t) each step; for common impurity mappings where the impurity is a two-level system with
Quantum impurity models often couple a local degree of freedom to a continuum bath with bandwidth
- the logarithmic discretization parameter
$\Lambda>1$ (smaller$\Lambda$ is more accurate but requires longer chains), - a
$z$ -shift (used for$z$ -averaging to reduce discretization artifacts), - the finite chain length
$L$ .
This repo contains MATLAB helpers for chain construction:
roklogdiscr.m: logarithmic discretization → star Hamiltonian → chain parameters (Appendix A of Phys. Rev. B 71, 045122 is referenced in the comments)lindiscr.m: linear discretization + Lanczos tridiagonalizationstar2tridiag.m: star → tridiagonal chain transformation
The Phys. Rev. B 79, 115137 (2009) study focuses on the (spinless) resonant-level model with a linearly swept impurity level,
with hybridization width
Under a standard mapping, this model is equivalent to an Ohmic spin-boson model at dissipation strength
separates fast (
Snake’s tDMRG engine supports this class of problems as long as you generate the corresponding model/ inputs (chain, two-site Hamiltonians, and time-dependent impurity gate).
Finite-temperature DMRG via purification and imaginary-time evolution is described in the thesis background, and this repo contains experimental/unfinished pieces (cppsrc/SuperChain_FTDMRG.cpp, cppsrc/FiniteTemperatureDMRG.*). The default Snake pipeline does not currently expose a complete FTDMRG workflow.
Snake’s C++ core assumes a nearest-neighbor chain built from local operators
From model/Hfac.dat it reads per-bond/per-site scalars:
-
hopT(i)(nearest-neighbor hopping strength on bond$i$ ) -
onesiteE(i)(on-site potential on site$i$ ) -
twositesV(i)(nearest-neighbor interaction on bond$i$ )
For real-valued chains, the built-in two-site term is (see cppsrc/Chain.cpp):
- hopping:
$t,(c_i^\dagger c_{i+1}+c_{i+1}^\dagger c_i)$ fromhopT - interaction (if nonzero):
$V,(n_i-\tfrac12)(n_{i+1}-\tfrac12)$ fromtwositesV - on-site term on the newly added site during growth/sweeps from
onesiteE
Anything more model-specific (impurity terms, explicitly time-dependent pieces, custom two-site Hamiltonians) is typically injected through the MATLAB-generated matrices:
-
model/HC.dat: two-site Hamiltonians$H_{i,i+1}$ used as the “middle bond” term during DMRG -
model/rt_T0.dat: real-time evolution gates for bonds 2…L−1 (bond 1 is handled separately; see below) -
model/rt_H1_T0.dat: real-time evolution gates for bond 1 (impurity + first bath site), one gate per time step
gen_model.m builds a chain with:
- site 1: a two-level “spin” site (MATLAB stores
sigmaxas the first operator andsigmazas the second) - sites 2…L: spinless fermion sites (first operator
c, second operatorn)
It prepares:
- a biased initial Hamiltonian (large
ed0field on the impurity) for the ground-state DMRG stage, and - a driven impurity Hamiltonian during real-time evolution (a term
$\Omega\cos(\omega t)\sigma_x$ plus static$\Delta\sigma_z/2$ and density couplings to the first fermion site).
Use it as a template for your own generator when implementing the resonant-level / Landau–Zener sweep (or other impurity models).
C++ (cppsrc/):
snake.cpp: hard-coded pipeline iDMRG → fDMRG → tDMRGInfiniteDMRG.*: infinite-system warmup (growing the chain)FiniteDMRG.*: finite-system sweeps and handoff to tDMRGAdaptiveTimeDependentDMRG.*: readsrt_*.datgates and runs time evolutionSuperChain.*: superblock object (wavefunction, bases, truncation matrices)SuperChain_iDMRG.cpp: ARPACK ground state solver and middle-bond mappingsSuperChain_tDMRG.cpp: Suzuki–Trotter time evolution + outputSuperChain_FTDMRG.cpp: experimental imaginary-time/temperature evolution
Chain.*,ChainHamiltonian.*: block growth and Hamiltonian assemblydtmat.*: reduced density matrix, truncation matrix, von Neumann entropysite.*: local operators and binary I/O from MATLAB-generated filesgqn*.{h,cpp},gqnmat.*: “good quantum number” bookkeeping and block-sparse matricespublic.*: binary readers/writers and small linear-algebra helpers
MATLAB (repo root):
gen_model.m,gen_model_multipara.m: example model generation (creates run folders)savepara.m,mat2file.m,sitebase2file.m: binary writers formodel/lindiscr.m,roklogdiscr.m,star2tridiag.m: bath discretization / chain mapping helpersplotsigmaz.m,plotspectral.m,plotvisibility.m: plotting scripts forresults/
Snake expects a model/ directory in the current working directory with fixed filenames. savepara.m writes them; the C++ code reads them via cppsrc/public.h and friends.
model/problemparmeters.dat(note the spelling; Snake expects it exactly)int32: chain lengthLint32: target good quantum numberTGQN(a particle-number-like constraint)
model/site_base.dat: local basis + good-quantum-number structure per site (sitebase2file.m)model/site_operators.dat: per-site local operators (savepara.m)- Snake reads two operators per site: “a” (first) and “n” (second)
model/Hfac.dat: real-valued scalar parameters (savepara.m)hopT(lengthL-1),onesiteE(lengthL),twositesV(lengthL-1)
model/HC.dat: two-site Hamiltonians for each bond (para.h0{i}ingen_model.m)model/rt_T0.dat: real-time evolution gates for all bondsint32: number of time stepst_num- then
L-1complex matrices: the half-step two-site gates for each bond - note: during evolution Snake uses
rt_H1_T0.datfor bond 1;rt_T0.datentry #1 is currently ignored (often set to identity), but must still be present for the reader
model/rt_H1_T0.dat: real-time evolution gates for bond 1t_numcomplex matrices: the half-step time-dependent impurity gate per step
Binary format (see mat2file.m and the C++ readers):
- integers are
int32 - real scalars are IEEE-754 doubles (
real*8) - matrices are written as:
int32 nrow,int32 ncol, then column-major data - complex matrices are written as interleaved
(real, imag)double pairs per element (column-major)
Snake creates and overwrites directories on each run:
results/: kept after exitdata/: intermediate blocks; deleted on exit and also cleared at startup
Common output files from tDMRG (cppsrc/SuperChain_tDMRG.cpp):
results/sigmaz_t.dat: one value per time step (no explicit time column)results/vonneumannentropy_t.dat: entanglement entropies during evolution (multiple cuts per step)results/rdm.dat: reduced-density-matrix coherence element(s), written as columnsRe Imresults/steperror_t.dat: created but not currently populated
Because time is not written explicitly, you must know the time step rt_tau, t_end in gen_model.m).
The main numerical controls are:
-
Kept states / bond dimension
$m$ ($\chi$ ): defaultm_KeptStatesNum = 50incppsrc/dmrg.cpp. -
Truncation cutoff:
Max_Truncate_Errorincppsrc/setting.h(minimum density-matrix eigenvalue kept during truncation). -
fDMRG sweep count:
m_sweepTimesincppsrc/FiniteDMRG.cpp. -
Real-time Trotter step
$\tau$ : set in your MATLAB generator (rt_tau); smaller$\tau$ reduces Trotter error but increases cost. -
Bath discretization and size:
$\Lambda$ ,$z$ , and chain length$L$ (set by your generator). For Wilson chains,$z$ -averaging (multiple runs with differentz) reduces discretization artifacts.
Practical notes:
- Real-time entanglement typically grows with time; long-time dynamics may require increasing
$m$ until truncation effects become acceptable. - Benchmark small systems against exact diagonalization or known limits before trusting large parameter sweeps.
-
Snakecallsrm -rf resultsandrm -rf dataviasystem(...); run it only inside a dedicated run folder.
Suggested general DMRG/MPS reference:
- U. Schollwöck, “The density-matrix renormalization group in the age of matrix product states” (Annals of Physics 2011)
- C. Guo "用密度矩阵重正化群方法研究一维量子系统及量子杂质系统的动力学性质" (My master these in Chinese)
If you use this code or build on these ideas, please cite:
@article{PhysRevB.79.115137,
title = {Density matrix renormalization group study of a quantum impurity model with Landau-Zener time-dependent Hamiltonian},
author = {Guo, Cheng and Weichselbaum, Andreas and Kehrein, Stefan and Xiang, Tao and von Delft, Jan},
journal = {Phys. Rev. B},
volume = {79},
issue = {11},
pages = {115137},
numpages = {6},
year = {2009},
month = {Mar},
publisher = {American Physical Society},
doi = {10.1103/PhysRevB.79.115137},
url = {https://link.aps.org/doi/10.1103/PhysRevB.79.115137}
}This program has been developed under the supervision of Prof. Tao Xiang and Prof. Jan von Delft. Thanks to Shijie Hu, Honggang Luo, Shaojing Qin, Jize Zhao, Hantao Lu, Zhihui Wang and Shuming Li for either direct contribution or discussions during the development of the program.
