<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 3.2 Final//EN">

<!--Converted with LaTeX2HTML 2008 (1.71)
original version by:  Nikos Drakos, CBLU, University of Leeds
* revised and updated by:  Marcus Hennecke, Ross Moore, Herb Swan
* with significant contributions from:
  Jens Lippmann, Marek Rouchal, Martin Wilck and others -->
<HTML>
<HEAD>
<TITLE>Subroutine mld_precset</TITLE>
<META NAME="description" CONTENT="Subroutine mld_precset">
<META NAME="keywords" CONTENT="userhtml">
<META NAME="resource-type" CONTENT="document">
<META NAME="distribution" CONTENT="global">

<META NAME="Generator" CONTENT="LaTeX2HTML v2008">
<META HTTP-EQUIV="Content-Style-Type" CONTENT="text/css">

<LINK REL="STYLESHEET" HREF="userhtml.css">

<LINK REL="next" HREF="node19.html">
<LINK REL="previous" HREF="node17.html">
<LINK REL="up" HREF="node16.html">
<LINK REL="next" HREF="node19.html">
</HEAD>

<BODY >
<!--Navigation Panel-->
<A NAME="tex2html284"
  HREF="node19.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html280"
  HREF="node16.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html274"
  HREF="node17.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html282"
  HREF="node2.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>  
<BR>
<B> Next:</B> <A NAME="tex2html285"
  HREF="node19.html">Subroutine mld_precbld</A>
<B> Up:</B> <A NAME="tex2html281"
  HREF="node16.html">User Interface</A>
<B> Previous:</B> <A NAME="tex2html275"
  HREF="node17.html">Subroutine mld_precinit</A>
 &nbsp; <B>  <A NAME="tex2html283"
  HREF="node2.html">Contents</A></B> 
<BR>
<BR>
<!--End of Navigation Panel-->

<H2><A NAME="SECTION00082000000000000000"></A><A NAME="sec:precset"></A>
<BR>
Subroutine mld_precset
</H2>

<P>
<DIV ALIGN="CENTER">
<code>mld_precset(p,what,val,info)</code>
<BR><code>mld_precset(p,smoother,info)</code>
<BR><code>mld_precset(p,solver,info)</code>
<BR>
</DIV>

<P>
This routine sets the parameters defining the preconditioner. More
precisely, the parameter identified by <code>what</code> is assigned the value
contained in <code>val</code>. The other two forms of this routine are
designed to allow extensions of the library by passing new smoothers
and  solvers to be employed in the preconditioner.

<P>
<FONT SIZE="+1"><B>Arguments</B></FONT>

<P>
<TABLE CELLPADDING=3>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>p</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>type(mld_</code><I>x</I><code>prec_type), intent(inout)</code>.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>The preconditioner data structure. Note that <I>x</I> must
                be chosen according to the real/complex, single/double precision
                 version of MLD2P4 under use.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>what</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>integer, intent(in)</code>.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>The number identifying the parameter to be set.
                A mnemonic constant has been associated to each of these
                numbers, as reported in Tables&nbsp;<A HREF="#tab:p_type">2</A>-<A HREF="#tab:p_coarse">5</A>.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>val </code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>integer</code> <I>or</I> <code>character(len=*)</code> <I>or</I>
                <code>real(psb_spk_)</code> <I>or</I> <code>real(psb_dpk_)</code>,
                <code>intent(in)</code>.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>The value of the parameter to be set. The list of allowed
                values and the corresponding data types is given in
                Tables&nbsp;<A HREF="#tab:p_type">2</A>-<A HREF="#tab:p_coarse">5</A>.
                When the value is of type <code>character(len=*)</code>,
                it is also treated as case insensitive.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>smoother</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>class(mld_x_base_smoother_type)</code></TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>The user-defined new smoother to be employed in the
                preconditioner.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>solver</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>class(mld_x_base_solver_type)</code></TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>The user-defined new solver to be employed in the
                preconditioner.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34><code>info</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340><code>integer, intent(out)</code>.</TD>
</TR>
<TR><TD ALIGN="LEFT" VALIGN="TOP" WIDTH=34>&nbsp;</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=340>Error code. If no error, 0 is returned. See Section&nbsp;<A HREF="node23.html#sec:errors">7</A>
                for details.</TD>
</TR>
</TABLE>

<P>
 
<BR>
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

