You cannot select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.
amg4psblas/docs/pdf/userinterface.tex

426 lines
21 KiB
TeX

\section{User Interface\label{sec:userinterface}}
\markboth{\textsc{MLD2P4 User's and Reference Guide}}
{\textsc{\ref{sec:userinterface} User Interface}}
The basic user interface of MLD2P4 consists of six routines. The four routines \verb|mld_precinit|,
\verb|mld_precset|, \verb|mld_precbld| and \verb|mld_precaply| encapsulate all the functionalities
for the setup and the application of any one-level and multi-level
preconditioner implemented in the package.
The routine \verb|mld_precfree| deallocates the preconditioner data structure, while
\verb|mld_precdescr| prints a description of the preconditioner setup by the user.
For each routine, the same user interface is overloaded with
respect to the real/complex case and the single/double precision;
arguments with appropriate data types must be passed to the routine,
i.e.
\begin{itemize}
\item the sparse matrix data structure, containing the matrix to be
preconditioned, must be of type \verb|mld_|\emph{x}\verb|spmat_type|
with \emph{x} = \verb|s| for real single precision, \emph{x} = \verb|d|
for real double precision, \emph{x} = \verb|c| for complex single precision,
\emph{x} = \verb|z| for complex double precision;
\item the preconditioner data structure must be of type
\verb|mld_|\emph{x}\verb|prec_type|, with \emph{x} =
\verb|s|, \verb|d|, \verb|c|, \verb|z|, according to the sparse
matrix data structure;
\item the arrays containing the vectors $v$ and $w$ involved in
the preconditioner application $w=M^{-1}v$ must be of type
\emph{type}\verb|(|\emph{kind\_parameter}\verb|)|, with \emph{type} =
\verb|real|, \verb|complex| and \emph{kind\_parameter} = \verb|kind(1.e0)|,
\verb|kind(1.d0)|, according to the sparse matrix and preconditioner
data structure; note that the PSBLAS module \verb|psb_base_mod|
provides the constants \verb|psb_spk_|
= \verb|kind(1.e0)| and \verb|psb_dpk_| = \verb|kind(1.d0)|;
\item real parameters defining the preconditioner must be declared
according to the precision of the sparse matrix and preconditioner
data structures (see Section~\ref{sec:precset}).
\end{itemize}
A description of each routine is given in the remainder of this section.
\subsection{Subroutine mld\_precinit\label{sec:precinit}}
\begin{center}
\verb|mld_precinit(p,ptype,info)| \\
\verb|mld_precinit(p,ptype,info,nlev)| \\
\end{center}
\noindent
This routine allocates and initializes the preconditioner data structure,
according to the preconditioner type chosen by the user.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.5cm}}
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(inout)|.\\
& The preconditioner data structure. Note that \emph{x}
must be chosen according to the real/complex, single/double
precision version of MLD2P4 under use.\\
\verb|ptype| & \verb|character(len=*), intent(in)|.\\
& The type of preconditioner. Its values are specified
in Table~\ref{tab:precinit}.\\
& Note that the strings are case insensitive.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\
\verb|nlev| & \verb|integer, optional, intent(in)|.\\
& The number of levels of the multilevel preconditioner.
If \verb|nlev| is not present and \verb|ptype|=\verb|'ML'|, \verb|'ml'|,
then \verb|nlev|=2 is assumed. Otherwise, \verb|nlev| is ignored.\\
\end{tabular}
\subsection{Subroutine mld\_precset\label{sec:precset}}
\begin{center}
\verb|mld_precset(p,what,val,info)|\\
\end{center}
\noindent
This routine sets the parameters defining the preconditioner. More
precisely, the parameter identified by \verb|what| is assigned the value
contained in \verb|val|.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.5cm}}
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(inout)|.\\
& The preconditioner data structure. Note that \emph{x} must
be chosen according to the real/complex, single/double precision
version of MLD2P4 under use.\\
\verb|what| & \verb|integer, intent(in)|. \\
& The number identifying the parameter to be set.
A mnemonic constant has been associated to each of these
numbers, as reported in Tables~\ref{tab:p_type}-\ref{tab:p_coarse}.\\
\verb|val | & \verb|integer| \emph{or} \verb|character(len=*)| \emph{or}
\verb|real(psb_spk_)| \emph{or} \verb|real(psb_dpk_)|,
\verb|intent(in)|.\\
& The value of the parameter to be set. The list of allowed
values and the corresponding data types is given in
Tables~\ref{tab:p_type}-\ref{tab:p_coarse}.
When the value is of type \verb|character(len=*)|,
it is also treated as case insensitive.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors}
for details.\\
%
%\verb|ilev| & \verb|integer, optional, intent(in)|.\\
% & For the multilevel preconditioner, the level at which the
% preconditioner parameter has to be set.
% The levels are numbered in increasing
% order starting from the finest one, i.e.\ level 1 is the finest level.
% If \verb|ilev| is not present, the parameter identified by \verb|what|
% is set at all the appropriate levels (see Table~\ref{tab:params}).
\end{tabular}
\ \\
A variety of (one-level and multi-level) preconditioners can be obtained
by a suitable setting of the preconditioner parameters. These parameters
can be logically divided into four groups, i.e.\ parameters defining
\begin{enumerate}
\item the type of multi-level preconditioner;
\item the one-level preconditioner used as smoother;
\item the aggregation algorithm;
\item the coarse-space correction at the coarsest level.
\end{enumerate}
A list of the parameters that can be set, along with their allowed and
default values, is given in Tables~\ref{tab:p_type}-\ref{tab:p_coarse}.
For a detailed description of the meaning of the parameters, please
refer to Section~\ref{sec:background}.
%
%Note that the routine allows to set different features of the
%preconditioner at each level through the use of \verb|ilev|.
%This should be done by users with experience in the field of
%multi-level preconditioners. Non-expert users are recommended
%to call \verb| mld_precset| without specifying \verb|ilev|.
\begin{sidewaystable}
\begin{center}
\begin{tabular}{|l|l|p{2cm}|l|p{7cm}|}
\hline
\verb|what| & \textsc{data type} & \verb|val| & \textsc{default} &
\textsc{comments} \\ \hline
%\multicolumn{5}{|c|}{\emph{type of the multi-level preconditioner}}\\ \hline
\verb|mld_ml_type_| & \verb|character(len=*)|
& \texttt{'ADD'} \ \ \ \texttt{'MULT'}
& \texttt{'MULT'}
& Basic multi-level framework: additive or multiplicative
among the levels (always additive inside a level). \\ \hline
\verb|mld_smoother_type_|& \verb|character(len=*)|
& \texttt{'DIAG'} \ \ \ \texttt{'BJAC'} \ \ \ \texttt{'AS'}
& \texttt{'AS'}
& Basic one-level preconditioner (i.e.\ smoother): diagonal,
block Jacobi, AS \\ \hline
\verb|mld_smoother_pos_| & \verb|character(len=*)|
& \texttt{'PRE'} \ \ \ \texttt{'POST'} \ \ \ \texttt{'TWOSIDE'}
& \texttt{'POST'}
& ``Position'' of the smoother: pre-smoother, post-smoother,
pre- and post-smoother. \\
\hline
\end{tabular}
\end{center}
\caption{Parameters defining the type of multi-level preconditioner.
\label{tab:p_type}}
\end{sidewaystable}
\begin{sidewaystable}
\begin{center}
\begin{tabular}{|l|l|p{3.2cm}|l|p{7cm}|}
\hline
\verb|what| & \textsc{data type} & \verb|val| & \textsc{default} &
\textsc{comments} \\ \hline
%\multicolumn{5}{|c|}{\emph{basic one-level preconditioner (smoother)}} \\ \hline
\verb|mld_sub_ovr_| & \verb|integer|
& any integer number $\ge 0$
& 1
& Number of overlap layers. \\ \hline
\verb|mld_sub_restr_| & \verb|character(len=*)|
& \texttt{'HALO'} \hspace{2.5cm} \texttt{'NONE'}
& \texttt{'HALO'}
& Type of restriction operator:
\texttt{'HALO'} for taking into account the overlap, \texttt{'NONE'}
for neglecting it. \\ \hline
\verb|mld_sub_prol_| & \verb|character(len=*)|
& \texttt{'SUM'} \hspace{2.5cm} \texttt{'NONE'}
& \texttt{'NONE'}
& Type of prolongator operator:
\texttt{'SUM'} for adding the contributions from the overlap, \texttt{'NONE'}
for neglecting them. \\ \hline
\verb|mld_sub_solve_| & \verb|character(len=*)|
& \texttt{'ILU'} \hspace{2.5cm} \texttt{'MILU'} \hspace{2.5cm} \texttt{'ILUT'}
\hspace{2.5cm} \texttt{'UMF'} \hspace{2.5cm} \texttt{'SLU'}
& \texttt{'UMF'}
& Local solver: ILU($p$), MILU($p$), ILU($p,t$), LU from UMFPACK, LU from SuperLU
(plus triangular solve). \\ \hline
\verb|mld_sub_fillin_| & \verb|integer|
& Any~int.~num.~$\ge 0$
& 0
& Fill-in level $p$ of the incomplete LU factorizations. \\ \hline
\verb|mld_sub_iluthrs_| & \verb|real(|\emph{kind\_parameter}\verb|)|
& Any~real~num.~$\ge 0$
& \texttt{0.e0} (or \texttt{0.d0})
& Drop tolerance $t$ in the ILU($p,t$) factorization. \\ \hline
\verb|mld_sub_ren_| & \verb|character(len=*)|
& \texttt{'RENUM\_NONE'} \texttt{'RENUM\_GLOBAL'} %, \texttt{'RENUM_GPS'}
& \texttt{'RENUM\_NONE'}
& Row and column reordering of the local submatrices: no reordering,
reordering according to the global numbering of the rows and columns of
the whole matrix. \\
\hline
\end{tabular}
\end{center}
\caption{Parameters defining the one-level preconditioner used as smoother.
\label{tab:p_smoother}}
\end{sidewaystable}
\begin{sidewaystable}
\begin{center}
\begin{tabular}{|l|l|p{2.3cm}|p{2.6cm}|p{7cm}|}
\hline
\verb|what| & \textsc{data type} & \verb|val| & \textsc{default} &
\textsc{comments} \\ \hline
%\multicolumn{5}{|c|}{\emph{aggregation algorithm}} \\ \hline
\verb|mld_aggr_alg_| & \verb|character(len=*)|
& \texttt{'DEC'}
& \texttt{'DEC'}
& Aggregation algorithm. Currently, only the decoupled aggregation is available. \\ \hline
\verb|mld_aggr_kind_| & \verb|character(len=*)|
& \texttt{'SMOOTH'} \hspace{2.5cm} \texttt{'RAW'}
& \texttt{'SMOOTH'}
& Type of aggregation: smoothed or raw, i.e.\ using the tentative prolongator. \\ \hline
\verb|mld_aggr_thresh_| & \verb|real(|\emph{kind\_parameter}\verb|)|
& Any~real~num. $\in [0, 1]$
& \texttt{0.e0} (or \texttt{0.d0})
& The threshold $\theta$ in the aggregation algorithm. \\ \hline
\verb|mld_aggr_eig_| & \verb|character(len=*)|
& \texttt{'A\_NORMI'}
& \texttt{'A\_NORMI'}
& Estimate of the maximum eigenvalue of $D^{-1}A$
for the smoothed aggregation. Currently, only the infinity norm of
the matrix is available. \\ \hline
\verb|mld_aggr_damp_| & \verb|real(|\emph{kind\_parameter}\verb|)|
& Any~real~num. $>0$
& \texttt{4.e0/3.e0} \hspace{1cm} (or \texttt{4.d0/3.d0})
& The damping parameter $\omega$ in the aggregation algorithm. \\
\hline
\end{tabular}
\end{center}
\caption{Parameters defining the aggregation algorithm.
\label{tab:p_aggregation}}
\end{sidewaystable}
\begin{sidewaystable}
\begin{center}
\begin{tabular}{|l|l|p{3.2cm}|l|p{7cm}|}
\hline
\verb|what| & \textsc{data type} & \verb|val| & \textsc{default} &
\textsc{comments} \\ \hline
%\multicolumn{5}{|c|}{\emph{coarse-space correction at the coarsest level}}\\ \hline
\verb|mld_coarse_mat_| & \verb|character(len=*)|
& \texttt{'DISTR'} \hspace{2.5cm} \texttt{'REPL'}
& \texttt{'DISTR'}
& Coarsest matrix: distributed among the processors or replicated on each of them. \\ \hline
\verb|mld_coarse_solve_| & \verb|character(len=*)|
& \texttt{'BJAC'} \hspace{2.5cm} \texttt{'UMF'} \hspace{2.5cm} \texttt{'SLU'}
\hspace{2.5cm} \texttt{'SLUDIST'}
& \texttt{'BJAC'}
& Solver used at the coarsest level: block Jacobi, sequential LU from UMFPACK,
sequential LU from SuperLU, distributed LU from SuperLU\_Dist.
If the coarsest matrix is distributed, only \texttt{'BJAC'} and \texttt{'SLUDIST'}
can be chosen; if it is replicated, only \texttt{'BJAC'} or \texttt{'SLUDIST'} can
be selected. \\ \hline
\verb|mld_coarse_subsolve_| & \verb|character(len=*)|
& \texttt{'ILU'} \hspace{2.5cm} \texttt{'MILU'} \hspace{2.5cm} \texttt{'ILUT'}
\hspace{2.5cm} \texttt{'UMF'} \hspace{2.5cm} \texttt{'SLU'}
& \texttt{'UMF'}
& Solver for the diagonal blocks of the coarse matrix, in case the block Jacobi solver
is chosen as coarsest-level solver: ILU($p$), MILU($p$), ILU($p,t$), LU from UMFPACK,
LU from SuperLU, plus triangular solve. \\ \hline
\verb|mld_coarse_sweeps_|& \verb|integer|
& Any~int.~num.~$> 0$
& \texttt{4}
& Number of Block-Jacobi sweeps when 'BJAC' is used as coarsest-level solver. \\
\verb|mld_coarse_fillin_| & \verb|integer|
& Any~int.~num.~$\ge 0$
& \texttt{0}
& Fill-in level $p$ of the incomplete LU factorizations. \\ \hline
\verb|mld_coarse_iluthrs_| & \verb|real(|\emph{kind\_parameter}\verb|)|
& Any~real.~num.~$\ge 0$
& \texttt{0.d0} (or \texttt{0.e0})
& Drop tolerance $t$ in the ILU($p,t$) factorization. \\
\hline
\end{tabular}
\end{center}
\caption{Parameters defining the coarse-space correction at the coarsest
level.\label{tab:p_coarse}}
\end{sidewaystable}
\clearpage
\subsection{Subroutine mld\_precbld\label{sec:precbld}}
\begin{center}
\verb|mld_precbld(a,desc_a,p,info)|\\
\end{center}
\noindent
This routine builds the preconditioner according to the requirements made by
the user through the routines \verb|mld_precinit| and \verb|mld_precset|.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.5cm}}
\verb|a| & \verb|type(psb_|\emph{x}\verb|spmat_type), intent(in)|. \\
& The sparse matrix structure containing the local part of the
matrix to be preconditioned. Note that \emph{x} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.
See the PSBLAS User's Guide for details \cite{PSBLASGUIDE}.\\
\verb|desc_a| & \verb|type(psb_desc_type), intent(in)|. \\
& The communication descriptor of \verb|a|. See the PSBLAS User's Guide for
details \cite{PSBLASGUIDE}.\\
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(inout)|.\\
& The preconditioner data structure. Note that \emph{x} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\
\end{tabular}
\clearpage
\subsection{Subroutine mld\_precaply\label{sec:precaply}}
\begin{center}
\verb|mld_precaply(p,x,y,desc_a,info)|\\
\verb|mld_precaply(p,x,y,desc_a,info,trans,work)|\\
\end{center}
\noindent
This routine computes $y = op(M^{-1})\, x$, where $M$ is a previously built
preconditioner, stored into \verb|p|, and $op$
denotes the preconditioner itself or its transpose, according to
the value of \verb|trans|.
Note that, when MLD2P4 is used with a Krylov solver from PSBLAS,
\verb|mld_precaply| is called within the PSBLAS routine \verb|mld_krylov|
and hence it is completely transparent to the user.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.5cm}}
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(inout)|.\\
& The preconditioner data structure, containing the local part of $M$.
Note that \emph{x} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|x| & \emph{type}\verb|(|\emph{kind\_parameter}\verb|), dimension(:), intent(in)|.\\
& The local part of the vector $x$. Note that \emph{type} and
\emph{kind\_parameter} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|y| & \emph{type}\verb|(|\emph{kind\_parameter}\verb|), dimension(:), intent(out)|.\\
& The local part of the vector $y$. Note that \emph{type} and
\emph{kind\_parameter} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|desc_a| & \verb|type(psb_desc_type), intent(in)|. \\
& The communication descriptor associated to the matrix to be
preconditioned.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\
\verb|trans| & \verb|character(len=1), optional, intent(in).|\\
& If \verb|trans| = \verb|'N','n'| then $op(M^{-1}) = M^{-1}$;
if \verb|trans| = \verb|'T','t'| then $op(M^{-1}) = M^{-T}$
(transpose of $M^{-1})$; if \verb|trans| = \verb|'C','c'| then $op(M^{-1}) = M^{-C}$
(conjugate transpose of $M^{-1})$.\\
\verb|work| & \emph{type}\verb|(|\emph{kind\_parameter}\verb|), dimension(:), optional, target|.\\
& Workspace. Its size should be at
least \verb|4 * psb_cd_get_local_cols(desc_a)| (see the PSBLAS User's Guide).
Note that \emph{type} and \emph{kind\_parameter} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\end{tabular}
\subsection{Subroutine mld\_precfree\label{sec:precfree}}
\begin{center}
\verb|mld_precfree(p,info)|\\
\end{center}
\noindent
This routine deallocates the preconditioner data structure.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.5cm}}
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(inout)|.\\
& The preconditioner data structure. Note that \emph{x} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\
\end{tabular}
\subsection{Subroutine mld\_precdescr\label{sec:precdescr}}
\begin{center}
\verb|mld_precdescr(p,info)|\\
\verb|mld_precdescr(p,info,iout)|\\
\end{center}
\noindent
This routine prints a description of the preconditioner
to the standard output or to a file.
\subsubsection*{Arguments}
\begin{tabular}{p{1.2cm}p{10.6cm}}
\verb|p| & \verb|type(mld_|\emph{x}\verb|prec_type), intent(in)|.\\
& The preconditioner data structure. Note that \emph{x} must be chosen according
to the real/complex, single/double precision version of MLD2P4 under use.\\
\verb|info| & \verb|integer, intent(out)|.\\
& Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\
\verb|iout| & \verb|integer, intent(in), optional|.\\
& The id of the file where the preconditioner description
will be printed; the default is the standard output.\\
\end{tabular}
%%% Local Variables:
%%% mode: latex
%%% TeX-master: "userguide"
%%% End: