|
|
|
\section{Getting Started\label{sec:started}}
|
|
|
|
\markboth{\textsc{MLD2P4 User's and Reference Guide}}
|
|
|
|
{\textsc{\ref{sec:started} Getting Started}}
|
|
|
|
|
|
|
|
We describe the basics for building and applying MLD2P4 one-level and multi-level
|
|
|
|
Schwarz preconditioners with the Krylov solvers included in PSBLAS \cite{PSBLASGUIDE}.
|
|
|
|
The following steps are required:
|
|
|
|
\begin{enumerate}
|
|
|
|
\item \emph{Declare the preconditioner data structure}. It is a derived data type,
|
|
|
|
\verb|mld_|\-\emph{x}\verb|prec_| \verb|type|, where \emph{x} may be \verb|s|, \verb|d|, \verb|c|
|
|
|
|
or \verb|z|, according to the basic data type of the sparse matrix
|
|
|
|
(\verb|s| = real single precision; \verb|d| = real double precision;
|
|
|
|
\verb|c| = complex single precision; \verb|z| = complex double precision).
|
|
|
|
This data structure is accessed by the user only through the MLD2P4 routines,
|
|
|
|
following an object-oriented approach.
|
|
|
|
\item \emph{Allocate and initialize the preconditioner data structure, according to
|
|
|
|
a preconditioner type chosen by the user}. This is performed by the routine
|
|
|
|
\verb|mld_precinit|, which also sets defaults for each preconditioner
|
|
|
|
type selected by the user. The defaults associated to each preconditioner
|
|
|
|
type are given in Table~\ref{tab:precinit}, where the strings used by
|
|
|
|
\verb|mld_precinit| to identify the preconditioner types are also given.
|
|
|
|
Note that these strings are valid also if uppercase letters are substituted by
|
|
|
|
corresponding lowercase ones.
|
|
|
|
\item \emph{Modify the selected preconditioner type, by properly setting
|
|
|
|
preconditioner parameters.} This is performed by the routine \verb|mld_precset|.
|
|
|
|
This routine must be called only if the user wants to modify the default values
|
|
|
|
of the parameters associated to the selected preconditioner type, to obtain a variant
|
|
|
|
of the preconditioner. Examples of use of \verb|mld_precset| are given in
|
|
|
|
Section~\ref{sec:examples}; a complete list of all the
|
|
|
|
preconditioner parameters and their allowed and default values is provided in
|
|
|
|
Section~\ref{sec:userinterface}, Tables~\ref{tab:p_type}-\ref{tab:p_coarse}.
|
|
|
|
\item \emph{Build the preconditioner for a given matrix.} This is performed by
|
|
|
|
the routine \verb|mld_precbld|.
|
|
|
|
\item \emph{Apply the preconditioner at each iteration of a Krylov solver.}
|
|
|
|
This is performed by the routine \verb|mld_precaply|. When using the PSBLAS Krylov solvers,
|
|
|
|
this step is completely transparent to the user, since \verb|mld_precaply| is called
|
|
|
|
by the PSBLAS routine implementing the Krylov solver (\verb|psb_krylov|).
|
|
|
|
\item \emph{Free the preconditioner data structure}. This is performed by
|
|
|
|
the routine \verb|mld_| \verb|precfree|. This step is complementary to step 1 and should
|
|
|
|
be performed when the preconditioner is no more used.
|
|
|
|
\end{enumerate}
|
|
|
|
A detailed description of the above routines is given in Section~\ref{sec:userinterface}.
|
|
|
|
Examples showing the basic use of MLD2P4 are reported in Section~\ref{sec:examples}.
|
|
|
|
|
|
|
|
Note that the Fortran 95 module \verb|mld_prec_mod|, containing the definition of the
|
|
|
|
preconditioner data type and the interfaces to the routines of MLD2P4,
|
|
|
|
must be used in any program calling such routines.
|
|
|
|
The modules \verb|psb_base_mod|, for the sparse matrix and communication descriptor
|
|
|
|
data types, and \verb|psb_krylov_mod|, for interfacing with the
|
|
|
|
Krylov solvers, must be also used (see Section~\ref{sec:examples}).
|
|
|
|
|
|
|
|
\ \\
|
|
|
|
\textbf{Remark 1.} The coarsest-level solver used by the default two-level
|
|
|
|
preconditioner has been chosen by taking into account that, on parallel
|
|
|
|
machines, it often leads to the smallest execution time when applied to
|
|
|
|
linear systems coming from finite-difference discretizations of basic
|
|
|
|
elliptic PDE problems, considered as standard tests for multi-level Schwarz
|
|
|
|
preconditioners \cite{aaecc_07,apnum_07}. However, this solver does
|
|
|
|
not necessarily correspond to the smallest number of iterations of the
|
|
|
|
preconditioned Krylov method, which is usually obtained by applying
|
|
|
|
a direct solver to the coarsest-level system, e.g.\ based on the LU
|
|
|
|
factorization (see Section~\ref{sec:userinterface}
|
|
|
|
for the coarsest-level solvers available in MLD2P4).
|
|
|
|
|
|
|
|
\ \\
|
|
|
|
\textbf{Remark 2.} The include path for MLD2P4 must override
|
|
|
|
those for PSBLAS, e.g.\ the latter must come first in the sequence
|
|
|
|
passed to the compiler, as the MLD2P4 version of the Krylov solver
|
|
|
|
interfaces must override that of PSBLAS. This will change in the future
|
|
|
|
when the support for the \verb|class| statement becomes widespread in Fortran
|
|
|
|
compilers.
|
|
|
|
|
|
|
|
|
|
|
|
\begin{table}[th]
|
|
|
|
\begin{center}
|
|
|
|
%{\small
|
|
|
|
\begin{tabular}{|l|l|p{7.8cm}|}
|
|
|
|
\hline
|
|
|
|
\textsc{type} & \textsc{string} & \textsc{default preconditioner} \\ \hline
|
|
|
|
No preconditioner &\verb|'NOPREC'|& Considered only to use the PSBLAS
|
|
|
|
Krylov solvers with no preconditioner. \\ \hline
|
|
|
|
Diagonal & \verb|'DIAG'| & --- \\ \hline
|
|
|
|
Block Jacobi & \verb|'BJAC'| & Block Jacobi with ILU(0) on the local blocks.\\ \hline
|
|
|
|
Additive Schwarz & \verb|'AS'| & Restricted Additive Schwarz (RAS),
|
|
|
|
with overlap 1 and ILU(0) on the local blocks. \\ \hline
|
|
|
|
Multilevel &\verb|'ML'| & Multi-level hybrid preconditioner (additive on the
|
|
|
|
same level and multiplicative through the levels),
|
|
|
|
with post-smoothing only.
|
|
|
|
Number of levels: 2.
|
|
|
|
Post-smoother: RAS with overlap 1 and ILU(0)
|
|
|
|
on the local blocks.
|
|
|
|
Aggregation: decoupled smoothed aggregation with
|
|
|
|
threshold $\theta = 0$.
|
|
|
|
Coarsest matrix: distributed among the processors.
|
|
|
|
Coarsest-level solver:
|
|
|
|
4 sweeps of the block-Jacobi solver,
|
|
|
|
with LU (or ILU) factorization of the blocks
|
|
|
|
(UMFPACK for the double precision versions and
|
|
|
|
SuperLU for the single precision ones, if the packages
|
|
|
|
have been installed; ILU(0), otherwise). \\
|
|
|
|
\hline
|
|
|
|
\end{tabular}
|
|
|
|
%}
|
|
|
|
\end{center}
|
|
|
|
|
|
|
|
\caption{Preconditioner types, corresponding strings and default choices.
|
|
|
|
\label{tab:precinit}}
|
|
|
|
\end{table}
|
|
|
|
|
|
|
|
\subsection{Examples\label{sec:examples}}
|
|
|
|
|
|
|
|
The code reported in Figure~\ref{fig:ex_default} shows how to set and apply the default
|
|
|
|
multi-level preconditioner available in the real double precision version
|
|
|
|
of MLD2P4 (see Table~\ref{tab:precinit}). This preconditioner is chosen
|
|
|
|
by simply specifying \verb|'ML'| as second argument of \verb|mld_precinit|
|
|
|
|
(a call to \verb|mld_precset| is not needed) and is applied with the BiCGSTAB
|
|
|
|
solver provided by PSBLAS. As previously observed, the modules \verb|psb_base_mod|,
|
|
|
|
\verb|mld_prec_mod| and \verb|psb_krylov_mod| must be used by the example program.
|
|
|
|
|
|
|
|
The part of the code concerning the
|
|
|
|
reading and assembling of the sparse matrix and the right-hand side vector, performed
|
|
|
|
through the PSBLAS routines for sparse matrix and vector management, is not reported
|
|
|
|
here for brevity; the statements concerning the deallocation of the PSBLAS
|
|
|
|
data structure are neglected too.
|
|
|
|
The complete code can be found in the example program file \verb|mld_dexample_ml.f90|,
|
|
|
|
in the directory \verb|examples/fileread| of the MLD2P4 tree (see
|
|
|
|
Section~\ref{sec:ex_and_test}).
|
|
|
|
For details on the use of the PSBLAS routines, see the PSBLAS User's
|
|
|
|
Guide \cite{PSBLASGUIDE}.
|
|
|
|
|
|
|
|
The setup and application of the default multi-level
|
|
|
|
preconditioners for the real single precision and the complex, single and double
|
|
|
|
precision, versions are obtained with straightforward modifications of the previous
|
|
|
|
example (see Section~\ref{sec:userinterface} for details). If these versions are installed,
|
|
|
|
the corresponding Fortran 95 codes are available in \verb|examples/fileread/|.
|
|
|
|
|
|
|
|
\begin{figure}[tbp]
|
|
|
|
\begin{center}
|
|
|
|
\begin{minipage}{.90\textwidth}
|
|
|
|
{\small
|
|
|
|
\begin{verbatim}
|
|
|
|
use psb_base_mod
|
|
|
|
use mld_prec_mod
|
|
|
|
use psb_krylov_mod
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! sparse matrix
|
|
|
|
type(psb_dspmat_type) :: A
|
|
|
|
! sparse matrix descriptor
|
|
|
|
type(psb_desc_type) :: desc_A
|
|
|
|
! preconditioner
|
|
|
|
type(mld_dprec_type) :: P
|
|
|
|
! right-hand side and solution vectors
|
|
|
|
real(kind(1.d0)) :: b(:), x(:)
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! initialize the parallel environment
|
|
|
|
call psb_init(ictxt)
|
|
|
|
call psb_info(ictxt,iam,np)
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! read and assemble the matrix A and the right-hand side b
|
|
|
|
! using PSBLAS routines for sparse matrix / vector management
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! initialize the default multi-level preconditioner, i.e. hybrid
|
|
|
|
! Schwarz, using RAS (with overlap 1 and ILU(0) on the blocks)
|
|
|
|
! as post-smoother and 4 block-Jacobi sweeps (with UMFPACK LU
|
|
|
|
! on the blocks) as distributed coarse-level solver
|
|
|
|
call mld_precinit(P,'ML',info)
|
|
|
|
!
|
|
|
|
! build the preconditioner
|
|
|
|
call mld_precbld(A,desc_A,P,info)
|
|
|
|
!
|
|
|
|
! set the solver parameters and the initial guess
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! solve Ax=b with preconditioned BiCGSTAB
|
|
|
|
call psb_krylov('BICGSTAB',A,P,b,x,tol,desc_A,info)
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! deallocate the preconditioner
|
|
|
|
call mld_precfree(P,info)
|
|
|
|
!
|
|
|
|
! deallocate other data structures
|
|
|
|
... ...
|
|
|
|
!
|
|
|
|
! exit the parallel environment
|
|
|
|
call psb_exit(ictxt)
|
|
|
|
stop
|
|
|
|
\end{verbatim}
|
|
|
|
}
|
|
|
|
\end{minipage}
|
|
|
|
\caption{Setup and application of the default multi-level Schwarz preconditioner.
|
|
|
|
\label{fig:ex_default}}
|
|
|
|
\end{center}
|
|
|
|
\end{figure}
|
|
|
|
|
|
|
|
Different versions of multi-level preconditioners can be obtained by changing
|
|
|
|
the default values of the preconditioner parameters. The code reported in
|
|
|
|
Figure~\ref{fig:ex_3lh} shows how to set a three-level hybrid Schwarz
|
|
|
|
preconditioner, which uses block Jacobi with ILU(0) on the
|
|
|
|
local blocks as post-smoother, has a coarsest matrix replicated on the processors,
|
|
|
|
and solves the coarsest-level system with the LU factorization from UMFPACK~\cite{UMFPACK}.
|
|
|
|
The number of levels is specified by using \verb|mld_precinit|; the other
|
|
|
|
preconditioner parameters are set by calling \verb|mld_precset|. Note that
|
|
|
|
the type of multilevel framework (i.e.\ multiplicative among the levels
|
|
|
|
with post-smoothing only) is not specified since it is the default
|
|
|
|
set by \verb|mld_precinit|.
|
|
|
|
|
|
|
|
Figure~\ref{fig:ex_3la} shows how to
|
|
|
|
set a three-level additive Schwarz preconditioner,
|
|
|
|
which uses RAS, with overlap 1 and ILU(0) on the blocks,
|
|
|
|
as pre- and post-smoother, and applies five block-Jacobi sweeps, with
|
|
|
|
the UMFPACK LU factorization on the blocks, as distributed coarsest-level
|
|
|
|
solver. Again, \verb|mld_precset| is used only to set
|
|
|
|
non-default values of the parameters (see Tables~\ref{tab:p_type}-\ref{tab:p_coarse}).
|
|
|
|
In both cases, the construction and the application of the preconditioner
|
|
|
|
are carried out as for the default multi-level preconditioner.
|
|
|
|
The code fragments shown in in Figures~\ref{fig:ex_3lh}-\ref{fig:ex_3la} are
|
|
|
|
included in the example program file \verb|mld_dexample_ml.f90| too.
|
|
|
|
|
|
|
|
Finally, Figure~\ref{fig:ex_1l} shows the setup of a one-level
|
|
|
|
additive Schwarz preconditioner, i.e.\ RAS with overlap 2. The corresponding
|
|
|
|
example program is available in \verb|mld_dexample_| \verb|1lev.f90|.
|
|
|
|
|
|
|
|
For all the previous preconditioners, example programs where the sparse matrix and
|
|
|
|
the right-hand side are generated by discretizing a PDE with Dirichlet
|
|
|
|
boundary conditions are also available in the directory \verb|examples/pdegen|.
|
|
|
|
|
|
|
|
\ \\
|
|
|
|
\textbf{Remark 3.} Any PSBLAS-based program using the basic preconditioners
|
|
|
|
implemented in PSBLAS 2.0, i.e.\ the diagonal and block-Jacobi ones,
|
|
|
|
can use the diagonal and block-Jacobi preconditioners
|
|
|
|
implemented in MLD2P4 without any change in the code.
|
|
|
|
The PSBLAS-based program must be only recompiled
|
|
|
|
and linked to the MLD2P4 library.
|
|
|
|
\\
|
|
|
|
|
|
|
|
|
|
|
|
\begin{figure}[tbh]
|
|
|
|
\begin{center}
|
|
|
|
\begin{minipage}{.90\textwidth}
|
|
|
|
{\small
|
|
|
|
\begin{verbatim}
|
|
|
|
... ...
|
|
|
|
! set a three-level hybrid Schwarz preconditioner, which uses
|
|
|
|
! block Jacobi (with ILU(0) on the blocks) as post-smoother,
|
|
|
|
! a coarsest matrix replicated on the processors, and the
|
|
|
|
! LU factorization from UMFPACK as coarse-level solver
|
|
|
|
call mld_precinit(P,'ML',info,nlev=3)
|
|
|
|
call_mld_precset(P,mld_smoother_type_,'BJAC',info)
|
|
|
|
call mld_precset(P,mld_coarse_mat_,'REPL',info)
|
|
|
|
call mld_precset(P,mld_coarse_solve_,'UMF',info)
|
|
|
|
... ...
|
|
|
|
\end{verbatim}
|
|
|
|
}
|
|
|
|
\end{minipage}
|
|
|
|
|
|
|
|
\caption{Setup of a hybrid three-level Schwarz preconditioner.\label{fig:ex_3lh}}
|
|
|
|
\end{center}
|
|
|
|
\end{figure}
|
|
|
|
|
|
|
|
\begin{figure}[tbh]
|
|
|
|
\begin{center}
|
|
|
|
\begin{minipage}{.90\textwidth}
|
|
|
|
{\small
|
|
|
|
\begin{verbatim}
|
|
|
|
... ...
|
|
|
|
! set a three-level additive Schwarz preconditioner, which uses
|
|
|
|
! RAS (with overlap 1 and ILU(0) on the blocks) as pre- and
|
|
|
|
! post-smoother, and 5 block-Jacobi sweeps (with UMFPACK LU
|
|
|
|
! on the blocks) as distributed coarsest-level solver
|
|
|
|
call mld_precinit(P,'ML',info,nlev=3)
|
|
|
|
call mld_precset(P,mld_ml_type_,'ADD',info)
|
|
|
|
call_mld_precset(P,mld_smoother_pos_,'TWOSIDE',info)
|
|
|
|
call mld_precset(P,mld_coarse_sweeps_,5,info)
|
|
|
|
... ...
|
|
|
|
\end{verbatim}
|
|
|
|
}
|
|
|
|
\end{minipage}
|
|
|
|
|
|
|
|
\caption{Setup of an additive three-level Schwarz preconditioner.\label{fig:ex_3la}}
|
|
|
|
\end{center}
|
|
|
|
\end{figure}
|
|
|
|
|
|
|
|
\begin{figure}[tbh]
|
|
|
|
\begin{center}
|
|
|
|
\begin{minipage}{.90\textwidth}
|
|
|
|
{\small
|
|
|
|
\begin{verbatim}
|
|
|
|
... ...
|
|
|
|
! set RAS with overlap 2 and ILU(0) on the local blocks
|
|
|
|
call mld_precinit(P,'AS',info)
|
|
|
|
call mld_precset(P,mld_sub_ovr_,2,info)
|
|
|
|
... ...
|
|
|
|
\end{verbatim}
|
|
|
|
}
|
|
|
|
\end{minipage}
|
|
|
|
\caption{Setup of a one-level Schwarz preconditioner.\label{fig:ex_1l}}
|
|
|
|
\end{center}
|
|
|
|
\end{figure}
|
|
|
|
|
|
|
|
|
|
|
|
%%% Local Variables:
|
|
|
|
%%% mode: latex
|
|
|
|
%%% TeX-master: "userguide"
|
|
|
|
%%% End:
|