ShellMatrix

The ShellMatrix class is a core feature of DanceQ. By providing an Operator, it allows a matrix-free matrix-vector multiplication using the Vector class. Not storing the matrix yields an increased runtime, but reduces the memory significantly allowing the access to much larger systems. In particular, acting on a Hilbert space spanned by \(L\) sites, utilizing the matrix-free routines of the ShellMatrix reduces the memory cost by \(L\) for ground-state searches. Computing the ground state for a system with \(L=42\) at half-filling (\(S=1/2\)) requires a large computing center with up to \(160\) TB. In contrast, the matrix-free computation reduces the memory cost to \(4\) TB which becomes now possible with a few computing nodes.

The routines are highly optimized and support distributed memory with MPI and shared memory with openMP. Note, that no third-party libraries like Eigen or Petsc/Slepc are required. If MPI has been found, it is used. If this is not the case and openMP is available, the program is parallelized using openMP. If no parallelization option is given, the code is executed serially.

/* Retrieving the ShellMatrix from the Operator class */
auto H = H_operator.create_ShellMatrix();

/* Retrieving the Vector from the ShellMatrix class */
auto v = H.create_Vector();

/* Random initialization of the Vector */
v.make_random();

/* Matrix-free multiplication */
auto w = H*v;

/* Output */
w.print();

/* Expectation value */
auto val0 = H.get_expectation_value(v,v);

/* val1 equals val0 */
auto val1 = w*v;

Note

The syntax is equivalent for MPI and openMP!

A full example for distributed and shared memory is available in the repository as discussed here. It computes multiple ground states.

Petsc wrapper

The following functions are only available if Petsc has been included. The ShellMatrix class is internally used by Petsc to create a MATSHELL object that behaves like a normal matrix and can be used for various functions within Petsc/Slepc. It works with the Petsc vector structure Vec which distributes the memory over all MPI ranks. See this exmaple for more information. The Petsc MATSHELL can be simply retrieved from an Operator via create_PetscShellMatrix():

/* MATHELL Petsc object using the ShellMatrix class */
auto H_MATSHELL = H.create_PetscShellMatrix();

The necassary functions are listed in:

Datatypes

The underlying scalar scalartype, its real part scalartype_real, and the State class state_class can be retrieved from an instance:

using danceq::internal::ShellMatrix::scalartype = typename Op::scalartype

Scalar.

The underlying scalar defined from the Operator class can be retrieved from the ShellMatrix class:

ShellMatrix::scalartype A;  

using danceq::internal::ShellMatrix::scalartype_real = decltype(std::real(std::declval<scalartype>()))

Real part of the scalar.

The real part of the underlying scalar defined from the Operator class can be retrieved from the ShellMatrix class:

ShellMatrix::scalartype_real A;  

This is equivalent to scalartype if it is real. If scalartype is complex, e.g., std::complex<double>, scalartype_real is double.

using danceq::internal::ShellMatrix::state_class = typename Op::state_class

State class.

State class can be retrieved from the ShellMatrix class:

// This equivalent to: State state;
ShellMatrix::state_class state;  

Class content