hypre/docs/usr_parasails.tex
2001-07-26 16:44:14 +00:00

147 lines
6.2 KiB
TeX

%==========================================================================
\section{ParaSails}
ParaSails is a parallel implementation of a sparse approximate inverse
preconditioner, using {\em a priori} sparsity patterns and least-squares
(Frobenius norm) minimization. Symmetric positive definite (SPD) problems
are handled using a factored SPD sparse approximate inverse. General
(nonsymmetric and/or indefinite) problems are handled with an
unfactored sparse approximate inverse. It is also possible to
precondition nonsymmetric but definite matrices with a factored, SPD
preconditioner.
ParaSails uses {\em a priori} sparsity patterns that are patterns of powers
of sparsified matrices. ParaSails also uses a post-filtering technique
to reduce the cost of applying the preconditioner.
In advanced usage not described here, the pattern of the
preconditioner can also be reused to generate preconditioners for different
matrices in a sequence of linear solves.
For more details about the ParaSails algorithm, see \cite{Chow:1999:APS}.
\subsection{Synopsis}
\begin{display}
\begin{verbatim}
#include "HYPRE_parcsr_ls.h"
int HYPRE_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver,
int symmetry);
int HYPRE_ParaSailsSetParams(HYPRE_Solver solver,
double thresh, int nlevel, double filter);
int HYPRE_ParaSailsSetup(HYPRE_Solver solver, HYPRE_ParCSRMatrix A,
HYPRE_ParVector b, HYPRE_ParVector x);
int HYPRE_ParaSailsSolve(HYPRE_Solver solver, HYPRE_ParCSRMatrix A,
HYPRE_ParVector b, HYPRE_ParVector x);
int HYPRE_ParaSailsStats(HYPRE_Solver solver);
int HYPRE_ParaSailsDestroy(HYPRE_Solver solver);
\end{verbatim}
\end{display}
The accuracy and cost of ParaSails are parameterized by the real {\em thresh}
and integer {\em nlevels} parameters,
$0 \le {\it thresh} \le 1$, $0 \le {\it nlevels}$.
Lower values of {\em thresh}
and higher values of {\em nlevels} lead to more accurate, but more expensive
preconditioners. More accurate preconditioners are also more expensive
per iteration. The default values are ${\it thresh} = 0.1$
and ${\it nlevels} = 1$. The parameters are set using
{\tt HYPRE\_ParaSailsSetParams}.
Mathematically, given a symmetric matrix $A$, the pattern of the
approximate inverse is the pattern of $\tilde{A}^m$ where $\tilde{A}$
is a matrix that has been sparsified from $A$. The sparsification
is performed by dropping all entries in a symmetrically diagonally scaled $A$
whose values are less than {\em thresh} in magnitude. The parameter
{\em nlevel} is equivalent to $m+1$.
Filtering is a post-thresholding procedure.
For more details about the algorithm, see \cite{Chow:1999:APS}.
The storage required for the ParaSails preconditioner depends on
the parameters {\em thresh} and {\em nlevels}. The default parameters
often produce a preconditioner that can be stored in less than the
space required to store the original matrix.
ParaSails does not need a large amount of intermediate storage in
order to construct the preconditioner.
\subsection{Interface functions}
A ParaSails solver {\tt solver} is returned with
\begin{display}
\begin{verbatim}
int HYPRE_ParaSailsCreate(MPI_Comm comm, HYPRE_Solver *solver,
int symmetry);
\end{verbatim}
\end{display}
where {\tt comm} is the MPI communicator.
The value of {\tt symmetry} has the following meanings, to indicate
the symmetry and definiteness of the problem, and to specify the
type of preconditioner to construct:
\begin{center}
\begin{tabular}{|c|l|} \hline
value & meaning \\ \hline
0 & nonsymmetric and/or indefinite problem, and nonsymmetric preconditioner \\
1 & SPD problem, and SPD (factored) preconditioner \\
2 & nonsymmetric, definite problem, and SPD (factored) preconditioner \\
\hline
\end{tabular}
\end{center}
For more information about the final case, see section \ref{nearly}.
Parameters for setting up the preconditioner are specified using
\begin{display}
\begin{verbatim}
int HYPRE_ParaSailsSetParams(HYPRE_Solver solver,
double thresh, int nlevel, double filter);
\end{verbatim}
\end{display}
The parameters are used to specify the sparsity pattern and filtering value
(see above), and are described with suggested values as follows:
\begin{center}
\begin{tabular}{|c|c|c|c|c|l|} \hline
parameter & type & range & sug. values & default & meaning \\ \hline
{\tt nlevel} & integer & ${\tt nlevel} \ge 0$ & 0, 1, 2 & 1 & $m={\tt nlevel}+1$\\
\hline
{\tt thresh} & real & ${\tt thresh} \ge 0$ & 0, 0.1, 0.01 & 0.1 & {\em thresh} $=$ {\tt thresh}\\
& & ${\tt thresh} < 0$ & -0.75, -0.90 & & {\em thresh} selected automatically\\
\hline
{\tt filter} & real & ${\tt filter} \ge 0$ & 0, 0.05, 0.001 & 0.05 & filter value $=$ {\tt filter}\\
& & ${\tt filter} < 0$ & -0.90 & & filter value selected automatically\\
\hline
\end{tabular}
\end{center}
When ${\tt thresh} < 0$, then a threshold is selected such that
$-{\tt thresh}$ represents the fraction of the nonzero elements
that are dropped. For example, if ${\tt thresh} = -0.9$ then
$\tilde{A}$ will contain approximately ten percent of the nonzeros
in $A$.
When ${\tt filter} < 0$, then a filter value is selected such that
$-{\tt filter}$ represents the fraction of the nonzero elements
that are dropped. For example, if ${\tt filter} = -0.9$ then
approximately 90 percent of the entries in the computed approximate
inverse are dropped.
\subsection{Preconditioning nearly symmetric matrices}
\label{nearly}
A nonsymmetric, but definite and nearly symmetric matrix $A$
may be preconditioned
with a symmetric preconditioner $M$. Using a symmetric preconditioner
has a few advantages, such as guaranteeing positive
definiteness of the preconditioner, as well as being less expensive
to construct.
The nonsymmetric matrix $A$ must be definite,
i.e., $(A+A^T)/2$ is SPD, and the {\em a priori} sparsity pattern to be used
must be symmetric. The latter may be guaranteed by 1)
constructing the sparsity pattern with a symmetric matrix, or 2) if the
matrix is structurally symmetric (has symmetric pattern), then
thresholding to construct the pattern is not used (i.e.,
zero value of the {\tt thresh} parameter is used).
%==========================================================================