ALPS Home Libraries License Support People ALPS Web Site

PrevUpHomeNext

Chapter 13. IETL

Robin Jäger

Distributed under the Boost Software License, Version 1.0. (See accompanying file LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt )

Table of Contents

IETL
IETL Concepts and Interfaces
Jacobi-Davidson
Reference
Header <alps/src/ietl/bandlanczos.h>
Header <alps/src/ietl/bicgstabl.h>
Header <alps/src/ietl/cg.h>
Header <alps/src/ietl/complex.h>
Header <alps/src/ietl/fmatrix.h>
Header <alps/src/ietl/gmres.h>
Header <alps/src/ietl/ietl2lapack.h>
Header <alps/src/ietl/interface/blas.h>
Header <alps/src/ietl/interface/blitz.h>
Header <alps/src/ietl/interface/mtl.h>
Header <alps/src/ietl/interface/ublas.h>
Header <alps/src/ietl/interface/valarray.h>
Header <alps/src/ietl/inverse.h>
Header <alps/src/ietl/iteration.h>
Header <alps/src/ietl/jacobi.h>
Header <alps/src/ietl/jd.h>
Header <alps/src/ietl/krylov_wrapper.h>
Header <alps/src/ietl/lanczos.h>
Header <alps/src/ietl/matrix.h>
Header <alps/src/ietl/power.h>
Header <alps/src/ietl/rayleigh.h>
Header <alps/src/ietl/simple_arnoldi.h>
Header <alps/src/ietl/tmatrix.h>
Header <alps/src/ietl/traits.h>
Header <alps/src/ietl/vectorspace.h>

IETL

IETL Concepts and Interfaces

Here we define the requirements on the types necesary to use the IETL algorithms:

The vector space concept

Linear eigensolvers mathematically rely on the concept of a vectorspace, which also defines the necessary functionality for the IETL. Going beyond the ITL we expand on the concept of a vectorspace and find the necessity to require a vector space class to be passed to the algorithms. Let

  • VS be a class modeling the concept of vectorspace
  • vs be an object of type VS

then the following types have to be defined:

Type

Documentation

ietl::vectorspace_traits<VS>::vector_type

the type of vectors in the vector space

ietl::vectorspace_traits<VS>::scalar_type

the type of scalars in the vector space

ietl::vectorspace_traits<VS>::magnitude_type

a scalar type appropriate for storing norms, usually the same as scalar_type for real types, or the corresponding real type for complex types.

ietl::vectorspace_traits<VS>::size_type

an integral type large enough to store the dimension of the vector space

and the following expression have to be defined for the vector space object:

Expression

return type

Documentation

new_vector(vs)

ietl::vectorspace_traits<VS>::vector_type

creates a new vector of the vector space

vec_dimension(vs)

ietl::vectorspace_traits<VS>::size_type

returns the dimension of the vector space

project(x,vs)

void

projects the vector x into the vector space

The new_vector function is necessary since we do not know the types and number of arguments needed for the constructor of a vector. The ITL assumes there is a constructor taking only the vector size as argument but this is insufficient in cases where more than the size is needed in the constructor, such as:

  • multi-dimensional arrays
  • distributed arrays taking additional information, requiring as the data layout in the constructor

The project function is required if the actual vector space is a subspace of the representation space. E.g. if special boundary conditions are imposed on an array. For the most common case, where the vector models the concept of an STL Container we provide two model implementations: vectorspace and wrapper_vectorspace.

Requirements on the matrix

Iterative algorithms require the matrix (linear operator) of which eigenvalues and eigenvectors are to be calculated only in the form of matrix-vector products. Hence only a single function is required. Let

  • VS be a class modeling the concept of vectorspace
  • x, y be an vectors in the vector space ( of type ietl::vectorspace_traits<VS>::vector_type )
  • a be the matrix

then the following expression needs to be defined

Expression

return type

Documentation

ietl::mult(a,x,y)

void

calculates the matrix-vector product y=a*x

For all uBlas matrix types, the ietl::mult function with boost::numeric:ublas::vector as vector type is implemented in the header ietl/interface/ublas.h .

Requirements on the vector

The vector type needs to fulfill the following requirements.

  • VS be a class modeling the concept of vectorspace
  • V be the type ietl::vectorspace_traits<VS>::vector_type
  • S be the type ietl::vectorspace_traits<VS>::scalar_type
  • M be the type ietl::vectorspace_traits<VS>::magnitude_type
  • x, y be a vectors in the vector space, of type V
  • t be a scalar, of type S
  • g be a model of a generator producing values of type S

Then the following expressions need to be defined for all IETL algorithms.

Expression

return type

Documentation

ietl::generate(x,g)

void

fills the vector x with numbers from the generator g. For a standard container this can be implemented as std::generate(x.begin(),x.end(),g);