<OL>
<LI>the type of multi-level preconditioner;
</LI>
<LI>the one-level preconditioner used as smoother;
</LI>
<LI>the aggregation algorithm;
</LI>
<LI>the coarse-space correction at the coarsest level.
</LI>
</OL>
A list of the parameters that can be set, along with their allowed and
default values, is given in Tables&nbsp;<A HREF="#tab:p_type">2</A>-<A HREF="#tab:p_coarse">5</A>.
For a detailed description  of the meaning of the parameters, please
refer to Section&nbsp;<A HREF="node11.html#sec:background">4</A>. 
The smoother and solver objects are arranged in a hierarchical manner;
when specifying a new smoother object, its parameters including the
contained solver are set to  default values, and when a new solver
object is specified its defaults are also set, overriding in both
cases any previous settings even if explicitly specified. Therefore if
the user specifies a new smoother, and whishes to use a new solver
which is not the default one, the call to set the solver must come
<I>after</I> the call to set the smoother. 

<P>
<BR><P></P>
<DIV ALIGN="CENTER"><A NAME="1258"></A>
<TABLE>
<CAPTION><STRONG>Table 2:</STRONG>
Parameters defining the type of multi-level preconditioner.
</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
<TABLE CELLPADDING=3 BORDER="1" ALIGN="CENTER">
<TR><TD ALIGN="LEFT"><code>what</code></TD>
<TD ALIGN="LEFT"><SMALL>DATA TYPE</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=57><code>val</code></TD>
<TD ALIGN="LEFT"><SMALL>DEFAULT</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198><SMALL>COMMENTS</SMALL></TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_ml_type_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=57><TT>'ADD'</TT>    <TT>'MULT'</TT></TD>
<TD ALIGN="LEFT"><TT>'MULT'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Basic multi-level framework: additive or multiplicative
                           among the levels (always additive inside a level).</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_smoother_type_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=57><TT>'DIAG'</TT>    <TT>'BJAC'</TT>    <TT>'AS'</TT></TD>
<TD ALIGN="LEFT"><TT>'AS'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Basic predefined one-level preconditioner (i.e. smoother): diagonal,
                           block Jacobi, AS.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_smoother_pos_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=57><TT>'PRE'</TT>    <TT>'POST'</TT>    <TT>'TWOSIDE'</TT></TD>
<TD ALIGN="LEFT"><TT>'POST'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>``Position'' of the smoother: pre-smoother, post-smoother, 
                           pre- and post-smoother.</TD>
</TR>
</TABLE>
</DIV>
</TD></TR>
</TABLE>
</DIV><P></P>
<BR>

<P>
<BR><P></P>
<DIV ALIGN="CENTER"><A NAME="1260"></A>
<TABLE>
<CAPTION><STRONG>Table 3:</STRONG>
Parameters defining the one-level preconditioner used as smoother.
</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
<TABLE CELLPADDING=3 BORDER="1" ALIGN="CENTER">
<TR><TD ALIGN="LEFT"><code>what</code></TD>
<TD ALIGN="LEFT"><SMALL>DATA TYPE</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><code>val</code></TD>
<TD ALIGN="LEFT"><SMALL>DEFAULT</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198><SMALL>COMMENTS</SMALL></TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_ovr_</code></TD>
<TD ALIGN="LEFT"><code>integer</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>any&nbsp;int.&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img90.png"
 ALT="$\ge 0$"></TD>
<TD ALIGN="LEFT">1</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Number of overlap layers.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_restr_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'HALO'</TT>   <TT>'NONE'</TT></TD>
<TD ALIGN="LEFT"><TT>'HALO'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Type of restriction operator:
                           <TT>'HALO'</TT> for taking into account the overlap, <TT>'NONE'</TT> 
                           for neglecting it.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_prol_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'SUM'</TT>   <TT>'NONE'</TT></TD>
<TD ALIGN="LEFT"><TT>'NONE'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Type of prolongation operator:
                           <TT>'SUM'</TT> for adding the contributions from the overlap, <TT>'NONE'</TT>
                           for neglecting them.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_solve_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'ILU'</TT>   <TT>'MILU'</TT>   <TT>'ILUT'</TT> 
                             <TT>'UMF'</TT>   <TT>'SLU'</TT></TD>
