Welcome to DanceQ!
DanceQ
Divide- And-conquer Number Conserving Exact diagonalization for Quantum systems
DanceQ provides an elegant and powerful approach to exploit one of the most common symmetries in quantum physics: The conservation of particle number or total magnetization. It is a flexible and modern – header-only – C++20 implementation to handle arbitrary tensor products in a given particle number sector. We provide multiple extensions to advanced linear algebra libraries like Eigen or Petsc/Slepc including matrix-free multiplications for scalable MPI computations. It is the bridge between an user-friendly interface and state-of-the-art libraries for exact diagonalization. Please find the source code on gitlab: https://gitlab.com/DanceQ/danceq.
Highlights
Arbitrary qudit and spin \(S\) Hamiltonians with and without number conservation
User-friendly interface similar to ITensor
Compatibility with advanced linear algebra libraries like Eigen and Petsc/Slepc
Matrix-free implementation for scalable MPI computations
It covers Hamiltonian and Lindbladian systems
The interface is user-friendly and straightforward to use:
/* Defining the operator (Lindbladian acting on the space of density matrices) with N sites and n particles */
Lindbladian_U1 L(N,n);
/* Simple XXX-Heisenberg model */
for(uint64_t i = 0; i < N; i++){
L.add_operator(.5, {i,(i+1)%N}, {"S+","S-"});
L.add_operator(.5, {i,(i+1)%N}, {"S-","S+"});
L.add_operator(1., {i,(i+1)%N}, {"Sz","Sz"});
}
/* Jump operator acting on the first site */
L.add_jump_operator(.1, {0}, {"Sz"});
/* Dense matrix using std::vector<std::vector<std::complex<double>>> */
auto L_dense = L.create_DenseMatrix();
/* Dense matrix with Eigen using the type Operator::EigenMatrixType */
auto L_eigen = L.create_EigenDenseMatrix();
/* Sparse matrix using the SparseMatrix class */
auto L_sparse = L.create_SparseMatrix();
/* Shell matrix using the ShellMatrix class for matrix-free multiplication */
auto L_shell = L.create_ShellMatrix();
/* Sparse matrix with Petsc MatType MATMPIAIJ */
auto L_sparse_petsc = L.create_PetscSparseMatrix();
/* Shell matrix with Petsc using the ShellMatrix class */
auto L_shell_petsc = L.create_PetscShellMatrix();
Many examples can be founds here.
The core
The base of the application is the State class that provides compressed storage of product states, similar to std::bitset
.
The Basis class implements efficient enumeration, iteration, and manipulations of all states in a given particle number sector.
A State \(\vert \Psi\rangle\) consists of \(L\) individual qudits \(\vert \sigma_i\rangle\in\{\vert0\rangle , .., \vert Q-1\rangle\}\) with a local Hilbert space dimension of \(Q\):
The total number of States with the same particle number \(n\) is (cf. proof in Appendix A)
The Basis class of DanceQ provides an elegant tool to uniquely map between states \(\vert \Psi_i\rangle\) and their zero based indices \(i=0,\dots,D_Q(L,n)-1\). Let’s have a look at a small example for \(L=4\) sites with \(Q=3\), and \(n=3\):
\(\vert\Psi_0 \rangle = \vert 0;0;1;2\rangle\) |
\(\vert\Psi_1 \rangle = \vert 0;0;2;1\rangle\) |
\(\vert\Psi_2 \rangle = \vert 0;1;0;2\rangle\) |
\(\vert\Psi_3 \rangle = \vert 0;1;1;1\rangle\) |
\(\vert\Psi_4 \rangle = \vert 0;1;2;0\rangle\) |
\(\vert\Psi_5 \rangle = \vert 0;2;0;1\rangle\) |
\(\vert\Psi_6 \rangle = \vert 0;2;1;0\rangle\) |
\(\vert\Psi_7 \rangle = \vert 1;0;0;2\rangle\) |
\(\vert\Psi_8 \rangle = \vert 1;0;1;1\rangle\) |
\(\vert\Psi_9 \rangle = \vert 1;0;2;0\rangle\) |
\(\vert\Psi_{10}\rangle = \vert 1;1;0;1\rangle\) |
\(\vert\Psi_{11}\rangle = \vert 1;1;1;0\rangle\) |
\(\vert\Psi_{12}\rangle = \vert 1;2;0;0\rangle\) |
\(\vert\Psi_{13}\rangle = \vert 2;0;0;1\rangle\) |
\(\vert\Psi_{14}\rangle = \vert 2;0;1;0\rangle\) |
\(\vert\Psi_{15}\rangle = \vert 2;1;0;0\rangle\) |
Why do you need it?
While the problem of enumerating the basis in particle number sector seems rather trivial at first glance, the exponentially growing complexity of quantum many-body problems quickly forbids any straightforward implementations such as simply storing states in a list. Our approach is an extension of Lin[1]’s original idea of splitting the system into two subsystems in the case \(Q=2\). We have extended this divide-and-conquer ansatz by splitting the model into arbitrarily many subsystems, allowing an optimal decomposition in a modern C++20 implementation. Our implementation is flexible and type-parametrized, eliminating overheads due to compile-time decision making.
Based on the Basis or State class you can define Operators from arbitrary tensor products. An Operator is a generalization of a matrix and encodes the concept of a matrix-vector product, typically without storing the matrix explicitly. The construction of local Operators like the Hamiltonian or local observables is straightforward. The interface is similar to ITensor and allows a user-friendly implementation of arbitrary Hamiltonians with many-body interactions. Besides native C++ containers, DanceQ comes with different extensions supporting shared and distributed memory.
In particular, Eigen can be used for dense matrix applications and provides access to LAPACK libraries such as Intel’s Math Kernel Library (MKL). For sparse matrix applications, the Operator class supports distributed memory computations with Petsc and Slepc. This is particularly powerful for sparse matrix applications such as determining ground states or real and imaginary time evolution, see here. Additionally, we provide matrix-free approaches based on openMP and MPI that can be used with and without Petsc. In this cases, the matrix is not stored and computed on-the-fly to reduce the memory costs.
Available matrix formats
create_DenseMatrix(): Dense matrix using two-dimensional
std::vector
create_EigenDenseMatrix(): Dense matrix using Eigen
create_SparseMatrix(): Sparse matrix using the SparseMatrix class
create_ShellMatrix(): Matrix-free multiplication using the ShellMatrix class
create_PetscSparseMatrix(): Sparse matrix using Petsc
create_PetscShellMatrix(): Matrix-free multiplication using Petsc with the ShellMatrix class
Content
A guide of how-to-install (remember it is header-only) and configurations for Eigen and Petsc/Slepc can be found in Installation and Configuration. After briefly outlining the basic usage in Usage, we discuss various examples for the different physical problems using the different extensions can be found here. A detailed explanation of the algorithm, an associated academic paper, and a bibtex entry for citations are listed in Paper.
There is a hierarchy of classes where the State class forms the algorithm’s core. Next, the Basis class is built on top of the State class which defines the Operator class. All classes are defined in the global namespace danceq. Different classes that are used by Operator class are collected in Miscellaneous.