Usage

This page introduces the basic usage with the predefined classes and the generic setup of the code.

Operator Types

The simplest way to include the DanceQ library is the following:

/* Maximal number that can be encoded in the basis class */
#define MaxSites 64

/* Local Hilbert space dimension */
#define Q 2

/* ScalarType for computation */
#define ScalarType std::complex<double>

#include "DanceQ.h"
using namespace danceq;

Including the header file and specifying the maximal number of sites that are potentially encoded, MaxSites (Default: 64), the local Hilbert space dimension, Q (Default: 2), and the scalar type, ScalarType (Default: std::complex<double>) defines different Operator class. That is all!

The header file:

defines Lindbladian and Hamiltonian versions of the Operator class as well as the corresponding matrix-free shells:

namespace danceq{

    /* Definition of the Hamiltonian class from the BasisU1 class */
    using Hamiltonian_U1=danceq::internal::Operator<Basis,ScalarType,danceq::Hamiltonian>;

    /* Definition of the Hamiltonian class from the State class */
    using Hamiltonian_NoU1=danceq::internal::Operator<State,ScalarType,danceq::Hamiltonian>;

    /* Definition of the Lindbladian class from the BasisU1 class */
    using Lindbladian_U1=danceq::internal::Operator<Basis,ScalarType,danceq::Lindbladian>;

    /* Definition of the Lindbladian class from the State class */
    using Lindbladian_NoU1=danceq::internal::Operator<State,ScalarType,danceq::Lindbladian>;

    /* Definition of the ShellMatrix class from the Operator class */
    using SellMatrix_Hamiltonian_NoU1=danceq::internal::ShellMatrix<Hamiltonian_NoU1>;

    /* Definition of the ShellMatrix class from the Operator class */
    using SellMatrix_Lindbladian_U1=danceq::internal::ShellMatrix<Lindbladian_U1>;

    /* Definition of the ShellMatrix class from the Operator class */
    using SellMatrix_Lindbladian_NoU1=danceq::internal::ShellMatrix<Lindbladian_NoU1>;

    /* Definition of the ShellMatrix class from the Operator class */
    using SellMatrix_Hamiltonian_U1=danceq::internal::ShellMatrix<Hamiltonian_U1>;

    /* Definition of the Vector class */
    using Vector=danceq::internal::Vector<ScalarType>;

    /* Definition of the SparseMatrix class */
    using SparseMatrix=danceq::internal::SparseMatrix<ScalarType>;

}

These classes can be used within the code by simply providing the system size and particle number (if number conservation is desired):

/* everything is defined in the namespace danceq */
using namespace danceq;

/* system size */
uint64_t L = 10;

/* particle number */
uint64_t n = 5;

/* Lindbladian with number conservation */
Lindbladian_U1 = L_U1(L,n)

/* Lindbladian without number conservation */
Lindbladian_NoU1 = L_NoU1(L)

/* Hamiltonian with number conservation */
Hamiltonian_U1 = H_U1(L,n)

/* Hamiltonian without number conservation */
Hamiltonian_NoU1 = H_NoU1(L)

Additionally, we have defined the ShellMatrix to perform the matrix-free multiplication. Each type has an associated ShellMatrix in the danceq namespace: SellMatrix_Lindbladian_U1, SellMatrix_Lindbladian_NoU1, SellMatrix_Hamiltonian_U1, and SellMatrix_Hamiltonian_NoU1. Furthermore, we provide a SparseMatrix container, and our own implementation of a Vector class that is compatible with openMP and MPI.

Hierarchy of classes

Before we look at the different examples we introduce the most important classes included in DanceQ:

The State class defines the Basis class which defines the Operator class.

Class definitions

\(\qquad\qquad\qquad\) State class \(\quad\stackrel{\mathrm{defines}}{\Longrightarrow}\quad\) Basis class \(\quad\stackrel{\mathrm{defines}}{\Longrightarrow}\quad\) Operator

There are two important template parameters that have to be set at compile time:

  • N_max_site: The maximal number of sites. (Default: 128)

  • Q: The local Hilbert space dimension. (Default: 2)

We use using to define the classes:

/* The operator class includes all files */
#include "Operator.h"

/* Setting the maximal number of sites and local Hilbert space Q */
#define MaxSites 256
#define Q 2

/* Definition of the State class */
using State=danceq::internal::State<MaxSites,Q>;

/* Definition of the Basis class */
using Basis=danceq::internal::BasisU1<State>;

/* Definition of the Operator class with double */
using MyOperator=danceq::internal::Operator<Basis,double,Hamiltonian> ;

Alternatively, you can define the Operator class straight from the State class:

using Hamiltonian_NoU1=danceq::internal::Operator<State,double>;

In this case, the operator is acting on the Hilbert space of dimension \(Q^L\).

Troubleshooting

Compiling cpp code can be annoying but a quick google search solves most problems. Remember to set all the paths correctly! Once the code is compiled but it the program crashes or returns non-sense, you can use the debugging option that is available in all classes:

#define dbug_level 10; // 10 maximal debugging, 0 is no debugging

Note

Many functions are returning an error_code which should be 0!

If you can not solve the problem, we are more than happy to help out! Please contact us via mail.