std::swap(x,y)

void

swaps the two vectors x and y

ietl::dot(x,y)

S

calculates the scalar product of the two vectors x and y.

ietl::two_norm(x)

M

calculates the 2-norm of the vector x. This is equivalent to sqrt(ietl::dot(x,x)).

ietl::copy(x,y)

void

a deep copy y = x. Modifications of x after the call are not allowed to modify y.

y = x

const V&

a (possibly shallow) copy

x *= t

const V&

x /= t

const V&

x += y

const V&

x += t*y

const V&

x -= t*y

const V&

x = t*y

const V&

The necessary functions are implemented in the header ietl/interface/ublas.h for ublas vectors, and in ietl/interface/blitz.h for Blitz++ arrays.

Requirements on the standard iteration control object

Iterative algorithms run until the desired accuracy is reached. To control the termination of the iteration, the IETL algorithms take an iteration control object. Let it be the iteration control object, lambda the best estimate for the eigenvalue and residual be the residual ||A v - lambda v|| then the following expressions need to be defined:

Expression

return type

Documentation

++it

void

is called at the start of a new iteration

it.finished(residual, lambda)

bool

returns true if the desired accuracy is reached and the iteration can be terminated, returns false otherwise. residual and lambda are of the magnitude_type and are the residual and and latest estimate of the eigenvalue respectively

A model implementation basic_iteration is provided in ietl/iteration.h.

Some algorithms require more complex iteration control objects. They will be discussed with the algorithms.

Requirements on the Lanczos iteration control object

Iterative algorithms run until the desired accuracy is reached. To control the termination of the iteration, the IETL algorithms take an iteration control object. Let it be the iteration control object, lambda the best estimate for the eigenvalue and residual be the residual ||A v - lambda v|| then the following expressions need to be defined:

Expression

return type

Documentation

++it

void

is called at the start of a new iteration

it.finished(tmatrix)

bool

returns true if the desired accuracy is reached and the iteration can be terminated, returns false otherwise. The T-matrix of the Lanczos iterations is passed as argument.

Model implementation are provided ietl/iteration.h.

Jacobi-Davidson

Generic implementation of the Jacobi-Davidson algorithm for hermitian hamiltonians.

Usage

This eigensolver can be used with any linear algebra library fullfilling the requirements of the IETL.

[Note] Note

Note that a function ietl::mult(A, x, b), calculating the matrix-vector product b = A x, must be defined and has to be declared before including the header ietl/jd.h.

If you are using the boost::ublas library, ietl/interface/ublas.h holds the proper function interface.

The constructor has the following syntax:

jd(const MATRIX& A, VS& vspace, size_t verbose = 0)

  • A is a matrix fullfilling the requirements.
  • VS is a vectorspace.
  • verbose is optional,
    • 0, for errors only
    • 1, print message when converged
    • 2, print current iteration and residual

The class ietl::jd basically provides two functions:

These are implementations of the algorithm templates jdqr and jdqrhar. The first can be used with a constant preconditioner.

Parameter

Description

Example

Info

IT& iter

Iterator object

ietl::jd_iteration(N)

GEN& gen

Random number generator for generating a starting vector

boost::lagged_fibonacci607()

size_t k

Number of eigenpairs to be calculated

3

PREC& K

Constant Preconditioner

example2

optional for eigensystem()

SOLV& solver

Krylov solver

ietl::gmres_wrapper()

bool search_highest = false

On default search for the lowest eigenpairs

true

only for eigensystem()

real_t tau

Target value near which eigenpairs are searched

0.

only for eigensystem_harmonic()

Iterator

The iterator jd_iteration contains all required parameters.

Parameter

Description

Default

reasonable constraints

size_t max_iter

max. number of iterations

~N

size_t m_min

size of the vectorset after a restart

10

>= 3

size_t m_max

max. size of the vectorset

20

> m_min

T reltol

relative tolerance

sqrt(std::numeric_limits<double>::epsilon())

>= default

T abstol

absolute tolerance

sqrt(std::numeric_limits<double>::epsilon())

>= default

Krylov Solver

For the Jacobi-Davidson algorithm one needs a Krylov solver to approximately solve the correction equation (I - QQ*)(A - \theta I)(I - QQ*) t = -r. The krylov-solver has to be wrapped in a function object. Already provided are GMRES and BiCGSTAB(L).

The function object must have the following syntax:

class krylov_wrapper {
   public:
       template <class VECTOR, class MATRIX, class REAL>
       //solve x from A x = b
       VECTOR operator() ( const MATRIX & A, // deflated matrix
                           const VECTOR & b,
                           const VECTOR & x0, // starting vector
                           REAL absolute_tolerance )
};

The class detail::deflated_matrix defines a deflated matrix-vector product.

