Amesos  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Amesos: An Interface to Serial and Parallel Direct Sparse Solver Packages.
AmesosLogo.jpeg

Table of Contents

Amesos provides an object-oriented interface to several direct sparse solvers, both sequential and parallel. Amesos is developed by M. Sala, K. Stanley, R. Hoekstra, T. Davis, and M. Heroux.

The main page of the Doxygen documentation is organized as follows:

Introduction

Amesos supports the following classes:

We refer to the Sandia report SAND-2004-2188 for a detailed overview of Amesos. A PDF version of this report can be found the Trilinos/packages/amesos/doc/AmesosReferenceGuide/AmesosReferenceGuide.pdf.

Copyright and licensing of the third party codes

Most of the Amesos classes are based on a third party code (that is, not distributed within Trilinos). Each third party code comes with its own copyright and/or licensing requirements. It is responsibility of the user to fulfill the requirements of each supported package' copyright.

Most of these third party codes are intended to be made available at no cost to users. Much of the copyright and licensing restrictions concern rights to modify, redistribute the code and generally include a request that credit be given in any papers which make use of their code. Please refer to the web page for the package that you are interested in for details.

Quick comparison of the Amesos classes

The following tables reports an overview of the Amesos classes.

Features Amesos_Lapack Amesos_Klu Amesos_Umfpack Amesos_Superlu Amesos_Scalapack
Default behavior Enabled Enabled Disabled Disabled Disabled
Package is written in FORTRAN77 C C C FORTRAN77
Communicator none none none none MPI
Processes used for Factorization/Solve Process 0 Process 0 Process 0 Process 0 Any
Distributed input matrix Yes Yes Yes Yes Yes
Unsymmetric matrices Yes Yes Yes Yes Yes
Requires BLAS Yes Yes Yes Yes Yes
Requires LAPACK Yes No No No No
Requires BLACS No No No No Yes
Requires ScaLAPACK No No No No Yes
Features Amesos_Pardiso Amesos_Taucs Amesos_Superludist Amesos_Mumps Amesos_Dscpack
Default behavior Disabled Disabled Disabled Disabled Disabled
Package is written in C C FORTRAN90 C
Communicator none none MPI FORTRAN MPI MPI
Processes used for Factorization/Solve Any Any Any Any (*) Any (**)
Distributed input matrix Yes Yes Yes Yes Yes
Unsymmetric matrices Yes No Yes Yes No
Requires BLAS Yes No Yes Yes Yes
Requires LAPACK Yes No No No No
Requires BLACS No No No Yes No
Requires ScaLAPACK No No No Yes No

(*) MUMPS required FORTRAN communicators. Some architectures (e.g., SGI) does not allow portable conversions from C/C++ communicator and FORTRAN. On these architectures, therefore, MPI_COMM_WORLD will be used by MUMPS. Instead, if C/C++ communicators can be converted to FORTRAN ones, then any number of prosesses can be used by Amesos_Mumps.

(**) DSCPACK requires a number of processes that is a power of 2.

Supported Matrix Formats.

The following table details the supported matrix formats for each class.

Amesos_Lapack Amesos_Klu Amesos_Umfpack Amesos_Pardiso Amesos_Taucs Amesos_Scalapack
Epetra_RowMatrix yes yes yes yes yes yes
Epetra_CrsMatrix yes yes yes yes yes yes
Epetra_VbrMatrix yes yes yes yes yes yes
Amesos_Superlu Amesos_Superludist Amesos_Mumps Amesos_Dscpack
Epetra_RowMatrix yes yes yes yes
Epetra_CrsMatrix yes yes yes yes
Epetra_VbrMatrix yes yes yes no

Installing Amesos

Each of the Amesos classes provides an interface to a third party direct sparse solver code. (Exception is KLU, whose sources are distributed within Amesos.) In order to install a particular class, you must first install the underlying direct sparse solver code.

Amesos is distributed through the Trilinos project, and can be downloaded from the web site http://trilinos.sandia.gov/download.

As all other Trilinos packages, Amesos is configured and built using the GNU autoconf and automake tools. To configure Amesos from the Trilinos top directory, a possible procedure is as follows. Let TRILINOS_HOME be a shell variable representing the location of the Trilinos source directory, and % the shell prompt sign. Let us suppose that we want to configure Amesos on a LINUX machine with MPI, with support for KLU and UMFPACK. Header files for UMFPACK are located in directory /usr/local/umfpack/include, while the library, called libumfpack.a is located in /usr/local/umfpack/lib. The configure like will look like:

% cd $TRILINOS_HOME
% mkdir LINUX_MPI
% cd LINUX_MPI
% ../configure \
--with-mpi-compilers \
--prefix=$TRILINOS_HOME/LINUX_MPI \
--enable-amesos \
--enable-amesos-klu \
--enable-amesos-umfpack \
--with-incdirs="-I/usr/local/umfpack/include" \
--with-ldflags="-L/usr/local/umfpack/lib" \
--with-libs="-lumfpack"
% make
% make install

Other flags may be required depending on the location of MPI, BLAS and LAPACK. The table below reports the architectures and compilers tested with each of the Amesos classes. "nightly" means that these tests are included in out routine (though not always nightly) testing. "yes" means that theses tests were run at least once.

This table may soon be replaced by a table which reports only the results of the routine, "nightly", tests. MANUAL TESTING will continue to list of architectures and compilers on which Amesos solvers have been tested manually.

