\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_| \verb|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|psb_|\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 \verb|psb_|\emph{x}\verb|vect_type| with \emph{x} = \verb|s|, \verb|d|, \verb|c|, \verb|z|, in a manner completely analogous to the sparse matrix type; \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. \clearpage \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. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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. This optional argument is deprecated, new codes should set the number of levels with \verb|mld_precset|.\\ % 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} \clearpage \subsection{Subroutine mld\_precset\label{sec:precset}} \begin{center} \verb|call 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|. The routine may also be invoked as a method of the preconditioner object as in the following: \begin{center} \verb|call p%set(what,val,info [,ilev, pos])|\\ \end{center} In this case it is also possible to specify an optional \verb|ilev| argument that restricts the effect of the call to the specified level. Finally, if the user has developed a new type of smoother and/or solver by extending one of the base MLD2P4 types, and has declared a variable of the new type in the main program, it is possible to pass the new smoother/solver variable to the setup routine as follows: \begin{center} \verb|call p%set(smoother,info [,ilev ,pos])|\\ \verb|call p%set(solver,info [,ilev ,pos])| \end{center} In this way, the variable will act as a \emph{mold} to which the preconditioner will conform, even though the MLD2P4 library is not modified, and thus has no direct knowledge about the new type. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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)| \emph{or} \verb|character(len=*)|. \\ & The parameter to be set. It can be specified by a predefined constant, or through its name; the string is case-insensitive. See also 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|smoother| & \verb|class(mld_x_base_smoother_type)| \\ & The user-defined new smoother to be employed in the preconditioner.\\ \verb|solver| & \verb|class(mld_x_base_solver_type)| \\ & The user-defined new solver to be employed in the preconditioner.\\ \verb|info| & \verb|integer, intent(out)|.\\ & Error code. If no error, 0 is returned. See Section~\ref{sec:errors} for details.\\ \verb|pos| & \verb|charater(len=*), intent(in)|.\\ & Whether the other arguments apply to the \verb|'PRE'| or to the \verb|'POST'| smoothers.\\ % %\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}. % The smoother and solver objects are arranged in a hierarchical manner; when specifying a smoother object, its parameters including the contained solver are set to default values, and when a solver object is specified its defaults are also set, overriding in both cases any previous settings even if explicitly specified. Therefore if the user sets a new smoother, and wishes to use a solver different from the default one, the call to set the solver must come \emph{after} the call to set the smoother. % The combination of a Jacobi smoother with a Diagonal Scaling local solver is equivalent to the strategy called Point Jacobi in the literature; similarly, having a Jacobi smoother with a Gauss-Seidel local solver is equivalent to a ``hybrid Gauss-Seidel'' solver. Completely new smoother and/or solver class derived from the base objects in the library may be used without recompiling the library itself. Once the new smoother/solver class has been developed, the user can declare a variable of that new type in the application, and pass that variable to the \verb|p%set(solver,info)| call; the new solver object is then dynamically included in the preconditioner structure. The \verb|what,val| pairs described here are those of the predefined smoother/solver objects; newly developed solvers may define new pairs according to their needs. \bsideways \begin{center} \begin{tabular}{|p{5cm}|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_| \break \verb|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_| \break \verb|SMOOTHER_TYPE| & \verb|character(len=*)| & \texttt{'JACOBI'} \ \ \ \texttt{'BJAC'} \ \ \ \texttt{'AS'} & \texttt{'AS'} & Basic predefined one-level preconditioner (i.e.\ smoother): Jacobi, block Jacobi, AS. \\ \hline \verb|mld_smoother_pos_| \break \verb|SMOOTHER_POS| & \verb|character(len=*)| & \texttt{'PRE'} \ \ \ \texttt{'POST'} \ \ \ \texttt{'TWOSIDE'} & \texttt{'TWOSIDE'} & ``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}} \esideways \bsideways \begin{center} \small \begin{tabular}{|p{3.5cm}|l|p{3.2cm}|l|p{5cm}|} \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_| \break \verb|SUB_OVR| & \verb|integer| & any~int.~num.~$\ge 0$ & 1 & Number of overlap layers. \\ \hline \verb|mld_sub_restr_| \break \verb|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_| \break \verb|SUB_PROL| & \verb|character(len=*)| & \texttt{'SUM'} \hspace{2.5cm} \texttt{'NONE'} & \texttt{'NONE'} & Type of prolongation operator: \texttt{'SUM'} for adding the contributions from the overlap, \texttt{'NONE'} for neglecting them. \\ \hline \verb|mld_sub_solve_| \break \verb|SUB_SOLVE| & \verb|character(len=*)| & \texttt{'DIAG'} \hspace{2.5cm} \texttt{'GS'} \hspace{2.5cm} \texttt{'BWGS'} \hspace{2.5cm} \texttt{'ILU'} \hspace{2.5cm} \texttt{'MILU'} \hspace{2.5cm} \texttt{'ILUT'} \hspace{2.5cm} \texttt{'UMF'} \hspace{2.5cm} \texttt{'SLU'} \hspace{2.5cm} \texttt{'MUMPS'} & \texttt{'ILU'} & Predefined local solver: pointwise Jacobi (diagonal scaling), (forward) Gauss-Seidel, Backward Gauss-Seidel, ILU($p$), MILU($p$), ILU($p,t$), LU from UMFPACK, LU from SuperLU (plus triangular solve), LU from MUMPS. \\ \hline \verb|mld_sub_fillin_| \break \verb|SUB_FILLIN| & \verb|integer| & Any~int.~num.~$\ge 0$ & 0 & Fill-in level $p$ of the incomplete LU factorizations. \\ \hline \verb|mld_sub_iluthrs_| \break \verb|SUB_ILUTHRS| & \verb|real(|\emph{kind\_parameter}\verb|)| & Any~real~num.~$\ge 0$ & 0 & Drop tolerance $t$ in the ILU($p,t$) factorization. \\ \hline \verb|mld_sub_ren_| \break \verb|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. \\ \verb|mld_solver_sweeps_| \break \verb|SOLVER_SWEEPS| & \verb|integer| & Any~int.~num.~$\ge 1$ & 1 & Number of sweeps for iterative local solver (currently only Gauss-Seidel). \\ \hline % \verb|mld_solver_eps_| \break \verb|SOLVER_EPS| & \verb|real| % & Any~real~num. % & 0 % & Stopping tolerance for iterative local solver % (currently only Gauss-Seidel); if $\le0$, then % perform prespecified number of iterations. \\ \hline \hline \end{tabular} \end{center} \caption{Parameters defining the one-level preconditioner used as smoother. \label{tab:p_smoother}} \esideways \bsideways \begin{center} \begin{tabular}{|p{5cm}|l|p{2.4cm}|p{2.4cm}|p{6cm}|} \hline \verb|what| & \textsc{data type} & \verb|val| & \textsc{default} & \textsc{comments} \\ \hline %\multicolumn{5}{|c|}{\emph{aggregation algorithm}} \\ \hline \verb|mld_coarse_aggr_size_| \break \verb|COARSE_AGGR_SIZE| & \verb|integer| & A positive number & The cubic root of the matrix size at the fine level. & Coarse size threshold. Aggregation will proceed until either the global number of variables is below this threshold, or the aggregation does not reduce the size any longer, or the maximum number of levels is reached. \\ \hline \verb|mld_min_aggr_ratio| \break \verb|MIN_AGGR_RATIO| & \verb|real| & A number greater than one & 1.5 & Minimum aggregation ratio. Aggregation will stop if the ratio between the matrix sizes at two consecutive levels drops below this threshold.\\ \hline \verb|mld_n_prec_levs_| \break \verb|N_PREC_LEVS| & \verb|integer| & A number greater than 1, or -1. & -1 & Number of levels; if set to a positive number greater than 1, it will force the aggregation algorithm to use this many levels (unless there is no reduction in the coarse matrix size). \\ \hline \verb|mld_max_prec_levs_| \break \verb|MAX_PREC_LEVS| & \verb|integer| & A positive number & 20 & Maximum number of levels: irrespective of the other settings, do not use more than this many levels. \\ \hline \verb|mld_aggr_alg_| \break \verb|AGGR_ALG| & \verb|character(len=*)| & \texttt{'DEC'}, \texttt{'SYMDEC'} & \texttt{'DEC'} & Aggregation algorithm. Currently, only the decoupled aggregation is available; the \verb|SYMDEC| option applies decoupled aggregation to the sparsity pattern of $A+A^T$.\\ \hline \verb|mld_aggr_ord_| \break \verb|AGGR_ORD| & \verb|character(len=*)| & \texttt{'NATURAL'} & \texttt{'DEGREE'} & Initial ordering of indices for aggregation algorithm: natural ordering or sorted by descending degree of the node in the matrix graph. Since aggregation is heuristics, results will be different. \\ \hline \hline \end{tabular} \end{center} \caption{Parameters defining the aggregation algorithm. \label{tab:p_aggregation}} \esideways \bsideways \begin{center} \begin{tabular}{|p{5cm}|l|p{2.4cm}|p{2.4cm}|p{6cm}|} \hline \verb|what| & \textsc{data type} & \verb|val| & \textsc{default} & \textsc{comments} \\ \hline \verb|mld_aggr_kind_| \break \verb|AGGR_KIND| & \verb|character(len=*)| & \texttt{'SMOOTHED'} \hspace{2.5cm} \texttt{'NONSMOOTHED'} & \texttt{'SMOOTHED'} & Type of aggregation: smoothed, nonsmoothed (i.e.\ using the tentative prolongator). \\ \hline \verb|mld_aggr_thresh_| \break \verb|AGGR_THRESH| & \verb|real(|\emph{kind\_parameter}\verb|)| & Any~real~num. $\in [0, 1]$ & 0.05 & Threshold $\theta$ in the aggregation algorithm. \\ \hline \verb|mld_aggr_scale_| \break \verb|AGGR_SCALE| & \verb|real(|\emph{kind\_parameter}\verb|)| & Any~real~num. $\in [0, 1]$ & 1.0 & Scale factor applied to the threshold going from level $ilev$ to level $ilev+1$. \\ \hline \verb|mld_aggr_omega_alg_| \break \verb|AGGR_OMEGA_ALG|& \verb|character(len=*)| & \texttt{'EIG\_EST'} \hspace{2.5cm} \texttt{'USER\_CHOICE'} & \texttt{'EIG\_EST'} & How the damping parameter $\omega$ in the smoothed aggregation should be computed: either via an estimate of the spectral radius of $D^{-1}A$, or explicily specified by the user. \\ \hline \verb|mld_aggr_eig_| \break \verb|AGGR_EIG| & \verb|character(len=*)| & \texttt{'A\_NORMI'} & \texttt{'A\_NORMI'} & How to estimate the spectral radius of $D^{-1}A$. Currently only the infinity norm estimate is available. \\ \hline \verb|mld_aggr_omega_val_| \break \verb|AGGR_OMEGA_VAL| & \verb|real(|\emph{kind\_parameter}\verb|)| & Any~nonnegative~real~num. & $4/(3\rho(D^{-1}A))$ & Damping parameter $\omega$ in the smoothed aggregation algorithm. It must be set by the user if \verb|USER_CHOICE| was specified for \verb|mld_aggr_omega_alg_|, otherwise it is computed by the library, using the selected estimate of the spectral radius $\rho(D^{-1}A)$ of $D^{-1}A$.\\ \hline \end{tabular} \end{center} \caption{Parameters defining the aggregation algorithm. \label{tab:p_aggregation_1}} \esideways \bsideways \begin{center} \begin{tabular}{|p{3.5cm}|l|p{3.2cm}|l|p{5cm}|} \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_| \break \verb|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_| \break \verb|COARSE_SOLVE| & \verb|character(len=*)| & \texttt{'BJAC'} \hspace{2.5cm} \texttt{'UMF'} \hspace{2.5cm} \texttt{'MUMPS'} \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 or MUMPS. \texttt{'SLUDIST'} and \texttt{'MUMPS'} require the coarsest matrix to be distributed, while \texttt{'UMF'} and \texttt{'SLU'} require it to be replicated. \\ \hline \verb|mld_coarse_subsolve_| \break \verb|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'} \hspace{2.5cm} \texttt{'MUMPS'} & See note & 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_| \break \verb|COARSE_SWEEPS| & \verb|integer| & Any~int.~num.~$> 0$ & 4 & Number of Block-Jacobi sweeps when 'BJAC' is used as coarsest-level solver. \\ \hline \verb|mld_coarse_fillin_| \break \verb|COARSE_FILLIN| & \verb|integer| & Any~int.~num.~$\ge 0$ & 0 & Fill-in level $p$ of the incomplete LU factorizations. \\ \hline \verb|mld_coarse_iluthrs_| \break \verb|COARSE_ILUTHRS| & \verb|real(|\emph{kind\_parameter}\verb|)| & Any~real.~num.~$\ge 0$ & 0 & Drop tolerance $t$ in the ILU($p,t$) factorization. \\ \hline \multicolumn{5}{|l|}{{\bfseries Note:} defaults for {\texttt mld\_coarse\_subsolve\_} are chosen as }\\ \multicolumn{5}{|l|}{single precision version: 'MUMPS' if installed, 'SLU' if installed, 'ILU' otherwise}\\ \multicolumn{5}{|l|}{double precision version: 'MUMPS' if installed, 'UMF' if installed, else 'SLU' if installed, 'ILU' otherwise}\\ \hline \end{tabular} \end{center} \caption{Parameters defining the coarse-space correction at the coarsest level.\label{tab:p_coarse}} \esideways % \par\noindent{\large\bfseries Note}\par\noindent % The defaults for parameter \verb|mld_coarse_subsolve_| in Table~\ref{tab:p_coarse} are determined % as follows: % \begin{description} % \item[Single precision data:] 'SLU' if installed, 'ILU' otherwise; % \item[Double precision data:] 'UMF' if installed, else 'SLU' if % installed, 'ILU' otherwise; % \end{description} \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|. For multilevel preconditioner this routine is supported for backward compatibility, but we recommend to use the routines of Sec.~\ref{sec:hier_bld} and~\ref{sec:smoothers_bld}. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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} \subsection{Subroutine mld\_hierarchy\_bld\label{sec:hier_bld}} \begin{center} \verb|mld_hierachy_bld(a,desc_a,p,info)|\\ \end{center} \noindent This routine builds the aggregation hierarchy according to the requirements made by the user through the routines \verb|mld_precinit| and \verb|mld_precset|. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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} \subsection{Subroutine mld\_smoothers\_bld\label{sec:smoothers_bld}} \begin{center} \verb|mld_smoothers_bld(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|, based on the aggregation hierahy produced by a previous call to \verb|mld_hierarchy_bld| (see Sec.~\ref{sec:hier_bld}). {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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|psb_krylov| and hence it is completely transparent to the user. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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_| \verb|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} \clearpage \subsection{Subroutine mld\_precfree\label{sec:precfree}} \begin{center} \verb|mld_precfree(p,info)|\\ \end{center} \noindent This routine deallocates the preconditioner data structure. {\vskip2\baselineskip\noindent\large\bfseries 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} \clearpage \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. It must be called after \verb|mld_precbld| has been called. {\vskip2\baselineskip\noindent\large\bfseries Arguments} \begin{tabular}{p{1.2cm}p{12cm}} \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: