.. _DanceQ: ################## Welcome to DanceQ! ################## .. admonition:: DanceQ **D**\ ivide- **A**\ nd-conquer **N**\ umber **C**\ onserving **E**\ xact diagonalization for **Q**\ uantum systems :ref:`DanceQ` is published on `SciPost `_:footcite:p:`DanceQ`. It 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 `_. .. admonition:: Highlights * Arbitrary qudit and spin :math:`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 :ref:`interface` is user-friendly and straightforward to use: .. code-block:: cpp /* 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>> */ 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 :ref:`here`. The core -------- The base of the application is the :ref:`State class` that provides compressed storage of product states, similar to :code:`std::bitset`. The :ref:`Basis class` implements efficient enumeration, iteration, and manipulations of all :ref:`states` in a given particle number sector. A :ref:`State` :math:`\vert \Psi\rangle` consists of :math:`L` individual qudits :math:`\vert \sigma_i\rangle\in\{\vert0\rangle , .., \vert Q-1\rangle\}` with a local Hilbert space dimension of :math:`Q`: .. math:: :name: eq:particle_number \vert\Psi\rangle = \vert\sigma_0\rangle \otimes \vert\sigma_1\rangle \dots \otimes \vert\sigma_{L-1}\rangle = \vert\sigma_0;\dots;\sigma_{L-1}\rangle \quad\mathrm{ with }\quad n = \sum_{i=0}^{L-1} \sigma_i\,. The total number of :ref:`States` with the same particle number :math:`n` is (cf. proof :ref:`in Appendix A`) .. math:: :name: eq:dim D_Q(L,n) = \sum_{k=0}^{\lfloor n / Q \rfloor} (-1)^k\binom{L}{k}\binom{L-1+n-qk}{L-1}\,. The :ref:`Basis class` of :ref:`DanceQ` provides an elegant tool to uniquely map between :ref:`states` :math:`\vert \Psi_i\rangle` and their zero based indices :math:`i=0,\dots,D_Q(L,n)-1`. Let's have a look at a small example for :math:`L=4` sites with :math:`Q=3`, and :math:`n=3`: +------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+ | :math:`\vert\Psi_0 \rangle = \vert 0;0;1;2\rangle` | :math:`\vert\Psi_1 \rangle = \vert 0;0;2;1\rangle` | :math:`\vert\Psi_2 \rangle = \vert 0;1;0;2\rangle` | :math:`\vert\Psi_3 \rangle = \vert 0;1;1;1\rangle` | +------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+ | :math:`\vert\Psi_4 \rangle = \vert 0;1;2;0\rangle` | :math:`\vert\Psi_5 \rangle = \vert 0;2;0;1\rangle` | :math:`\vert\Psi_6 \rangle = \vert 0;2;1;0\rangle` | :math:`\vert\Psi_7 \rangle = \vert 1;0;0;2\rangle` | +------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+ | :math:`\vert\Psi_8 \rangle = \vert 1;0;1;1\rangle` | :math:`\vert\Psi_9 \rangle = \vert 1;0;2;0\rangle` | :math:`\vert\Psi_{10}\rangle = \vert 1;1;0;1\rangle` | :math:`\vert\Psi_{11}\rangle = \vert 1;1;1;0\rangle` | +------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+------------------------------------------------------+ | :math:`\vert\Psi_{12}\rangle = \vert 1;2;0;0\rangle` | :math:`\vert\Psi_{13}\rangle = \vert 2;0;0;1\rangle` | :math:`\vert\Psi_{14}\rangle = \vert 2;0;1;0\rangle` | :math:`\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 :ref:`states` in a list. Our approach is an extension of :footcite:t:`Lin`'s original idea of splitting the system into two subsystems in the case :math:`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 :ref:`Basis` or :ref:`State class` you can define :ref:`Operators` from arbitrary tensor products. An :ref:`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 :ref:`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, :ref:`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 :ref:`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 :ref:`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. .. admonition:: Available matrix formats * :ref:`create_DenseMatrix()`: Dense matrix using two-dimensional ``std::vector`` * :ref:`create_EigenDenseMatrix()`: Dense matrix using `Eigen `_ * :ref:`create_SparseMatrix()`: Sparse matrix using the :ref:`SparseMatrix class` * :ref:`create_ShellMatrix()`: **Matrix-free** multiplication using the :ref:`ShellMatrix class` * :ref:`create_PetscSparseMatrix()`: Sparse matrix using `Petsc `_ * :ref:`create_PetscShellMatrix()`: **Matrix-free** multiplication using `Petsc `_ with the :ref:`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 :ref:`Installation and Configuration`. After briefly outlining the basic usage in :ref:`Usage`, we discuss various examples for the different physical problems using the different extensions can be found :ref:`here`. A detailed explanation of the algorithm, an associated academic paper, and a **bibtex** entry for citations are listed in :ref:`Paper`. .. toctree:: :maxdepth: 2 Installation.rst Usage.rst Examples.rst Paper.rst There is a :ref:`hierarchy of classes` where the :ref:`State class` forms the algorithm's core. Next, the :ref:`Basis class` is built on top of the :ref:`State class` which defines the :ref:`Operator class`. All classes are defined in the global namespace **danceq**. Different classes that are used by :ref:`Operator class` are collected in :ref:`Miscellaneous`. .. toctree:: :maxdepth: 1 :titlesonly: State.rst Basis.rst Operator.rst Miscellaneous.rst .. * :ref:`genindex` .. * :ref:`modindex` .. * :ref:`search` .. footbibliography::