architecture Amesos_Lapack Amesos_Klu Amesos_Umfpack Amesos_Pardiso Amesos_Taucs Amesos_Scalapack
LINUX (GNU), SERIAL yes yes yes yes yes no
LINUX (Intel), MPI yes yes yes no no no
LINUX (GNU), MPI nightly nightly yes yes yes nightly
SGI 64, MPI nightly nightly nightly no no nightly
DEC/Alpha MPI nightly nightly nightly no no no
MAC OSX, SERIAL no yes yes no yes no
MAC OSX, MPI nightly nightly no no no no
Sun OS, MPI nightly yes yes no no no
Sandia Cplant, MPI no yes yes no no no
architecture Amesos_Superlu Amesos_Superludist Amesos_Mumps Amesos_Dscpack
LINUX (GNU), SERIAL yes yes no yes
LINUX (Intel), MPI no no yes no
LINUX (GNU), MPI nightly nightly no nightly
SGI 64, MPI yes no nightly nightly
DEC/Alpha MPI no no no no
MAC OSX, SERIAL yes no no no
MAC OSX, MPI no no no no
Sun OS, MPI nightly no no nightly
Sandia Cplant, MPI yes yes no no

Once Amesos has been compiled, it can be tested using the example compare_solvers.exe, as follows:

% cd amesos/example
% make
% ./compare_solvers.exe

or, using MPI,

% mpirun -np 4 ./compare_solvers.exe

This will test all the supported solvers on an example matrix.

Example of Use of Amesos Classes

A simple fragment of code using Amesos can read as follows. First, we need to include the header files for Amesos:

#include "Amesos.h"

Note that these header files will not include the header files for the supported libraries (which are of course needed to compile the Amesos library itself). Now, let define the linear system matrix, the vector that will contain the solution, and the right-hand side as:

Epetra_RowMatrix* A; // linear system matrix
Epetra_MultiVector* LHS; // vector that will contain the solution
Epetra_MultiVector* RHS; // right-hand side

All Amesos object (derived from pure virtual class Amesos_BaseSolver) can be created using the function class Amesos, as follows:

Amesos_BaseSolver * Solver;
Amesos Factory;
char* SolverType = "Amesos_Klu";
Solver = Factory.Create(SolverType, *Problem);
Note
It is important to note that all available solvers can be selected simply by changing an input parameter. Some solvers are serial, other parallel; generally, each solver has its own matrix format. However, the user will still have the same interface to all of them.

The Factory object will create an Amesos_Klu object (if Amesos has been configure to support this solver). Factory.Create() returns 0 if the requested solver is not available. Parameter names are case-sensitive; misspelled parameters will not be recognized. Method Factory.Query() can be used to query the factory about the availability of a given solver:

char* SolverType = "Amesos_Klu";
bool IsAvailable = Factory.Query(SolverType);

Here, we simply recall that the parameters list can be created as

and parameters can be set as

List.set(ParameterName,ParameterValue);

Here, ParameterName is a string containing the parameter name, and ParameterValue is any valid C++ object that specifies the parameter value (for instance, an integer, a pointer to an array or to an object). Please consult the Amesos Reference Guide for more details. The Doxygen documentation of each class can also be of help.

Note that Problem is still empty. After setting the pointer to the linear system matrix, we can perform the symbolic factorization of the linear system matrix:

AMESOS_CHK_ERR(Solver->SymbolicFactorization());

(AMESOS_CHK_ERR is a macro, that checks the return error and, if not null, prints out a message and returns.) This phase does not require the numerical values of A, that can therefore be changed after the call to SymbolicFactorization(). However, the nonzero pattern of A cannot be changed.

The numeric factorization is performed by

AMESOS_CHK_ERR(Solver->NumericFactorization());

The values of RHS must be set before solving the linear system, which simply reads

Problem.SetLHS(LHS);
Problem.SetRHS(RHS);
AMESOS_CHK_ERR(Solver->Solve());

Should the user need to re-factorize the matrix, he or she must call NumericFactorization(). If the structure of the matrix is changed, he or she must call SymbolicFactorization(). However, it is supposed that the linear system matrix and the solution and right-hand side vectors are still defined with the same Epetra_Map.

Please consult the examples reported in the amesos/example subdirectory:

Browse all of Amesos as a single doxygen collection

You can browse all of Amesos as a single doxygen collection. Warning, this is not the recommended way to learn about Amesos software. However, this is a good way to browse the directory structure of amesos, to locate files, etc.

Python support

Amesos has an interface to PyTrilinos, and therefore all the enabled third-party libraries can be used within Python to solve problems defined using Python/Epetra objects. Please refer to the PyTrilinos documentation http://trilinos.sandia.gov/packages/pytrilinos/overview.html for more details.

Thyra support

The Amesos to Thyra Linear Solver Adapters in Stratimikos take Thyra::LinearOpBase objects that wrap Epetra_Operator objects and turn them into Thyra::LinearOpWithSolveBase objects which can then be used to solve linear systems using Amesos_BaseSolver implementations. Please refer to the Stratimikos documentation http://trilinos.sandia.gov/packages/stratimikos/documentation.html for more details.

Known bugs and future developments

Incompatibility between old versions of SuperLU

Versions of Superlu and Superludist prior to August 2005 exhibit some mutual incompatibility as well as incompatbility with KLU and UMFPACK.

Depending on the version of SuperLU and SuperLU_DIST, you may or may not configure Amesos with support for both. This is due to conflicts in the header files of SuperLU/SuperLU_DIST.

bugs

<p> Amesos_Dscpack has memory leaks, prints to standard out and may 

crash if SymbolicFactorization() is called multiple times on the same Amesos object.

Amesos packages other than Amesos_Klu do not accept non-standard maps (e.g. one with indices {1,2,5,8,9, ... } )

Amesos packages other than Amesos_Klu do not accept range and domain maps which differ from the matrix map

Superludist fails on some matrices

Error handling on singular and near singular matrices is inconsistent

Final Remarks

This page is maintained by Marzio Sala, SNL 9214 and Ken Stanley.

This page last updated 26-Aug-05.