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

393 lines
20 KiB
HTML

<!DOCTYPE html PUBLIC "-//W3C//DTD HTML 4.01 Transitional//EN"
"http://www.w3.org/TR/html4/loose.dtd">
<html >
<head><title>Examples</title>
<meta http-equiv="Content-Type" content="text/html; charset=iso-8859-1">
<meta name="generator" content="TeX4ht (http://www.tug.org/tex4ht/)">
<meta name="originator" content="TeX4ht (http://www.tug.org/tex4ht/)">
<!-- html,3 -->
<meta name="src" content="userhtml.tex">
<link rel="stylesheet" type="text/css" href="userhtml.css">
</head><body
>
<!--l. 109--><div class="crosslinks"><p class="noindent"><span
class="cmbx-12">[</span><a
href="#tailuserhtmlsu9.html"><span
class="cmbx-12">tail</span></a><span
class="cmbx-12">] [</span><a
href="userhtmlse5.html#userhtmlsu9.html" ><span
class="cmbx-12">up</span></a><span
class="cmbx-12">] </span></p></div>
<h4 class="subsectionHead"><span class="titlemark"><span
class="cmbx-12">5.1 </span></span> <a
id="x18-170005.1"></a><span
class="cmbx-12">Examples</span></h4>
<!--l. 111--><p class="noindent" ><span
class="cmbx-12">The code reported in Figure</span><span
class="cmbx-12">&#x00A0;</span><a
href="#x18-170012"><span
class="cmbx-12">2</span><!--tex4ht:ref: fig:ex1 --></a> <span
class="cmbx-12">shows how to set and apply the default</span>
<span
class="cmbx-12">multilevel preconditioner available in the real double precision version of</span>
<span
class="cmbx-12">AMG4PSBLAS (see Table</span><span
class="cmbx-12">&#x00A0;</span><a
href="userhtmlse5.html#x17-160151"><span
class="cmbx-12">1</span><!--tex4ht:ref: tab:precinit --></a><span
class="cmbx-12">). This preconditioner is chosen by simply</span>
<span
class="cmbx-12">specifying </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">&#8217;ML&#8217;</span></span></span> <span
class="cmbx-12">as the second argument of </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">P%init</span></span></span> <span
class="cmbx-12">(a call to </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">P%set</span></span></span> <span
class="cmbx-12">is not</span>
<span
class="cmbx-12">needed) and is applied with the CG solver provided by PSBLAS (the</span>
<span
class="cmbx-12">matrix of the system to be solved is assumed to be positive definite). As</span>
<span
class="cmbx-12">previously observed, the modules </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">psb_base_mod</span></span></span><span
class="cmbx-12">, </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">amg_prec_mod</span></span></span> <span
class="cmbx-12">and</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">psb_krylov_mod</span></span></span> <span
class="cmbx-12">must be used by the example program.</span>
<!--l. 121--><p class="indent" > <span
class="cmbx-12">The part of the code concerning the reading and assembling of the</span>
<span
class="cmbx-12">sparse matrix and the right-hand side vector, performed through the</span>
<span
class="cmbx-12">PSBLAS routines for sparse matrix and vector management, is not</span>
<span
class="cmbx-12">reported here for brevity; the statements concerning the deallocation of the</span>
<span
class="cmbx-12">PSBLAS data structure are neglected too. The complete code can be found</span>
<span
class="cmbx-12">in the example program file </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">amg_dexample_ml.f90</span></span></span><span
class="cmbx-12">, in the directory</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">examples/fileread</span></span></span> <span
class="cmbx-12">of the AMG4PSBLAS implementation (see Section</span><span
class="cmbx-12">&#x00A0;</span><a
href="userhtmlsu5.html#x12-110003.5"><span
class="cmbx-12">3.5</span><!--tex4ht:ref: sec:ex_and_test --></a><span
class="cmbx-12">).</span>
<span
class="cmbx-12">A sample test problem along with the relevant input data is available in</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">examples/fileread/runs</span></span></span><span
class="cmbx-12">. For details on the use of the PSBLAS routines, see</span>
<span
class="cmbx-12">the PSBLAS User&#8217;s Guide</span><span
class="cmbx-12">&#x00A0;</span><span class="cite"><span
class="cmbx-12">[</span><a
href="userhtmlli4.html#XPSBLASGUIDE"><span
class="cmbx-12">15</span></a><span
class="cmbx-12">]</span></span><span
class="cmbx-12">.</span>
<!--l. 133--><p class="indent" > <span
class="cmbx-12">The setup and application of the default multilevel preconditioner for</span>
<span
class="cmbx-12">the real single precision and the complex, single and double precision,</span>
<span
class="cmbx-12">versions are obtained with straightforward modifications of the previous</span>
<span
class="cmbx-12">example (see Section</span><span
class="cmbx-12">&#x00A0;</span><a
href="userhtmlse6.html#x19-180006"><span
class="cmbx-12">6</span><!--tex4ht:ref: sec:userinterface --></a> <span
class="cmbx-12">for details). If these versions are installed, the</span>
<span
class="cmbx-12">corresponding codes are available in </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">examples/fileread/</span></span></span><span
class="cmbx-12">.</span>
<!--l. 139--><p class="indent" > <hr class="figure"><div class="figure"
>
<a
id="x18-170012"></a>
<div class="center"
>
<!--l. 140--><p class="noindent" >
<div class="minipage"><div class="verbatim" id="verbatim-6">
&#x00A0;&#x00A0;use&#x00A0;psb_base_mod
&#x00A0;<br />&#x00A0;&#x00A0;use&#x00A0;amg_prec_mod
&#x00A0;<br />&#x00A0;&#x00A0;use&#x00A0;psb_krylov_mod
&#x00A0;<br />...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;sparse&#x00A0;matrix
&#x00A0;<br />&#x00A0;&#x00A0;type(psb_dspmat_type)&#x00A0;::&#x00A0;A
&#x00A0;<br />!&#x00A0;sparse&#x00A0;matrix&#x00A0;descriptor
&#x00A0;<br />&#x00A0;&#x00A0;type(psb_desc_type)&#x00A0;&#x00A0;&#x00A0;::&#x00A0;desc_A
&#x00A0;<br />!&#x00A0;preconditioner
&#x00A0;<br />&#x00A0;&#x00A0;type(amg_dprec_type)&#x00A0;&#x00A0;::&#x00A0;P
&#x00A0;<br />!&#x00A0;right-hand&#x00A0;side&#x00A0;and&#x00A0;solution&#x00A0;vectors
&#x00A0;<br />&#x00A0;&#x00A0;type(psb_d_vect_type)&#x00A0;::&#x00A0;b,&#x00A0;x
&#x00A0;<br />...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;initialize&#x00A0;the&#x00A0;parallel&#x00A0;environment
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;psb_init(ictxt)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;psb_info(ictxt,iam,np)
&#x00A0;<br />...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;read&#x00A0;and&#x00A0;assemble&#x00A0;the&#x00A0;spd&#x00A0;matrix&#x00A0;A&#x00A0;and&#x00A0;the&#x00A0;right-hand&#x00A0;side&#x00A0;b
&#x00A0;<br />!&#x00A0;using&#x00A0;PSBLAS&#x00A0;routines&#x00A0;for&#x00A0;sparse&#x00A0;matrix&#x00A0;/&#x00A0;vector&#x00A0;management
&#x00A0;<br />...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;initialize&#x00A0;the&#x00A0;default&#x00A0;multilevel&#x00A0;preconditioner,&#x00A0;i.e.&#x00A0;V-cycle
&#x00A0;<br />!&#x00A0;with&#x00A0;basic&#x00A0;smoothed&#x00A0;aggregation,&#x00A0;1&#x00A0;hybrid&#x00A0;forward/backward
&#x00A0;<br />!&#x00A0;GS&#x00A0;sweep&#x00A0;as&#x00A0;pre/post-smoother&#x00A0;and&#x00A0;UMFPACK&#x00A0;as&#x00A0;coarsest-level
&#x00A0;<br />!&#x00A0;solver
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;build&#x00A0;the&#x00A0;preconditioner
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
&#x00A0;<br />
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;set&#x00A0;the&#x00A0;solver&#x00A0;parameters&#x00A0;and&#x00A0;the&#x00A0;initial&#x00A0;guess
&#x00A0;<br />&#x00A0;&#x00A0;...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;solve&#x00A0;Ax=b&#x00A0;with&#x00A0;preconditioned&#x00A0;CG
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;psb_krylov(&#8217;CG&#8217;,A,P,b,x,tol,desc_A,info)
&#x00A0;<br />&#x00A0;&#x00A0;...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;deallocate&#x00A0;the&#x00A0;preconditioner
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%free(info)
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;deallocate&#x00A0;other&#x00A0;data&#x00A0;structures
&#x00A0;<br />&#x00A0;&#x00A0;...&#x00A0;...
&#x00A0;<br />!
&#x00A0;<br />!&#x00A0;exit&#x00A0;the&#x00A0;parallel&#x00A0;environment
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;psb_exit(ictxt)
&#x00A0;<br />&#x00A0;&#x00A0;stop
</div>
<!--l. 195--><p class="nopar" >
</div>
<br /> <div class="caption"
><span class="id">Figure&#x00A0;2: </span><span
class="content">setup and application of the default multilevel preconditioner (example 1). </span></div><!--tex4ht:label?: x18-170012 -->
</div>
<!--l. 201--><p class="indent" > </div><hr class="endfigure">
<!--l. 203--><p class="indent" > <span
class="cmbx-12">Different versions of the multilevel preconditioner can be obtained by</span>
<span
class="cmbx-12">changing the default values of the preconditioner parameters. The code</span>
<span
class="cmbx-12">reported in Figure</span><span
class="cmbx-12">&#x00A0;</span><a
href="#x18-170023"><span
class="cmbx-12">3</span><!--tex4ht:ref: fig:ex2 --></a> <span
class="cmbx-12">shows how to set a V-cycle preconditioner which</span>
<span
class="cmbx-12">applies 1 block-Jacobi sweep as pre- and post-smoother, and solves the</span>
<span
class="cmbx-12">coarsest-level system with 8 block-Jacobi sweeps. Note that the ILU(0)</span>
<span
class="cmbx-12">factorization (plus triangular solve) is used as local solver for the</span>
<span
class="cmbx-12">block-Jacobi sweeps, since this is the default associated with block-Jacobi</span>
<span
class="cmbx-12">and set by</span><span
class="cmbx-12">&#x00A0;</span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">P%init</span></span></span><span
class="cmbx-12">. Furthermore, specifying block-Jacobi as coarsest-level</span>
<span
class="cmbx-12">solver implies that the coarsest-level matrix is distributed among the</span>
<span
class="cmbx-12">processes. Figure</span><span
class="cmbx-12">&#x00A0;</span><a
href="#x18-170034"><span
class="cmbx-12">4</span><!--tex4ht:ref: fig:ex3 --></a> <span
class="cmbx-12">shows how to set a W-cycle preconditioner which</span>
<span
class="cmbx-12">applies 2 hybrid Gauss-Seidel sweeps as pre- and post-smoother, and solves</span>
<span
class="cmbx-12">the coarsest-level system with the multifrontal LU factorization</span>
<span
class="cmbx-12">implemented in MUMPS. It is specified that the coarsest-level matrix is</span>
<span
class="cmbx-12">distributed, since MUMPS can be used on both replicated and distributed</span>
<span
class="cmbx-12">matrices, and by default it is used on replicated ones. The code fragments</span>
<span
class="cmbx-12">shown in Figures</span><span
class="cmbx-12">&#x00A0;</span><a
href="#x18-170023"><span
class="cmbx-12">3</span><!--tex4ht:ref: fig:ex2 --></a> <span
class="cmbx-12">and </span><a
href="#x18-170034"><span
class="cmbx-12">4</span><!--tex4ht:ref: fig:ex3 --></a> <span
class="cmbx-12">are included in the example program file</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">amg_dexample_ml.f90</span></span></span> <span
class="cmbx-12">too.</span>
<!--l. 227--><p class="indent" > <span
class="cmbx-12">Finally, Figure</span><span
class="cmbx-12">&#x00A0;</span><a
href="#x18-170045"><span
class="cmbx-12">5</span><!--tex4ht:ref: fig:ex4 --></a> <span
class="cmbx-12">shows the setup of a one-level additive Schwarz</span>
<span
class="cmbx-12">preconditioner, i.e., RAS with overlap 2. Note also that a Krylov method</span>
<span
class="cmbx-12">different from CG must be used to solve the preconditioned system, since</span>
<span
class="cmbx-12">the preconditione in nonsymmetric. The corresponding example program is</span>
<span
class="cmbx-12">available in the file </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">amg_dexample_1lev.f90</span></span></span><span
class="cmbx-12">.</span>
<!--l. 234--><p class="indent" > <span
class="cmbx-12">For all the previous preconditioners, example programs where the sparse</span>
<span
class="cmbx-12">matrix and the right-hand side are generated by discretizing a PDE with</span>
<span
class="cmbx-12">Dirichlet boundary conditions are also available in the directory</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">examples/pdegen</span></span></span><span
class="cmbx-12">.</span>
<!--l. 238--><p class="indent" > <hr class="figure"><div class="figure"
>
<a
id="x18-170023"></a>
<div class="center"
>
<!--l. 239--><p class="noindent" >
<div class="minipage"><div class="verbatim" id="verbatim-7">
...&#x00A0;...
&#x00A0;<br />!&#x00A0;build&#x00A0;a&#x00A0;V-cycle&#x00A0;preconditioner&#x00A0;with&#x00A0;1&#x00A0;block-Jacobi&#x00A0;sweep&#x00A0;(with
&#x00A0;<br />!&#x00A0;ILU(0)&#x00A0;on&#x00A0;the&#x00A0;blocks)&#x00A0;as&#x00A0;pre-&#x00A0;and&#x00A0;post-smoother,&#x00A0;and&#x00A0;8&#x00A0;&#x00A0;block-Jacobi
&#x00A0;<br />!&#x00A0;sweeps&#x00A0;(with&#x00A0;ILU(0)&#x00A0;on&#x00A0;the&#x00A0;blocks)&#x00A0;as&#x00A0;coarsest-level&#x00A0;solver
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call_P%set(&#8217;SMOOTHER_TYPE&#8217;,&#8217;BJAC&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SOLVE&#8217;,&#8217;BJAC&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SWEEPS&#8217;,8,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
&#x00A0;<br />...&#x00A0;...
</div>
<!--l. 254--><p class="nopar" ></div>
<br /> <div class="caption"
><span class="id">Figure&#x00A0;3: </span><span
class="content">setup of a multilevel preconditioner</span></div><!--tex4ht:label?: x18-170023 -->
</div>
<!--l. 260--><p class="indent" > </div><hr class="endfigure">
<!--l. 262--><p class="indent" > <hr class="figure"><div class="figure"
>
<a
id="x18-170034"></a>
<div class="center"
>
<!--l. 263--><p class="noindent" >
<div class="minipage"><div class="verbatim" id="verbatim-8">
...&#x00A0;...
&#x00A0;<br />!&#x00A0;build&#x00A0;a&#x00A0;W-cycle&#x00A0;preconditioner&#x00A0;with&#x00A0;2&#x00A0;hybrid&#x00A0;Gauss-Seidel&#x00A0;sweeps
&#x00A0;<br />!&#x00A0;as&#x00A0;pre-&#x00A0;and&#x00A0;post-smoother,&#x00A0;a&#x00A0;distributed&#x00A0;coarsest
&#x00A0;<br />!&#x00A0;matrix,&#x00A0;and&#x00A0;MUMPS&#x00A0;as&#x00A0;coarsest-level&#x00A0;solver
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;ML_CYCLE&#8217;,&#8217;WCYCLE&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SMOOTHER_TYPE&#8217;,&#8217;FBGS&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SMOOTHER_SWEEPS&#8217;,2,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SOLVE&#8217;,&#8217;MUMPS&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_MAT&#8217;,&#8217;DIST&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
&#x00A0;<br />...&#x00A0;...
</div>
<!--l. 280--><p class="nopar" ></div>
<br /> <div class="caption"
><span class="id">Figure&#x00A0;4: </span><span
class="content">setup of a multilevel preconditioner</span></div><!--tex4ht:label?: x18-170034 -->
</div>
<!--l. 285--><p class="indent" > </div><hr class="endfigure">
<!--l. 287--><p class="indent" > <hr class="figure"><div class="figure"
>
<a
id="x18-170045"></a>
<div class="center"
>
<!--l. 288--><p class="noindent" >
<div class="minipage"><div class="verbatim" id="verbatim-9">
...&#x00A0;...
&#x00A0;<br />!&#x00A0;set&#x00A0;RAS&#x00A0;with&#x00A0;overlap&#x00A0;2&#x00A0;and&#x00A0;ILU(0)&#x00A0;on&#x00A0;the&#x00A0;local&#x00A0;blocks
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;AS&#8217;,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SUB_OVR&#8217;,2,info)
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;P%bld(A,desc_A,info)
&#x00A0;<br />...&#x00A0;...
&#x00A0;<br />!&#x00A0;solve&#x00A0;Ax=b&#x00A0;with&#x00A0;preconditioned&#x00A0;BiCGSTAB
&#x00A0;<br />&#x00A0;&#x00A0;call&#x00A0;psb_krylov(&#8217;BICGSTAB&#8217;,A,P,b,x,tol,desc_A,info)
</div>
<!--l. 300--><p class="nopar" ></div>
<br /> <div class="caption"
><span class="id">Figure&#x00A0;5: </span><span
class="content">setup of a one-level Schwarz preconditioner.</span></div><!--tex4ht:label?: x18-170045 -->
</div>
<!--l. 305--><p class="indent" > </div><hr class="endfigure">
<!--l. 1--><div class="crosslinks"><p class="noindent"><span
class="cmbx-12">[</span><a
href="userhtmlsu9.html" ><span
class="cmbx-12">front</span></a><span
class="cmbx-12">] [</span><a
href="userhtmlse5.html#userhtmlsu9.html" ><span
class="cmbx-12">up</span></a><span
class="cmbx-12">] </span></p></div>
<!--l. 1--><p class="indent" > <a
id="tailuserhtmlsu9.html"></a>
</body></html>