Example 1 A simple Hamiltonian

In this example a hamiltonian of a simple one particle system with periodic boundary conditions is shown, it does not store a matrix. The lowest eigenvalue is -2 and all interior eigenvalues are two-fold degenerate. Here is part of the code. To see the full sourcecode, look at jacobidavidson1.cpp.

namespace ietl {
    template<class Vector>
    void mult(Hamiltonian<Vector> const & H, Vector& x, Vector& y)
    {
        H.mult(x, y);
    }
}

#include <ietl/jd.h>

typedef boost::numeric::ublas::vector<double> vector_t;
typedef ietl::vectorspace<vector_t> vecspace_t;
typedef boost::lagged_fibonacci607 gen_t;
typedef boost::lagged_fibonacci44497 gen2_t;
typedef Hamiltonian<vector_t> ham_t;

int main(int argc, char **argv)
{
    int L = 100;    // dimension of the Hamiltonian
    
    vecspace_t vs(L);
    ham_t H(L);

    // random generator to create a starting vector x0
    // to find an eigenvector v, we require dot(x0,v) != 0
    gen_t gen;
    gen2_t gen2;
    
    int n_evals = 10;   // number of eigenpairs to be calculated
    int max_iter = L*10; // maximal number of iterations

    int m_max = 40;
    int m_min = 20;
 
    // tolerance
    double rel_tol = sqrt(std::numeric_limits<double>::epsilon());
    double abs_tol = rel_tol;

    // maximal iterations for the correction equation solver
    unsigned max_cor_iter = 10;
    // on default 5 steps are used

    ietl::jd_iteration<double> iter(max_iter, m_min, m_max, rel_tol, abs_tol);
    ietl::jd<Hamiltonian<vector_t>, vecspace_t> jd(H, vs /*, 2 for verbose mode */ );

    // the correction equation solver must be an function object
    ietl::gmres_wrapper gmres(max_cor_iter);
    
    // to find degenerated (or not well seperated) eigenvalues in the right order
    // the correction equation has to be solved exactly, this is very expensive. 
    // here we use different starting vectors instead.
    try {
            jd.eigensystem(iter, gen, n_evals/2, gmres);
            jd.eigensystem(iter, gen2, n_evals-n_evals/2, gmres);
    } catch (std::runtime_error& e) {
            cerr << "Error in eigenpair calculation: " << e.what() << "\n";
            return 1;
    }
    // if you use std=c++0x gmres is used as default

    std::vector<double> evals = jd.eigenvalues(); // copy eigenvalues
    // cout << "First eigenvector: \n" << jd.eigenvector(0) << "\n";
    std::cout.precision(10);
    std::sort(evals.begin(), evals.end());
    cout << "Sorted Eigenvalues: \n";
    std::copy(evals.begin(), evals.end(), std::ostream_iterator<double>(cout, "\n"));
}
    

Example 2 Jacobi preconditioning for a random matrix

Here we calculate the eigenvalues of a diagonal dominant random matrix. First we do this without preconditioning, then with preconditioning comparing the iterations and time used until convergence is reached. Again here is only the code skeleton of the example jacobidavidson2.cpp. The preconditioner has to define a function ietl::mult(K, x, r) which approximates r ~= K x.

The preconditioner:

    template <class MATRIX, class SCALAR>
    class jacobi_prec {
        public:
            jacobi_prec( const MATRIX& a, SCALAR l )
             : A_(a), lambda_(l) {}

            template <class VECTOR>
            void mult( const VECTOR & x, VECTOR & y ) const
            {
                assert(A_.size2() == x.size());
                if(y.size() != x.size())
                    y = VECTOR(x.size());

                for(unsigned i = 0; i < x.size(); ++i)
                    y(i) = x(i) / ( A_(i,i) - lambda_);
            }
        private:
            const MATRIX& A_;
            SCALAR lambda_;
    };

namespace ietl {
    template <class MATRIX, class SCALAR, class VECTOR>
    void mult( jacobi_prec<MATRIX,SCALAR> const & K, const VECTOR & x, VECTOR & y )
    {
        K.mult(x, y);
    }
}

The function call:

// create jacobi preconditioner
jacobi_prec<matrix_t,double> K(A, lambda);

jd_test.reset();
// correction equation solver steps are more expensive with preconditioning,
// but fewer are needed to have 'good' convergence
ietl::gmres_wrapper solver2(3);

std::cout << "solve with jacobi preconditioning...";
std::cout.flush();
clock.restart();
try{
    jd_test.eigensystem(iter3, gen, k, K, solver2);
}
catch (std::runtime_error& e) {
    std::cerr << "Something went wrong: " << e.what() << "\n";
}

Last revised: April 27, 2012 at 09:04:07 GMT


PrevUpHomeNext