Functions

We collect different Functions that are not part of any class. They belong to the name space danceq::.

Lanczos

This function uses the Lanczos algorithm to compute the ground state and excitations using the matrix-free multiplication provided by ShellMatrix class. It is compatible with MPI and openMP.

template<class Op>
std::vector<typename Op::scalartype_real> danceq::lanczos(const Op &H, uint64_t number_of_states = 1, std::vector<danceq::internal::Vector<typename Op::scalartype>> *states = nullptr, uint32_t seed = -1, double overhead_in_GB_per_core = 2., double conv = 1e-9, uint64_t iter_max = 500, bool quiet = false)

A primitive Krylov subspace method to compute the ground state and a few excited states using the ShellMatrix.

This function uses the matrix-free multiplication in the ShellMatrix class to construct the Krylov space and computes the ground state and few excited states. The procedure is known as the Lanczos algorithm. Note that the Operator H has to be hermitian. The generalization to non-hermitian operators is known as the Arnoldi algorithm.

Only the Operator H is required as an input. All other parameters are optional and control the performance and output:

  • H: Pointer to the Operator class. It has to be hermitian.

  • number_of_states (Default: 1): Number of lowest energy states that should be computed.

    Tip

    This is a rather primitive implementation of the algorithm without reorthogonalization or other tricks. If you require multiple excited states (or states within the spectrum) you might want to use Petsc and Slepc as demonstrated here.

  • states (Default: nullptr): Pointer to computed states in form the of Vector class. If it is not defined (set to nullptr), states are not computed and not returned.

  • seed (Default: -1): Seed to initialize a random state. If it is set to -1 a random seed is generated.

  • overhead_in_GB_per_core (Default: 2.): This is only relevant if MPI is enabled. Memory is buffered during the communication between different processes and can be limited. The overhead is set in GB and refers to each process.

  • conv (Default: 1e-9): Convergence criteria to control the precession. The subspace construction is stopped if the difference in lowest energies or beta is below conv. beta refers to the norm of the new basis state within the algorithm.

  • iter_max (Default: 200): Maximal subspace dimension.

The ShellMatrix class uses the Vector class which spans the subspace. Returned states are also members of the Vector class. This function minimizes the number of necessary vectors and, therefore, the memory consumption. The number of vectors used is 1 + number_of_states. An additional vector is required if the eigen states are requested.

Depending on availability different modules are used. First, it utilizes the distributed memory structure if MPI has been found. Second, if possible, matrix-vector multiplication is parallelized with openMP. Third, if neither MPI or openMP is available, the multiplication is serial and uses a single thread.

Note

The functions is only available if Eigen is included to diagonalize the projected Operator in the Krylov subspace.

Parameters:
  • H – Pointer to a hermitian operator

  • number_of_states – The number lowest energy levels that are computed

  • states – Output of states if set (unequal to nullptr)

  • seed – Seed to initialize random vector (-1 will generate a random seed)

  • overhead_in_GB_per_core – Overhead in memory (GB) per rank due to the MPI communication (irrelevant if MPI is not used)

  • conv – Convergence criteria

  • iter_max – Maximal subspace dimension

  • quiet – No output of true

Returns:

Ground-state and excited energies (if requested)

Basic functions

The other functions are generic and used by different classes.

template<class T>
constexpr bool danceq::is_complex(void)

A simple function that evaluates a given primitive datatype at compile time.

Returns true if the T has the form std::complex<…>, otherwise false:

ScalarType c;
if constexpr (danceq::is_complex<ScalarType>()){
    c.real(1.);
    c.imag(0.);
} else {
    c = 1.;
}

Returns:

true it T is complex otherwise false

template<class T>
constexpr MPI_Datatype danceq::get_mpi_type(void) noexcept

A simple function that returns the MPI_Datatype at compile time.

The code is inspired from stackoverflow. The functions is only available if MPI is included.

MPI_Datatype dtype = danceq::get_mpi_type<ScalarType>()

Returns:

MPI_Datatype

template<class T>
std::vector<T> danceq::uniform_dist_mpi(T dim)

Function create a uniform partition of size dim across all MPI ranks.

The size of the returned vector, say v, equals the number of MPI ranks. Each rank is in charge of v[myrank] to v[myrank+1]. Note that the returned vector does not have the entry v[myrank+1] for the last rank which would equal the input dim. The functions is only available if MPI is included.

Parameters:

dim – Total dimension

Returns:

Vector with a uniform partition

Petsc functions

The following functions are only available if Petsc is included.

double danceq::von_Neumann_entropy_from_PetscVec(const Vec &vec, const int32_t root = 0)

Returns the von Neumann entropy from a vectorized density matrix given as a Petsc Vec.

Given a density matrix \(\rho\), the entropy is defined by

\[S = - \text{Tr}\left(\rho \log \rho\right) = - \sum_i e_i\log(e_i)\,.\]

\(e_i\in[0,1]\) are the eigenvalues of the density matrix.

The density matrix will be only stored and diagonalized using the root rank. The function is only available if Petsc and Eigen are included.

Parameters:
  • vec – density matrix using Vec from Petsc

  • root – Designated MPI rank (Default: 0)

Returns:

von Neumann entropy S

Eigen::MatrixXcd danceq::get_density_matrix_in_eigen_from_PetscVec(const Vec &vec, const int32_t root = 0)

Returns the density matrix in an Eigen format from a Petsc Vec generated by the Lindbladian.

We represent the density matrix \(\rho\) of dimension \(d\) as a large vector \(v\) in Petsc with dimension \(d^2\) when we work with the Lindbladian OperatorType. This function returns the density matrix using Eigen:

\[\rho[i,j] = v[i\cdot d+j]\]

The density matrix will be only returned to the root rank. The function is only available if Petsc and Eigen are included.

Parameters:
  • vec – density matrix using Vec from Petsc

  • root – Designated MPI rank (Default: 0)

Returns:

density matrix

uint64_t danceq::get_dim_from_PetscVec(const Vec &v)

Returns dimension of a Vec object from Petsc.

This is only available if Petsc is included.

Returns:

dimension

std::vector<uint64_t> danceq::get_ownership_per_rank_PetscVec(const Vec &v)

Returns memory layout (MPI ownership) of a Vec object from Petsc.

This is only available if Petsc is included.

Returns:

ownership_per_rank

uint64_t danceq::get_dim_from_PetscMat(const Mat &mat)

Returns dimension of a Mat object from Petsc.

This is only available if Petsc is included.

Returns:

dimension

std::vector<uint64_t> danceq::get_ownership_per_rank_PetscMat(const Mat &mat)

Returns memory layout (MPI ownership) of a Mat object from Petsc.

This is only available if Petsc is included.

Returns:

ownership_per_rank