< !DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN"
"http://www.w3.org/TR/html4/loose.dtd">
< html >
< head > < title > Preconditioner routines< / title >
< meta http-equiv = "Content-Type" content = "text/html; charset=iso-8859-1" >
< meta name = "generator" content = "TeX4ht (https://tug.org/tex4ht/)" >
< meta name = "originator" content = "TeX4ht (https://tug.org/tex4ht/)" >
<!-- html,3 -->
< meta name = "src" content = "userhtml.tex" >
< link rel = "stylesheet" type = "text/css" href = "userhtml.css" >
< / head > < body
>
<!-- l. 1 --> < div class = "crosslinks" > < p class = "noindent" > [< a
href="userhtmlse9.html" >prev< / a > ] [< a
href="userhtmlse9.html#tailuserhtmlse9.html" >prev-tail< / a > ] [< a
href="userhtmlse7.html#tailuserhtmlse10.html">tail< / a > ] [< a
href="userhtml.html#userhtmlse13.html" >up< / a > ] < / p > < / div >
< h3 class = "sectionHead" > < span class = "titlemark" > 10 < / span > < a
id="x15-13500010">< / a > Preconditioner routines< / h3 >
<!-- l. 6 --> < p class = "noindent" > The base PSBLAS library contains the implementation of two simple preconditioning
techniques:
< ul class = "itemize1" >
< li class = "itemize" >
<!-- l. 9 --> < p class = "noindent" > Diagonal Scaling
< / li >
< li class = "itemize" >
<!-- l. 10 --> < p class = "noindent" > Block Jacobi with ILU(0) factorization< / li > < / ul >
<!-- l. 14 --> < p class = "noindent" > The supporting data type and subroutine interfaces are defined in the module
< span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">psb_prec_mod< / span > < / span > < / span > . The old interfaces < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">psb_precinit< / span > < / span > < / span > and < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">psb_precbld< / span > < / span > < / span > are still supported
for backward compatibility
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.1 < / span > < a
id="x15-13600010.1">< / a > init — Initialize a preconditioner< / h4 >
< pre class = "verbatim" id = "verbatim-97" >
call  prec%init(icontxt,ptype,  info)
< / pre >
<!-- l. 30 --> < p class = "nopar" >
<!-- l. 32 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 33 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 33 --> < p class = "noindent" > Asynchronous.
< / dd > < dt class = "description" >
<!-- l. 34 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 34 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 35 --> < p class = "noindent" >
< span
class="pplb7t-">icontxt< / span > < / dt > < dd
class="description">
<!-- l. 35 --> < p class = "noindent" > the communication context.< br
class="newline" />Scope:< span
class="pplb7t-">global< / span > .< br
class="newline" />Type:< span
class="pplb7t-">required< / span > .< br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an integer value.
< / dd > < dt class = "description" >
<!-- l. 40 --> < p class = "noindent" >
< span
class="pplb7t-">ptype< / span > < / dt > < dd
class="description">
<!-- l. 40 --> < p class = "noindent" > the type of preconditioner. Scope: < span
class="pplb7t-">global < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: a character string, see usage notes.
< / dd > < dt class = "description" >
<!-- l. 53 --> < p class = "noindent" >
< span
class="pplb7t-">On Exit< / span > < / dt > < dd
class="description">
<!-- l. 53 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 55 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 55 --> < p class = "noindent" > Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a preconditioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 60 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 60 --> < p class = "noindent" > Scope: < span
class="pplb7t-">global < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">out< / span > .< br
class="newline" />Error code: if no error, 0 is returned.< / dd > < / dl >
<!-- l. 66 --> < p class = "noindent" > < span
class="pplb7t-x-x-120">Notes < / span > Legal inputs to this subroutine are interpreted depending on the < span
class="zplmr7m-">ptype < / span > string as
follows< span class = "footnote-mark" > < a
href="userhtml16.html#fn4x0">< sup class = "textsuperscript" > 4< / sup > < / a > < / span > < a
id="x15-136001f4">< / a > :
< dl class = "description" > < dt class = "description" >
<!-- l. 74 --> < p class = "noindent" >
< span
class="pplb7t-">NONE< / span > < / dt > < dd
class="description">
<!-- l. 74 --> < p class = "noindent" > No preconditioning, i.e. the preconditioner is just a copy operator.
< / dd > < dt class = "description" >
<!-- l. 76 --> < p class = "noindent" >
< span
class="pplb7t-">DIAG< / span > < / dt > < dd
class="description">
<!-- l. 76 --> < p class = "noindent" > Diagonal scaling; each entry of the input vector is multiplied by the
reciprocal of the sum of the absolute values of the coefficients in the
corresponding row of matrix < span
class="zplmr7m-">A< / span > ;
< / dd > < dt class = "description" >
<!-- l. 79 --> < p class = "noindent" >
< span
class="pplb7t-">BJAC< / span > < / dt > < dd
class="description">
<!-- l. 79 --> < p class = "noindent" > Precondition by a factorization of the block-diagonal of matrix < span
class="zplmr7m-">A< / span > , where
block boundaries are determined by the data allocation boundaries
for each process; requires no communication. Only the incomplete
factorization < span
class="zplmr7m-">ILU< / span > < span
class="zplmr7t-">(< / span > 0< span
class="zplmr7t-">) < / span > is currently implemented.< / dd > < / dl >
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.2 < / span > < a
id="x15-13700010.2">< / a > Set — set preconditioner parameters< / h4 >
< div class = "center"
>
<!-- l. 92 --> < p class = "noindent" >
<!-- l. 93 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > call< / span > < span style = "color:#000000" > < / span > < span style = "color:#000000" > p< / span > < span style = "color:#000000" > %< / span > < span style = "color:#000000" > set< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > what< / span > < span style = "color:#000000" > ,< / span > < span style = "color:#000000" > val< / span > < span style = "color:#000000" > ,< / span > < span style = "color:#000000" > info< / span > < span style = "color:#000000" > )< / span > < / code > < / div >
<!-- l. 96 --> < p class = "noindent" > This method sets the parameters defining the subdomain solver when the
preconditioner type is < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">BJAC< / span > < / span > < / span > . More precisely, the parameter identified by < code class = "lstinline" > < span style = "color:#000000" > what< / span > < / code > is
assigned the value contained in < code class = "lstinline" > < span style = "color:#000000" > val< / span > < / code > .
<!-- l. 102 --> < p class = "noindent" > < span
class="pplb7t-x-x-120">Arguments< / span >
< div class = "tabular" > < table id = "TBL-23" class = "tabular"
>< colgroup id = "TBL-23-1g" > < col
id="TBL-23-1">< col
id="TBL-23-2">< / colgroup > < tr
style="vertical-align:baseline;" id="TBL-23-1-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-1-1"
class="td11"><!-- l. 105 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > what< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-1-2"
class="td11"><!-- l. 105 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > character< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > len< / span > < span style = "color:#000000" > =*)< / span > < / code > . < / td >
< / tr > < tr
style="vertical-align:baseline;" id="TBL-23-2-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-2-1"
class="td11"><!-- l. 106 --> < p class = "noindent" > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-2-2"
class="td11"><!-- l. 106 --> < p class = "noindent" > The parameter to be set. It can be specified through its name; the string is
case-insensitive. See Tables  < span
class="pplb7t-">??< / span > -< a
href="#x15-137002r22">22<!-- tex4ht:ref: tab:p_smoother_1 --> < / a > . < / td >
< / tr > < tr
style="vertical-align:baseline;" id="TBL-23-3-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-3-1"
class="td11"><!-- l. 109 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > val< / span > < span style = "color:#000000" > < / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-3-2"
class="td11"><!-- l. 109 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > integer< / span > < / code > < span
class="pplri7t-">or < / span > < code class = "lstinline" > < span style = "color:#000000" > character< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > len< / span > < span style = "color:#000000" > =*)< / span > < / code > < span
class="pplri7t-">or < / span > < code class = "lstinline" > < span style = "color:#000000" > real< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > psb_spk_< / span > < span style = "color:#000000" > )< / span > < / code > < span
class="pplri7t-">or < / span > < code class = "lstinline" > < span style = "color:#000000" > real< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > psb_dpk_< / span > < span style = "color:#000000" > )< / span > < / code > ,
< code class = "lstinline" > < span style = "color:#000000" > intent< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > in< / span > < span style = "color:#000000" > )< / span > < / code > . < / td >
< / tr > < tr
style="vertical-align:baseline;" id="TBL-23-4-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-4-1"
class="td11"><!-- l. 112 --> < p class = "noindent" > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-4-2"
class="td11"><!-- l. 112 --> < p class = "noindent" > The value of the parameter to be set. The list of allowed values and the
corresponding data types is given in Tables  < span
class="pplb7t-">??< / span > -< a
href="#x15-137002r22">22<!-- tex4ht:ref: tab:p_smoother_1 --> < / a > . When the value is of type
< code class = "lstinline" > < span style = "color:#000000" > character< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > len< / span > < span style = "color:#000000" > =*)< / span > < / code > , it is also treated as case insensitive. < / td >
< / tr > < tr
style="vertical-align:baseline;" id="TBL-23-5-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-5-1"
class="td11"><!-- l. 117 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > info< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-5-2"
class="td11"><!-- l. 117 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > integer< / span > < span style = "color:#000000" > ,< / span > < span style = "color:#000000" > < / span > < span style = "color:#000000" > intent< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > out< / span > < span style = "color:#000000" > )< / span > < / code > . < / td >
< / tr > < tr
style="vertical-align:baseline;" id="TBL-23-6-">< td style = "white-space:normal; text-align:left;" id = "TBL-23-6-1"
class="td11"><!-- l. 118 --> < p class = "noindent" > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-23-6-2"
class="td11"><!-- l. 118 --> < p class = "noindent" > Error code. If no error, 0 is returned. See Section  < span
class="pplb7t-">??< / span > for details. < / td > < / tr > < / table > < / div >
<!-- l. 123 --> < p class = "noindent" > A number of subdomain solvers can be chosen with this method; a list of the
parameters that can be set, along with their allowed and default values, is given in
Tables  < span
class="pplb7t-">??< / span > -< a
href="#x15-137002r22">22<!-- tex4ht:ref: tab:p_smoother_1 --> < / a > .< br
class="newline" />
< div class = "table" >
<!-- l. 130 --> < p class = "indent" > < a
id="x15-137001r21">< / a > < hr class = "float" > < div class = "float"
>
< div class = "center"
>
<!-- l. 130 --> < p class = "noindent" >
< div class = "tabular" > < table id = "TBL-24" class = "tabular"
>< colgroup id = "TBL-24-1g" > < col
id="TBL-24-1">< / colgroup > < colgroup id = "TBL-24-2g" > < col
id="TBL-24-2">< / colgroup > < colgroup id = "TBL-24-3g" > < col
id="TBL-24-3">< / colgroup > < colgroup id = "TBL-24-4g" > < col
id="TBL-24-4">< / colgroup > < colgroup id = "TBL-24-5g" > < col
id="TBL-24-5">< / colgroup > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-24-1-">< td style = "white-space:normal; text-align:left;" id = "TBL-24-1-1"
class="td11"><!-- l. 134 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > what< / span > < / code > < / td > < td style = "white-space:nowrap; text-align:left;" id = "TBL-24-1-2"
class="td11">< span
class="pplrc7t-x-x-90">< span
class="small-caps">d< / span > < span
class="small-caps">a< / span > < span
class="small-caps">t< / span > < span
class="small-caps">a< / span > < span
class="small-caps">t< / span > < span
class="small-caps">y< / span > < span
class="small-caps">p< / span > < span
class="small-caps">e< / span > < / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-1-3"
class="td11"><!-- l. 134 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > val< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-1-4"
class="td11"><!-- l. 134 --> < p class = "noindent" > < span
class="pplrc7t-x-x-90">< span
class="small-caps">d< / span > < span
class="small-caps">e< / span > < span
class="small-caps">f< / span > < span
class="small-caps">a< / span > < span
class="small-caps">u< / span > < span
class="small-caps">l< / span > < span
class="small-caps">t< / span > < / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-1-5"
class="td11"><!-- l. 135 --> < p class = "noindent" > < span
class="pplrc7t-x-x-90">< span
class="small-caps">c< / span > < span
class="small-caps">o< / span > < span
class="small-caps">m< / span > < span
class="small-caps">m< / span > < span
class="small-caps">e< / span > < span
class="small-caps">n< / span > < span
class="small-caps">t< / span > < span
class="small-caps">s< / span > < / span > < / td > < / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-24-2-">< td style = "white-space:normal; text-align:left;" id = "TBL-24-2-1"
class="td
11"><!-- l. 137 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > SUB_SOLVE< / span > < span style = "color:#000000" > ’ < / span > < / code > < / td > < td style = "white-space:nowrap; text-align:left;" id = "TBL-24-2-2"
class="td11">< code class = "lstinline" > <!-- l. 137 --> < p class = "noindent" > < span style = "color:#000000" > character< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > len< / span > < span style = "color:#000000" > =*)< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-2-3"
class="td11"><!-- l. 138 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > ILU< / span > < span style = "color:#000000" > ’ < / span > < / code >
<!-- l. 139 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > ILUT< / span > < span style = "color:#000000" > ’ < / span > < / code >
<!-- l. 140 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > INVT< / span > < span style = "color:#000000" > ’ < / span > < / code >
<!-- l. 140 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > INVK< / span > < span style = "color:#000000" > ’ < / span > < / code >
<!-- l. 140 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > AINV< / span > < span style = "color:#000000" > ’ < / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-2-4"
class="td11"><!-- l. 142 --> < p class = "noindent" > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-24-2-5"
class="td11"><!-- l. 142 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">The local solver to be used with the smoother< / span >
< span
class="pplr7t-x-x-90">or one-level preconditioner ILU(< / span > < span
class="zplmr7m-x-x-90">p< / span > < span
class="pplr7t-x-x-90">), ILU(< / span > < span
class="zplmr7m-x-x-90">p< / span > < span
class="pplr7t-x-x-90">, < / span > < span
class="zplmr7m-x-x-90">t< / span > < span
class="pplr7t-x-x-90">),< / span >
< span
class="pplr7t-x-x-90">Approximate Inverses< / span >
< span
class="pplr7t-x-x-90">INVK(< / span > < span
class="zplmr7m-x-x-90">p< / span > < span
class="pplr7t-x-x-90">, < / span > < span
class="zplmr7m-x-x-90">q< / span > < span
class="pplr7t-x-x-90">), INVT(< / span > < span
class="zplmr7m-x-x-90">p< / span > < sub > < span
class="pplr7t-x-x-70">1< / span > < / sub > < span
class="pplr7t-x-x-90">, < / span > < span
class="zplmr7m-x-x-90">p< / span > < span
class="pplr7t-x-x-90">2, < / span > < span
class="zplmr7m-x-x-90">t< / span > < sub > < span
class="pplr7t-x-x-70">1< / span > < / sub > < span
class="pplr7t-x-x-90">, < / span > < span
class="zplmr7m-x-x-90">t< / span > < sub > < span
class="pplr7t-x-x-70">2< / span > < / sub > < span
class="pplr7t-x-x-90">) and AINV(< / span > < span
class="zplmr7m-x-x-90">t< / span > < span
class="pplr7t-x-x-90">);< / span >
< span
class="pplr7t-x-x-90">note that approximate inverses are specifically< / span >
< span
class="pplr7t-x-x-90">suited for GPUs since they do not employ< / span >
< span
class="pplr7t-x-x-90">triangular system solve kernels, see< / span > < span
class="pplr7t-x-x-90">  < / span > < span class = "cite" > < span
class="pplr7t-x-x-90">[< / span > < span
class="pplb7t-x-x-90">?< / span > < span
class="pplr7t-x-x-90">]< / span > < / span > < span
class="pplr7t-x-x-90">.< / span > < / td >
< / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-24-3-">< td style = "white-space:normal; text-align:left;" id = "TBL-24-3-1"
class="td11"> < / td > < / tr > < / table > < / div > < / div >
< br / > < div class = "caption"
>< span class = "id" > Table  21: < / span > < span
class="content">Parameters defining the solver or the details of the one-level
preconditioner. < / span > < / div > <!-- tex4ht:label?: x15 - 137001r21 -->
< / div > < hr class = "endfloat" / >
< / div >
< div class = "table" >
<!-- l. 159 --> < p class = "indent" > < a
id="x15-137002r22">< / a > < hr class = "float" > < div class = "float"
>
< div class = "center"
>
<!-- l. 159 --> < p class = "noindent" >
< div class = "tabular" > < table id = "TBL-25" class = "tabular"
>< colgroup id = "TBL-25-1g" > < col
id="TBL-25-1">< / colgroup > < colgroup id = "TBL-25-2g" > < col
id="TBL-25-2">< / colgroup > < colgroup id = "TBL-25-3g" > < col
id="TBL-25-3">< / colgroup > < colgroup id = "TBL-25-4g" > < col
id="TBL-25-4">< / colgroup > < colgroup id = "TBL-25-5g" > < col
id="TBL-25-5">< / colgroup > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-25-1-">< td style = "white-space:normal; text-align:left;" id = "TBL-25-1-1"
class="td11"><!-- l. 163 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > what< / span > < / code > < / td > < td style = "white-space:nowrap; text-align:left;" id = "TBL-25-1-2"
class="td11">< span
class="pplrc7t-x-x-90">< span
class="small-caps">d< / span > < span
class="small-caps">a< / span > < span
class="small-caps">t< / span > < span
class="small-caps">a< / span > < span
class="small-caps">t< / span > < span
class="small-caps">y< / span > < span
class="small-caps">p< / span > < span
class="small-caps">e< / span > < / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-1-3"
class="td11"><!-- l. 163 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > val< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-1-4"
class="td11"><!-- l. 163 --> < p class = "noindent" > < span
class="pplrc7t-x-x-90">< span
class="small-caps">d< / span > < span
class="small-caps">e< / span > < span
class="small-caps">f< / span > < span
class="small-caps">a< / span > < span
class="small-caps">u< / span > < span
class="small-caps">l< / span > < span
class="small-caps">t< / span > < / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-1-5"
class="td11"><!-- l. 164 --> < p class = "noindent" > < span
class="pplrc7t-x-x-90">< span
class="small-caps">c< / span > < span
class="small-caps">o< / span > < span
class="small-caps">m< / span > < span
class="small-caps">m< / span > < span
class="small-caps">e< / span > < span
class="small-caps">n< / span > < span
class="small-caps">t< / span > < span
class="small-caps">s< / span > < / span > < / td > < / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-25-2-">< td style = "white-space:normal; text-align:left;" id = "TBL-25-2-1"
class="td
11"><!-- l. 165 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > SUB_FILLIN< / span > < span style = "color:#000000" > ’ < / span > < / code > < / td > < td style = "white-space:nowrap; text-align:left;" id = "TBL-25-2-2"
class="td11">< code class = "lstinline" > <!-- l. 165 --> < p class = "noindent" > < span style = "color:#000000" > integer< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-2-3"
class="td11"><!-- l. 166 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">Any integer< / span >
<!-- l. 166 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">number< / span > < span
class="pplr7t-x-x-90">  < / span > < span
class="zplmr7y-x-x-90">≥ < / span > < span
class="pplr7t-x-x-90">0< / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-2-4"
class="td11"><!-- l. 167 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">0< / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-2-5"
class="td11"><!-- l. 168 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">Fill-in level < / span > < span
class="zplmr7m-x-x-90">p < / span > < span
class="pplr7t-x-x-90">of the incomplete LU< / span >
< span
class="pplr7t-x-x-90">factorizations.< / span > < / td >
< / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-25-3-">< td style = "white-space:normal; text-align:left;" id = "TBL-25-3-1"
class="td11"><!-- l. 169 --> < p class = "noindent" > < code class = "lstinline" > < span style = "color:#000000" > ’ < / span > < span style = "color:#000000" > SUB_ILUTHRS< / span > < span style = "color:#000000" > ’ < / span > < / code > < / td > < td style = "white-space:nowrap; text-align:left;" id = "TBL-25-3-2"
class="td11">< code class = "lstinline" > <!-- l. 169 --> < p class = "noindent" > < span style = "color:#000000" > real< / span > < span style = "color:#000000" > (< / span > < span style = "color:#000000" > kind_parameter< / span > < span style = "color:#000000" > )< / span > < / code > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-3-3"
class="td11"><!-- l. 170 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">Any real< / span >
< span
class="pplr7t-x-x-90">number< / span > < span
class="pplr7t-x-x-90">  < / span > < span
class="zplmr7y-x-x-90">≥ < / span > < span
class="pplr7t-x-x-90">0< / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-3-4"
class="td11"><!-- l. 171 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">0< / span > < / td > < td style = "white-space:normal; text-align:left;" id = "TBL-25-3-5"
class="td11"><!-- l. 172 --> < p class = "noindent" > < span
class="pplr7t-x-x-90">Drop tolerance < / span > < span
class="zplmr7m-x-x-90">t < / span > < span
class="pplr7t-x-x-90">in the ILU(< / span > < span
class="zplmr7m-x-x-90">p< / span > < span
class="pplr7t-x-x-90">, < / span > < span
class="zplmr7m-x-x-90">t< / span > < span
class="pplr7t-x-x-90">) factorization.< / span > < / td >
< / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
class="hline">< td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < td > < hr > < / td > < / tr > < tr
style="vertical-align:baseline;" id="TBL-25-4-">< td style = "white-space:normal; text-align:left;" id = "TBL-25-4-1"
class="td11"> < / td > < / tr > < / table > < / div > < / div >
< br / > < div class = "caption"
>< span class = "id" > Table  22: < / span > < span
class="content">Parameters defining the smoother or the details of the one-level
preconditioner (continued).< / span > < / div > <!-- tex4ht:label?: x15 - 137002r22 -->
< / div > < hr class = "endfloat" / >
< / div >
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.3 < / span > < a
id="x15-13800010.3">< / a > build — Builds a preconditioner< / h4 >
< pre class = "verbatim" id = "verbatim-98" >
call  prec%build(a,  desc_a,  info[,amold,vmold,imold])
< / pre >
<!-- l. 187 --> < p class = "nopar" >
<!-- l. 189 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 190 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 190 --> < p class = "noindent" > Synchronous.
< / dd > < dt class = "description" >
<!-- l. 191 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 191 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 192 --> < p class = "noindent" >
< span
class="pplb7t-">a< / span > < / dt > < dd
class="description">
<!-- l. 192 --> < p class = "noindent" > the system sparse matrix. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > , target.< br
class="newline" />Specified as: a sparse matrix data structure < a
href="userhtmlse3.html#spdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_Tspmat< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 197 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 197 --> < p class = "noindent" > the preconditioner.< br
class="newline" />Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: an already initialized precondtioner data structure
< a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > < br
class="newline" />
< / dd > < dt class = "description" >
<!-- l. 202 --> < p class = "noindent" >
< span
class="pplb7t-">desc< / span > < span
class="pplb7t-">_a< / span > < / dt > < dd
class="description">
<!-- l. 202 --> < p class = "noindent" > the problem communication descriptor. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > , target.< br
class="newline" />Specified as: a communication descriptor data structure < a
href="userhtmlse3.html#descdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_desc< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 212 --> < p class = "noindent" >
< span
class="pplb7t-">amold< / span > < / dt > < dd
class="description">
<!-- l. 212 --> < p class = "noindent" > The desired dynamic type for the internal matrix storage.< br
class="newline" />Scope: < span
class="pplb7t-">local< / span > .< br
class="newline" />Type: < span
class="pplb7t-">optional< / span > .< br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an object of a class derived from < a
id="spbasedata">< / a > < span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_T< / span > < span
class="cmtt-10">_base< / span > < span
class="cmtt-10">_sparse< / span > < span
class="cmtt-10">_mat< / span > .
< / dd > < dt class = "description" >
<!-- l. 217 --> < p class = "noindent" >
< span
class="pplb7t-">vmold< / span > < / dt > < dd
class="description">
<!-- l. 217 --> < p class = "noindent" > The desired dynamic type for the internal vector storage.< br
class="newline" />Scope: < span
class="pplb7t-">local< / span > .< br
class="newline" />Type: < span
class="pplb7t-">optional< / span > .< br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an object of a class derived from < a
id="vbasedata">< / a > < span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_T< / span > < span
class="cmtt-10">_base< / span > < span
class="cmtt-10">_vect< / span > < span
class="cmtt-10">_type< / span > .
< / dd > < dt class = "description" >
<!-- l. 222 --> < p class = "noindent" >
< span
class="pplb7t-">imold< / span > < / dt > < dd
class="description">
<!-- l. 222 --> < p class = "noindent" > The desired dynamic type for the internal integer vector storage.< br
class="newline" />Scope: < span
class="pplb7t-">local< / span > .< br
class="newline" />Type: < span
class="pplb7t-">optional< / span > .< br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an object of a class derived from (integer)
< a
id="vbasedata">< / a > < span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_T< / span > < span
class="cmtt-10">_base< / span > < span
class="cmtt-10">_vect< / span > < span
class="cmtt-10">_type< / span > .< / dd > < / dl >
<!-- l. 229 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 230 --> < p class = "noindent" >
< span
class="pplb7t-">On Return< / span > < / dt > < dd
class="description">
<!-- l. 230 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 231 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 231 --> < p class = "noindent" > the preconditioner.< br
class="newline" />Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a precondtioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > < br
class="newline" />
< / dd > < dt class = "description" >
<!-- l. 236 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 236 --> < p class = "noindent" > Error code.< br
class="newline" />Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required < / span > < br
class="newline" />Intent: < span
class="pplb7t-">out< / span > .< br
class="newline" />An integer value; 0 means no error has been detected.< / dd > < / dl >
<!-- l. 242 --> < p class = "noindent" > The < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">amold< / span > < / span > < / span > , < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">vmold< / span > < / span > < / span > and < span class = "obeylines-h" > < span class = "verb" > < span
class="cmtt-10">imold< / span > < / span > < / span > arguments may be employed to interface with special
devices, such as GPUs and other accelerators.
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.4 < / span > < a
id="x15-13900010.4">< / a > apply — Preconditioner application routine< / h4 >
< pre class = "verbatim" id = "verbatim-99" >
call  prec%apply(x,y,desc_a,info,trans,work)
call  prec%apply(x,desc_a,info,trans)
< / pre >
<!-- l. 254 --> < p class = "nopar" >
<!-- l. 256 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 257 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 257 --> < p class = "noindent" > Synchronous.
< / dd > < dt class = "description" >
<!-- l. 258 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 258 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 259 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 259 --> < p class = "noindent" > the preconditioner. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: a preconditioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 264 --> < p class = "noindent" >
< span
class="pplb7t-">x< / span > < / dt > < dd
class="description">
<!-- l. 264 --> < p class = "noindent" > the source vector. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a rank one array or an object of type < a
href="userhtmlse3.html#vdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_T< / span > < span
class="cmtt-10">_vect< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 269 --> < p class = "noindent" >
< span
class="pplb7t-">desc< / span > < span
class="pplb7t-">_a< / span > < / dt > < dd
class="description">
<!-- l. 269 --> < p class = "noindent" > the problem communication descriptor. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: a communication data structure < a
href="userhtmlse3.html#descdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_desc< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 274 --> < p class = "noindent" >
< span
class="pplb7t-">trans< / span > < / dt > < dd
class="description">
<!-- l. 274 --> < p class = "noindent" > Scope: < br
class="newline" />Type: < span
class="pplb7t-">optional< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: a character.
< / dd > < dt class = "description" >
<!-- l. 279 --> < p class = "noindent" >
< span
class="pplb7t-">work< / span > < / dt > < dd
class="description">
<!-- l. 279 --> < p class = "noindent" > an optional work space Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">optional< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a double precision array.< / dd > < / dl >
<!-- l. 286 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 287 --> < p class = "noindent" >
< span
class="pplb7t-">On Return< / span > < / dt > < dd
class="description">
<!-- l. 287 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 288 --> < p class = "noindent" >
< span
class="pplb7t-">y< / span > < / dt > < dd
class="description">
<!-- l. 288 --> < p class = "noindent" > the destination vector. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a rank one array or an object of type < a
href="userhtmlse3.html#vdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_T< / span > < span
class="cmtt-10">_vect< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 293 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 293 --> < p class = "noindent" > Error code.< br
class="newline" />Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required < / span > < br
class="newline" />Intent: < span
class="pplb7t-">out< / span > .< br
class="newline" />An integer value; 0 means no error has been detected.< / dd > < / dl >
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.5 < / span > < a
id="x15-14000010.5">< / a > descr — Prints a description of current preconditioner< / h4 >
< pre class = "verbatim" id = "verbatim-100" >
call  prec%descr(info)
call  prec%descr(info,iout,  root)
< / pre >
<!-- l. 308 --> < p class = "nopar" >
<!-- l. 310 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 311 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 311 --> < p class = "noindent" > Asynchronous.
< / dd > < dt class = "description" >
<!-- l. 312 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 312 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 313 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 313 --> < p class = "noindent" > the preconditioner. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: a preconditioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 318 --> < p class = "noindent" >
< span
class="pplb7t-">iout< / span > < / dt > < dd
class="description">
<!-- l. 318 --> < p class = "noindent" > output unit. Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">optional< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an integer number. Default: default output unit.
< / dd > < dt class = "description" >
<!-- l. 323 --> < p class = "noindent" >
< span
class="pplb7t-">root< / span > < / dt > < dd
class="description">
<!-- l. 323 --> < p class = "noindent" > Process from which to print Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">optional< / span > < br
class="newline" />Intent: < span
class="pplb7t-">in< / span > .< br
class="newline" />Specified as: an integer number between 0 and < span
class="zplmr7m-">np< / span > < span
class="zplmr7y-">- < / span > 1, in which case
the specified process will print the description, or < span
class="zplmr7y-">-< / span > 1, in which case all
processes will print. Default: 0.
< / dd > < dt class = "description" >
<!-- l. 330 --> < p class = "noindent" >
< span
class="pplb7t-">On Return< / span > < / dt > < dd
class="description">
<!-- l. 330 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 331 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 331 --> < p class = "noindent" > Error code.< br
class="newline" />Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required < / span > < br
class="newline" />Intent: < span
class="pplb7t-">out< / span > .< br
class="newline" />An integer value; 0 means no error has been detected.< / dd > < / dl >
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.6 < / span > < a
id="x15-14100010.6">< / a > clone — clone current preconditioner< / h4 >
< pre class = "verbatim" id = "verbatim-101" >
call    prec%clone(precout,info)
< / pre >
<!-- l. 344 --> < p class = "nopar" >
<!-- l. 346 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 347 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 347 --> < p class = "noindent" > Asynchronous.
< / dd > < dt class = "description" >
<!-- l. 348 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 348 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 349 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 349 --> < p class = "noindent" > the preconditioner.< br
class="newline" />Scope: < span
class="pplb7t-">local< / span > .< br
class="newline" />< / dd > < / dl >
<!-- l. 356 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 357 --> < p class = "noindent" >
< span
class="pplb7t-">On Return< / span > < / dt > < dd
class="description">
<!-- l. 357 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 358 --> < p class = "noindent" >
< span
class="pplb7t-">precout< / span > < / dt > < dd
class="description">
<!-- l. 358 --> < p class = "noindent" > A copy of the input object.
< / dd > < dt class = "description" >
<!-- l. 359 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 359 --> < p class = "noindent" > Return code.< / dd > < / dl >
< h4 class = "subsectionHead" > < span class = "titlemark" > 10.7 < / span > < a
id="x15-14200010.7">< / a > free — Free a preconditioner< / h4 >
< pre class = "verbatim" id = "verbatim-102" >
call  prec%free(info)
< / pre >
<!-- l. 367 --> < p class = "nopar" >
<!-- l. 369 --> < p class = "indent" >
< dl class = "description" > < dt class = "description" >
<!-- l. 370 --> < p class = "noindent" >
< span
class="pplb7t-">Type:< / span > < / dt > < dd
class="description">
<!-- l. 370 --> < p class = "noindent" > Asynchronous.
< / dd > < dt class = "description" >
<!-- l. 371 --> < p class = "noindent" >
< span
class="pplb7t-">On Entry< / span > < / dt > < dd
class="description">
<!-- l. 371 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 372 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 372 --> < p class = "noindent" > the preconditioner.< br
class="newline" />Scope: < span
class="pplb7t-">local< / span > .< br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a preconditioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 385 --> < p class = "noindent" >
< span
class="pplb7t-">On Exit< / span > < / dt > < dd
class="description">
<!-- l. 385 --> < p class = "noindent" >
< / dd > < dt class = "description" >
<!-- l. 387 --> < p class = "noindent" >
< span
class="pplb7t-">prec< / span > < / dt > < dd
class="description">
<!-- l. 387 --> < p class = "noindent" > Scope: < span
class="pplb7t-">local < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">inout< / span > .< br
class="newline" />Specified as: a preconditioner data structure < a
href="userhtmlse3.html#precdata">< span
class="cmtt-10">psb< / span > < span
class="cmtt-10">_prec< / span > < span
class="cmtt-10">_type< / span > < / a > .
< / dd > < dt class = "description" >
<!-- l. 392 --> < p class = "noindent" >
< span
class="pplb7t-">info< / span > < / dt > < dd
class="description">
<!-- l. 392 --> < p class = "noindent" > Scope: < span
class="pplb7t-">global < / span > < br
class="newline" />Type: < span
class="pplb7t-">required< / span > < br
class="newline" />Intent: < span
class="pplb7t-">out< / span > .< br
class="newline" />Error code: if no error, 0 is returned.< / dd > < / dl >
<!-- l. 398 --> < p class = "noindent" > < span
class="pplb7t-x-x-120">Notes < / span > Releases all internal storage.
<!-- l. 1 --> < div class = "crosslinks" > < p class = "noindent" > [< a
href="userhtmlse9.html" >prev< / a > ] [< a
href="userhtmlse9.html#tailuserhtmlse9.html" >prev-tail< / a > ] [< a
href="userhtmlse10.html" >front< / a > ] [< a
href="userhtml.html#userhtmlse13.html" >up< / a > ] < / p > < / div >
<!-- l. 1 --> < p class = "indent" > < a
id="tailuserhtmlse10.html">< / a >
< / body > < / html >