IFPACK  Development
 All Classes Namespaces Files Functions Variables Enumerations Friends Pages
Trilinos/IFPACK: Object-oriented Algebraic Preconditioner Package
IFPACKLogo.gif

Table of Contents

Introduction

IFPACK, developed by Marzio Sala (ETHZ/D-INFK) and Micheal Heroux (SNL 9214), provides a suite of object-oriented algebraic preconditioners for the solution of preconditioned iterative solvers. IFPACK constructors expect an Epetra_RowMatrix object for construction. IFPACK is part of the Trilinos Solver Project and IFPACK object interacts well with other Trilinos classes. In particular, IFPACK can be used as a preconditioner for AztecOO.

Most IFPACK preconditioners can be rewritten as additive Schwarz methods, of overlapping domain decomposition type. The user can adopt a minimal-overlap (that is, zero-row overlap), or ask IFPACK to extend the overlap. The resulting preconditioner reads:

\[ P_{IFPACK}^{-1} = \sum_{i = 0}^{NumProcs-1} P_i A_i^{-1} R_i \]

where \(R_i\) is the restriction operator from the global vector, to the overlapping subdomain \(i\), and \(P_i\) is the prolongator operator. \(P_i\) is generally the transpose of \(R_i\) (in which case the resulting preconditioner is symmetic). It is assumed that each subdomain is assigned to a different processor.

A key component of the previous formula is a strategy to apply \(A_i^{-1}\). Using IFPACK, this can be defined by one of the following:

More precisely, the choices for \(A_i^{-1}\) are:

Note
IFPACK preconditioners can be used as smoothers for multilevel solvers, like ML.

Download IFPACK

IFPACK can be downloaded from the web page http://trilinos.sandia.gov/download/

Configuring IFPACK

IFPACK is configured with autotools. IFPACK is enabled by default if configured at the Trilinos level (please refer to the Trilinos documentation for more details). This section briefly recalls some of the configure parameters that affect the IFPACK compilation. For a complete list of parameters, please type

% $TRILINOS_HOME/packages/ifpack/configure --help

An (incomplete) list of parameters recognized by the IFPACK configure is reported below.

Note
IFPACK cannot be build without Epetra. Most of IFPACK preconditioners cannot be used without the Teuchos library. Most of the examples requires triutils and AztecOO.

Examples of Usage

This section details how to use the IFPACK Factory class.

Probably, the easiest way to use IFPACK is through the function class Ifpack. This is a factory class, that contains only one method, Create(). A call to Create() returns a pointer to a newly-created IFPACK object. The user is responsible of deleting this object.

An example of usage is reported by the following code:

#include "Ifpack.h"
...
Ifpack Factory;
Epetra_RowMatrix* A; // A is already FillComplete()'d
// but its values can still change (not its structure)
string PrecType = "Amesos"; // exact solve on each subdomain
int OverlapLevel = 1; // one row of overlap among the processes
Ifpack_Preconditioner* Prec = Factory.Create(PrecType, A, OverlapLevel);
assert (Prec != 0);
IFPACK_CHK_ERR(Prec->SetParameters(List));
IFPACK_CHK_ERR(Prec->Initialize());
IFPACK_CHK_ERR(Prec->Compute());
// we can query the status of the preconditioner
assert (Prec->IsInitialized() == true);
assert (Prec->IsComputed() == true);
// information about the preconditioner
// 1.- number of initialization phases
cout << Prec->NumInitialize() << endl;
// 2.- number of computation phases
cout << Prec->NumCompute() << endl;
// 3.- number of applications of the prconditioner
cout << Prec->NumApplyInverse() << endl;
// we can query for flops and CPU time as well, consult the
// documentation of Ifpack_Preconditioner, or simply write
cout << *Prec;
// now destroy the preconditioner
delete Prec;

Although the Ifpack factory is appropriate for most users, it is possible to create IFPACK preconditioners by handling directly an Ifpack_AdditiveSchwarz object. The preconditioner of the previous example is equivalent to the one created in the following fragment of code:

#include "Ifpack_AdditiveSchwarz.h"
...
Epetra_RowMatrix* A;
int OverlapLevel = 1;
Ifpack_AdditiveSchwarz Prec<Ifpack_Amesos> (A, OverlapLevel);

Another example is the following, where the local solve is a block Jacobi method, which uses dense containers (LAPACK).

#include "Ifpack_AdditiveSchwarz.h"
#include "Ifpack_DenseContainer.h"
#include "Ifpack_BlockJacobi.h"
...
Epetra_RowMatrix* A;
int OverlapLevel = 1;
Ifpack_AdditiveSchwarz Prec<Ifpack_BlockJacobi<Ifpack_DenseContainer> > (A, OverlapLevel);

Using sparse containers and Amesos the code may look as follows.

#include "Ifpack_AdditiveSchwarz.h"
#include "Ifpack_SparseContainer.h"
#include "Ifpack_BlockJacobi.h"
#include "Ifpack_Amesos.h"
...
Epetra_RowMatrix* A;
int OverlapLevel = 1;
Ifpack_AdditiveSchwarz Prec<Ifpack_BlockJacobi<Ifpack_SparseContainer<Ifpack_Amesos> > > (A, OverlapLevel);

Other examples are reported in the example subdirectory with several comments:

List of Supported Parameters

The complete list of supported parameters is reported below.

The table below reports the C++ type and the default value of each parameter.

Parameter nameC++ typeDefaut value
"relaxation: type" string "Jacobi"
"relaxation: sweeps" int "1"
"relaxation: damping factor" double "1.0"
"relaxation: min diagonal value" double "0.0"
"relaxation: zero starting solution" bool "true"
"relaxation: backward mode" bool "false"
"relaxation: use l1" bool "false"
"relaxation: l1 eta" double "1.5"
"relaxation: number of local smoothing indices" int "A->NumMyRows()"
"relaxation: local smoothing indices" int* "0"
"partitioner: type" string "greedy"
"partitioner: overlap" int "0"
"partitioner: local parts" int "1"
"partitioner: root node" int "0"
"partitioner: use symmetric graph" bool "true"
"amesos: solver type" string "Amesos_Klu"
"fact: absolute threshold" double "0.0"
"fact: relative threshold" double "0.0"
"fact: relax value" double "0.0"
"schwarz: combine mode" Epetra_CombineMode "Zero"
"schwarz: compute condest" bool "true"
"schwarz: reordering type" string "none"
"schwarz: filter singletons" string "false"
Note
Whenever possible, the Teuchos::ParameterList given in input to SetParameters() is only parsed. There are cases, however, where the solver needs to store a copy of this list, because it will be used in a later stage. This happens, for example, in block relaxation schemes that use Amesos as local solver. As the the list given in input to SetParameters() is parsed or copied, this list can go out of scope after this function returns.

Brief Description IFPACK Classes

IFPACK classes are divided into classes for users and classes for developers.

Users' Tools

Developers' Tools

Error Table

The following table defines the class of IFPACK errors. Return code should be checked using the IFPACK_CHK_ERR() macro. In general terms, we follow this convention:

Error CodeMeaning
-1Generic Error (called method or function returned an error)
-2Input data not valid (wrong parameter, out-of-bounds, wrong dimensions, matrix is not square,...)
-3Data has not been correctly pre-processed
-4Problem encountered during application of the algorithm (division by zero, out-of-bounds, ...)
-5Memory allocation error
-98Feature is not supported
-99Feature is not implemented yet (check Known Bugs and Future Developments, or submit a bug)

Known Bugs and Future Developments

Copyright

Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
license for use of this work by or on behalf of the U.S. Government.

This library is free software; you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as
published by the Free Software Foundation; either version 2.1 of the
License, or (at your option) any later version.

This library is distributed in the hope that it will be useful, but
WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
Lesser General Public License for more details.

You should have received a copy of the GNU Lesser General Public
License along with this library; if not, write to the Free Software
Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
USA

Browse all of Ifpack as a single doxygen collection

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