<TD ALIGN="LEFT"><TT>'ILU'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Predefined local solver: ILU(<IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$">),
                         MILU(<IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$">), ILU(<IMG
 WIDTH="27" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img36.png"
 ALT="$p,t$">), LU from UMFPACK, LU  from SuperLU 
                           (plus triangular solve).</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_fillin_</code></TD>
<TD ALIGN="LEFT"><code>integer</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>Any&nbsp;int.&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img90.png"
 ALT="$\ge 0$"></TD>
<TD ALIGN="LEFT">0</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Fill-in level <IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$"> of the incomplete LU factorizations.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_iluthrs_</code></TD>
<TD ALIGN="LEFT"><code>real(</code><I>kind_parameter</I><code>)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>Any&nbsp;real&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img90.png"
 ALT="$\ge 0$"></TD>
<TD ALIGN="LEFT">0</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Drop tolerance <IMG
 WIDTH="11" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
 SRC="img91.png"
 ALT="$t$"> in the ILU(<IMG
 WIDTH="27" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img36.png"
 ALT="$p,t$">) factorization.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_sub_ren_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'RENUM_NONE'</TT>  <TT>'RENUM_GLOBAL'</TT> </TD>
<TD ALIGN="LEFT"><TT>'RENUM_NONE'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>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.</TD>
</TR>
</TABLE>
</DIV>
</TD></TR>
</TABLE>
</DIV><P></P>
<BR>

<P>
<BR><P></P>
<DIV ALIGN="CENTER"><A NAME="1262"></A>
<TABLE>
<CAPTION><STRONG>Table 4:</STRONG>
Parameters defining the aggregation algorithm.
</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
<TABLE CELLPADDING=3 BORDER="1" ALIGN="CENTER">
<TR><TD ALIGN="LEFT"><code>what</code></TD>
<TD ALIGN="LEFT"><SMALL>DATA TYPE</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><code>val</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><SMALL>DEFAULT</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198><SMALL>COMMENTS</SMALL></TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_aggr_size_</code></TD>
<TD ALIGN="LEFT"><code>integer</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68>A positive number</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68>0, meaning that the size is fixed at
                         <code>precinit</code> time</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Coarse size threshold. Disregard the
                         original specification of number of levels in
                         <code>precinit</code> and continue aggregation
                         until either the global number of variables
                         is below this threshold, or the aggregation
                         does not reduce the size any longer.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_alg_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'DEC'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'DEC'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Aggregation algorithm. Currently, only the
                         decoupled aggregation is available.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_kind_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'SMOOTHED'</TT>   <TT>'NONSMOOTHED'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'SMOOTHED'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Type of aggregation: smoothed, nonsmoothed
                         (i.e. using the tentative prolongator).</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_thresh_</code></TD>
<TD ALIGN="LEFT"><code>real(</code><I>kind_parameter</I><code>)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68>Any&nbsp;real&nbsp;num. <IMG
 WIDTH="56" HEIGHT="36" ALIGN="MIDDLE" BORDER="0"
 SRC="img92.png"
 ALT="$\in [0, 1]$"></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68>0</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Threshold <IMG
 WIDTH="13" HEIGHT="15" ALIGN="BOTTOM" BORDER="0"
 SRC="img93.png"
 ALT="$\theta$"> in the aggregation algorithm.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_omega_alg_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'EIG_EST'</TT>   <TT>'USER_CHOICE'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'EIG_EST'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>How the damping parameter <IMG
 WIDTH="16" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
 SRC="img86.png"
 ALT="$\omega$"> in the
                         smoothed aggregation should be computed:
                         either via an estimate of the spectral radius of
                         <IMG
 WIDTH="50" HEIGHT="21" ALIGN="BOTTOM" BORDER="0"
 SRC="img87.png"
 ALT="$D^{-1}A$">, or explicily
                         specified by the user.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_eig_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'A_NORMI'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><TT>'A_NORMI'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>How to estimate the spectral radius of <IMG
 WIDTH="50" HEIGHT="21" ALIGN="BOTTOM" BORDER="0"
 SRC="img87.png"
 ALT="$D^{-1}A$">.
                           Currently only the infinity norm estimate
                           is available.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_aggr_omega_val_</code></TD>
<TD ALIGN="LEFT"><code>real(</code><I>kind_parameter</I><code>)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68>Any&nbsp;nonnegative&nbsp;real&nbsp;num.</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=68><!-- MATH
 $4/(3\rho(D^{-1}A))$
 -->
<IMG
 WIDTH="113" HEIGHT="39" ALIGN="MIDDLE" BORDER="0"
 SRC="img94.png"
 ALT="$4/(3\rho(D^{-1}A))$"></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Damping parameter <IMG
 WIDTH="16" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
 SRC="img86.png"
 ALT="$\omega$"> in the smoothed aggregation algorithm. 
                           It must be set by the user if
                           <code>USER_CHOICE</code> was specified for 
                           <code>mld_aggr_omega_alg_</code>,
                           otherwise it is computed by the library, using the
                           selected estimate of the spectral radius <IMG
 WIDTH="73" HEIGHT="39" ALIGN="MIDDLE" BORDER="0"
 SRC="img95.png"
 ALT="$\rho(D^{-1}A)$"> of
                           <IMG
 WIDTH="50" HEIGHT="21" ALIGN="BOTTOM" BORDER="0"
 SRC="img87.png"
 ALT="$D^{-1}A$">.</TD>
</TR>
</TABLE>
</DIV>
</TD></TR>
</TABLE>
</DIV><P></P>
<BR>

<P>
<BR><P></P>
<DIV ALIGN="CENTER"><A NAME="1265"></A>
<TABLE>
<CAPTION><STRONG>Table 5:</STRONG>
Parameters defining the coarse-space correction at the coarsest
level.</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
<TABLE CELLPADDING=3 BORDER="1" ALIGN="CENTER">
<TR><TD ALIGN="LEFT"><code>what</code></TD>
<TD ALIGN="LEFT"><SMALL>DATA TYPE</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><code>val</code></TD>
<TD ALIGN="LEFT"><SMALL>DEFAULT</SMALL></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198><SMALL>COMMENTS</SMALL></TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_mat_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'DISTR'</TT>   <TT>'REPL'</TT></TD>
<TD ALIGN="LEFT"><TT>'DISTR'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Coarsest matrix: distributed among the processors or
                           replicated on each of them.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_solve_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'BJAC'</TT>   <TT>'UMF'</TT>  
                           <TT>'SLU'</TT>   <TT>'SLUDIST'</TT></TD>
<TD ALIGN="LEFT"><TT>'BJAC'</TT></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Solver used at the coarsest level: block Jacobi, sequential
                           LU from UMFPACK, sequential LU from SuperLU, 
                           distributed LU from SuperLU_Dist.
                           <TT>'SLUDIST'</TT> requires the coarsest 
                           matrix to be distributed, while <TT>'UMF'</TT> and
                           <TT>'SLU'</TT> require it to be replicated.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_subsolve_</code></TD>
<TD ALIGN="LEFT"><code>character(len=*)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91><TT>'ILU'</TT>   <TT>'MILU'</TT>
                             <TT>'ILUT'</TT>
                             <TT>'UMF'</TT>   <TT>'SLU'</TT></TD>
<TD ALIGN="LEFT">See note</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Solver for the diagonal blocks of the coarse matrix,
                           in case the block Jacobi solver
                           is chosen as coarsest-level solver: ILU(<IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$">), MILU(<IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$">),
                           ILU(<IMG
 WIDTH="27" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img36.png"
 ALT="$p,t$">), LU from UMFPACK,
                           LU from SuperLU, plus triangular solve.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_sweeps_</code></TD>
<TD ALIGN="LEFT"><code>integer</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>Any&nbsp;int.&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img96.png"
 ALT="$&gt; 0$"></TD>
<TD ALIGN="LEFT">4</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Number of Block-Jacobi sweeps when 'BJAC' is used as
                           coarsest-level solver.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_fillin_</code></TD>
<TD ALIGN="LEFT"><code>integer</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>Any&nbsp;int.&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img90.png"
 ALT="$\ge 0$"></TD>
<TD ALIGN="LEFT">0</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Fill-in level <IMG
 WIDTH="13" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img35.png"
 ALT="$p$"> of the incomplete LU factorizations.</TD>
</TR>
<TR><TD ALIGN="LEFT"><code>mld_coarse_iluthrs_</code></TD>
<TD ALIGN="LEFT"><code>real(</code><I>kind_parameter</I><code>)</code></TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=91>Any&nbsp;real.&nbsp;num.&nbsp;<IMG
 WIDTH="32" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img90.png"
 ALT="$\ge 0$"></TD>
<TD ALIGN="LEFT">0</TD>
<TD ALIGN="LEFT" VALIGN="TOP" WIDTH=198>Drop tolerance <IMG
 WIDTH="11" HEIGHT="18" ALIGN="BOTTOM" BORDER="0"
 SRC="img91.png"
 ALT="$t$"> in the ILU(<IMG
 WIDTH="27" HEIGHT="31" ALIGN="MIDDLE" BORDER="0"
 SRC="img36.png"
 ALT="$p,t$">) factorization.</TD>
</TR>
<TR><TD ALIGN="LEFT" COLSPAN=5><B>Note:</B> defaults for
  <TT>m</TT>ld_coarse_subsolve_ are chosen as</TD>
</TR>
<TR><TD ALIGN="LEFT" COLSPAN=5>single precision version: 'SLU' if installed, 'ILU' otherwise</TD>
</TR>
<TR><TD ALIGN="LEFT" COLSPAN=5>double precision version: 'UMF' if installed,
  else 'SLU' if installed, 'ILU' otherwise</TD>
</TR>
</TABLE>
</DIV>
</TD></TR>
</TABLE>
</DIV><P></P>
<BR>

<P>

<P>
<HR>
<!--Navigation Panel-->
<A NAME="tex2html284"
  HREF="node19.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html280"
  HREF="node16.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html274"
  HREF="node17.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html282"
  HREF="node2.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>  
<BR>
<B> Next:</B> <A NAME="tex2html285"
  HREF="node19.html">Subroutine mld_precbld</A>
<B> Up:</B> <A NAME="tex2html281"
  HREF="node16.html">User Interface</A>
<B> Previous:</B> <A NAME="tex2html275"
  HREF="node17.html">Subroutine mld_precinit</A>
 &nbsp; <B>  <A NAME="tex2html283"
  HREF="node2.html">Contents</A></B> 
<!--End of Navigation Panel-->

</BODY>
</HTML>