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/userhtmlsu6.html

400 lines
19 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 (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. 114--><div class="crosslinks"><p class="noindent"><span
class="cmr-12">[</span><a
href="userhtmlsu7.html" ><span
class="cmr-12">next</span></a><span
class="cmr-12">] [</span><a
href="#tailuserhtmlsu6.html"><span
class="cmr-12">tail</span></a><span
class="cmr-12">] [</span><a
href="userhtmlse4.html#userhtmlsu6.html" ><span
class="cmr-12">up</span></a><span
class="cmr-12">] </span></p></div>
<h4 class="subsectionHead"><span class="titlemark"><span
class="cmr-12">4.1 </span></span> <a
id="x15-140004.1"></a><span
class="cmr-12">Examples</span></h4>
<!--l. 116--><p class="noindent" ><span
class="cmr-12">The code reported in Figure</span><span
class="cmr-12">&#x00A0;</span><a
href="#x15-14001r1"><span
class="cmr-12">1</span><!--tex4ht:ref: fig:ex1 --></a> <span
class="cmr-12">shows how to set and apply the default multilevel</span>
<span
class="cmr-12">preconditioner available in the real double precision version of AMG4PSBLAS</span>
<span
class="cmr-12">(see Table</span><span
class="cmr-12">&#x00A0;</span><a
href="userhtmlse4.html#x14-13015r1"><span
class="cmr-12">1</span><!--tex4ht:ref: tab:precinit --></a><span
class="cmr-12">). This preconditioner is chosen by simply specifying </span><span class="lstinline"></span><span
class="cmtt-12">&#8217;</span><span
class="cmtt-12">ML</span><span
class="cmtt-12">&#8217;</span> <span
class="cmr-12">as the</span>
<span
class="cmr-12">second argument of </span><span class="lstinline"></span><span
class="cmtt-12">P</span><span
class="cmtt-12">%</span><span
class="cmtt-12">init</span> <span
class="cmr-12">(a call to </span><span class="lstinline"></span><span
class="cmtt-12">P</span><span
class="cmtt-12">%</span><span
class="cmtt-12">set</span> <span
class="cmr-12">is not needed) and is applied</span>
<span
class="cmr-12">with the CG solver provided by PSBLAS (the matrix of the system to be</span>
<span
class="cmr-12">solved is assumed to be positive definite). As previously observed, the modules</span>
<span class="lstinline"></span><span
class="cmtt-12">psb_base_mod</span><span
class="cmr-12">, </span><span class="lstinline"></span><span
class="cmtt-12">amg_prec_mod</span> <span
class="cmr-12">and </span><span class="lstinline"></span><span
class="cmtt-12">psb_krylov_mod</span> <span
class="cmr-12">must be used by the example</span>
<span
class="cmr-12">program.</span>
<!--l. 126--><p class="indent" > <span
class="cmr-12">The part of the code dealing with reading and assembling the sparse matrix and the</span>
<span
class="cmr-12">right-hand side vector and the deallocation of the relevant data structures, performed</span>
<span
class="cmr-12">through the PSBLAS routines for sparse matrix and vector management,</span>
<span
class="cmr-12">is not reported here for the sake of conciseness. The complete code can be</span>
<span
class="cmr-12">found 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="cmr-12">, in the directory</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">samples/simple/fileread</span></span></span> <span
class="cmr-12">of the AMG4PSBLAS implementation (see Section</span><span
class="cmr-12">&#x00A0;</span><a
href="userhtmlsu5.html#x13-120003.5"><span
class="cmr-12">3.5</span><!--tex4ht:ref: sec:ex_and_test --></a><span
class="cmr-12">). A</span>
<span
class="cmr-12">sample test problem along with the relevant input data is available in</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">samples/simple/fileread/runs</span></span></span><span
class="cmr-12">. For details on the use of the PSBLAS routines, see</span>
<span
class="cmr-12">the PSBLAS User&#8217;s Guide</span><span
class="cmr-12">&#x00A0;</span><span class="cite"><span
class="cmr-12">[</span><a
href="userhtmlli5.html#XPSBLASGUIDE"><span
class="cmr-12">20</span></a><span
class="cmr-12">]</span></span><span
class="cmr-12">.</span>
<!--l. 138--><p class="indent" > <span
class="cmr-12">The setup and application of the default multilevel preconditioner for the real single</span>
<span
class="cmr-12">precision and the complex, single and double precision, versions are obtained</span>
<span
class="cmr-12">with straightforward modifications of the previous example (see Section</span><span
class="cmr-12">&#x00A0;</span><a
href="userhtmlse5.html#x17-160005"><span
class="cmr-12">5</span><!--tex4ht:ref: sec:userinterface --></a> <span
class="cmr-12">for</span>
<span
class="cmr-12">details). If these versions are installed, the corresponding codes are available in</span>
<span class="obeylines-h"><span class="verb"><span
class="cmtt-12">samples/simple/fileread/</span></span></span><span
class="cmr-12">.</span>
<!--l. 144--><p class="indent" > <a
id="x15-14001r1"></a><hr class="float"><div class="float"
>
<div class="center"
>
<!--l. 145--><p class="noindent" >
<div class="minipage"><pre class="verbatim" id="verbatim-6">
&#x00A0;&#x00A0;use&#x00A0;psb_base_mod
&#x00A0;&#x00A0;use&#x00A0;amg_prec_mod
&#x00A0;&#x00A0;use&#x00A0;psb_krylov_mod
...&#x00A0;...
!
!&#x00A0;sparse&#x00A0;matrix
&#x00A0;&#x00A0;type(psb_dspmat_type)&#x00A0;::&#x00A0;A
!&#x00A0;sparse&#x00A0;matrix&#x00A0;descriptor
&#x00A0;&#x00A0;type(psb_desc_type)&#x00A0;&#x00A0;&#x00A0;::&#x00A0;desc_A
!&#x00A0;preconditioner
&#x00A0;&#x00A0;type(amg_dprec_type)&#x00A0;&#x00A0;::&#x00A0;P
!&#x00A0;right-hand&#x00A0;side&#x00A0;and&#x00A0;solution&#x00A0;vectors
&#x00A0;&#x00A0;type(psb_d_vect_type)&#x00A0;::&#x00A0;b,&#x00A0;x
...&#x00A0;...
!
!&#x00A0;initialize&#x00A0;the&#x00A0;parallel&#x00A0;environment
&#x00A0;&#x00A0;call&#x00A0;psb_init(ctxt)
&#x00A0;&#x00A0;call&#x00A0;psb_info(ctxt,iam,np)
...&#x00A0;...
!
!&#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;using&#x00A0;PSBLAS&#x00A0;routines&#x00A0;for&#x00A0;sparse&#x00A0;matrix&#x00A0;/&#x00A0;vector&#x00A0;management
...&#x00A0;...
!
!&#x00A0;initialize&#x00A0;the&#x00A0;default&#x00A0;multilevel&#x00A0;preconditioner,&#x00A0;i.e.&#x00A0;V-cycle
!&#x00A0;with&#x00A0;basic&#x00A0;smoothed&#x00A0;aggregation,&#x00A0;1&#x00A0;hybrid&#x00A0;forward/backward
!&#x00A0;GS&#x00A0;sweep&#x00A0;as&#x00A0;pre/post-smoother&#x00A0;and&#x00A0;UMFPACK&#x00A0;as&#x00A0;coarsest-level
!&#x00A0;solver
&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
!
!&#x00A0;build&#x00A0;the&#x00A0;preconditioner
&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
!
!&#x00A0;set&#x00A0;the&#x00A0;solver&#x00A0;parameters&#x00A0;and&#x00A0;the&#x00A0;initial&#x00A0;guess
&#x00A0;&#x00A0;...&#x00A0;...
!
!&#x00A0;solve&#x00A0;Ax=b&#x00A0;with&#x00A0;preconditioned&#x00A0;CG
&#x00A0;&#x00A0;call&#x00A0;psb_krylov(&#8217;CG&#8217;,A,P,b,x,tol,desc_A,info)
&#x00A0;&#x00A0;...&#x00A0;...
!
!&#x00A0;deallocate&#x00A0;the&#x00A0;preconditioner
&#x00A0;&#x00A0;call&#x00A0;P%free(info)
!
!&#x00A0;deallocate&#x00A0;other&#x00A0;data&#x00A0;structures
&#x00A0;&#x00A0;...&#x00A0;...
!
!&#x00A0;exit&#x00A0;the&#x00A0;parallel&#x00A0;environment
&#x00A0;&#x00A0;call&#x00A0;psb_exit(ctxt)
&#x00A0;&#x00A0;stop
</pre>
<!--l. 255--><p class="nopar" >
</div>
<br /> <div class="caption"
><span class="id">Listing 1: </span><span
class="content">setup and application of the default multilevel preconditioner (example 1).
</span></div><!--tex4ht:label?: x15-14001r1 -->
</div>
</div><hr class="endfloat" />
<!--l. 264--><p class="indent" > <span
class="cmr-12">Different versions of the multilevel preconditioner can be obtained by changing the</span>
<span
class="cmr-12">default values of the preconditioner parameters. The code reported in Figure</span><span
class="cmr-12">&#x00A0;</span><a
href="#x15-14002r2"><span
class="cmr-12">2</span><!--tex4ht:ref: fig:ex2 --></a> <span
class="cmr-12">shows</span>
<span
class="cmr-12">how to set a V-cycle preconditioner which applies 1 block-Jacobi sweep as pre-</span>
<span
class="cmr-12">and post-smoother, and solves the coarsest-level system with 8 block-Jacobi</span>
<span
class="cmr-12">sweeps. Note that the ILU(0) factorization (plus triangular solve) is used as</span>
<span
class="cmr-12">local solver for the block-Jacobi sweeps, since this is the default associated</span>
<span
class="cmr-12">with block-Jacobi and set by</span><span
class="cmr-12">&#x00A0;</span><span class="lstinline"></span><span
class="cmtt-12">P</span><span
class="cmtt-12">%</span><span
class="cmtt-12">init</span><span
class="cmr-12">. Furthermore, specifying block-Jacobi as</span>
<span
class="cmr-12">coarsest-level solver implies that the coarsest-level matrix is distributed among</span>
<span
class="cmr-12">the processes. Figure</span><span
class="cmr-12">&#x00A0;</span><a
href="#x15-14003r3"><span
class="cmr-12">3</span><!--tex4ht:ref: fig:ex3 --></a> <span
class="cmr-12">shows how to set a W-cycle preconditioner using the</span>
<span
class="cmr-12">Coarsening based on Compatible Weighted Matching, aggregates of size at</span>
<span
class="cmr-12">most 8 and smoothed prolongators. It applies 2 hybrid Gauss-Seidel sweeps as</span>
<span
class="cmr-12">pre- and post-smoother, and solves the coarsest-level system with the parallel</span>
<span
class="cmr-12">flexible Conjugate Gradient method (KRM) coupled with the block-Jacobi</span>
<span
class="cmr-12">preconditioner having ILU(0) on the blocks. Default parameters are used for stopping</span>
<span
class="cmr-12">criterion of the coarsest solver. Note that, also in this case, specifying KRM as</span>
<span
class="cmr-12">coarsest-level solver implies that the coarsest-level matrix is distributed among the</span>
<span
class="cmr-12">processes.</span>
<!--l. 291--><p class="indent" > <span
class="cmr-12">The code fragments shown in Figures</span><span
class="cmr-12">&#x00A0;</span><a
href="#x15-14002r2"><span
class="cmr-12">2</span><!--tex4ht:ref: fig:ex2 --></a> <span
class="cmr-12">and </span><a
href="#x15-14003r3"><span
class="cmr-12">3</span><!--tex4ht:ref: fig:ex3 --></a> <span
class="cmr-12">are included in the example program</span>
<span
class="cmr-12">file </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">amg_dexample_ml.f90</span></span></span> <span
class="cmr-12">too.</span>
<!--l. 294--><p class="indent" > <span
class="cmr-12">Finally, Figure</span><span
class="cmr-12">&#x00A0;</span><a
href="#x15-14004r4"><span
class="cmr-12">4</span><!--tex4ht:ref: fig:ex4 --></a> <span
class="cmr-12">shows the setup of a one-level additive Schwarz preconditioner,</span>
<span
class="cmr-12">i.e., RAS with overlap 2. Note also that a Krylov method different from CG</span>
<span
class="cmr-12">must be used to solve the preconditioned system, since the preconditione in</span>
<span
class="cmr-12">nonsymmetric. The corresponding example program is 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="cmr-12">.</span>
<!--l. 301--><p class="indent" > <span
class="cmr-12">For all the previous preconditioners, example programs where the sparse matrix</span>
<span
class="cmr-12">and the right-hand side are generated by discretizing a PDE with Dirichlet</span>
<span
class="cmr-12">boundary conditions are also available in the directory </span><span class="obeylines-h"><span class="verb"><span
class="cmtt-12">samples/simple/pdegen</span></span></span><span
class="cmr-12">.</span>
<!--l. 304--><p class="indent" > <a
id="x15-14002r2"></a><hr class="float"><div class="float"
>
<div class="center"
>
<!--l. 318--><p class="noindent" >
<div class="minipage"><pre class="verbatim" id="verbatim-7">
...&#x00A0;...
!&#x00A0;build&#x00A0;a&#x00A0;V-cycle&#x00A0;preconditioner&#x00A0;with&#x00A0;1&#x00A0;block-Jacobi&#x00A0;sweep&#x00A0;(with
!&#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;sweeps&#x00A0;(with&#x00A0;ILU(0)&#x00A0;on&#x00A0;the&#x00A0;blocks)&#x00A0;as&#x00A0;coarsest-level&#x00A0;solver
&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SMOOTHER_TYPE&#8217;,&#8217;BJAC&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SOLVE&#8217;,&#8217;BJAC&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SWEEPS&#8217;,8,info)
&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
...&#x00A0;...
</pre>
<!--l. 333--><p class="nopar" ></div></div>
<br /><div class="caption"
><span class="id">Listing 2: </span><span
class="content">setup of a multilevel preconditioner based on the default decoupled coarsening</span></div><!--tex4ht:label?: x15-14002r2 -->
</div><hr class="endfloat" />
<!--l. 340--><p class="indent" > <a
id="x15-14003r3"></a><hr class="float"><div class="float"
>
<div class="center"
>
<!--l. 360--><p class="noindent" >
<div class="minipage"><pre class="verbatim" id="verbatim-8">
...&#x00A0;...
!&#x00A0;build&#x00A0;a&#x00A0;W-cycle&#x00A0;preconditioner&#x00A0;with&#x00A0;2&#x00A0;hybrid&#x00A0;Gauss-Seidel&#x00A0;sweeps
!&#x00A0;as&#x00A0;pre-&#x00A0;and&#x00A0;post-smoother,&#x00A0;a&#x00A0;distributed&#x00A0;coarsest
!&#x00A0;matrix,&#x00A0;and&#x00A0;MUMPS&#x00A0;as&#x00A0;coarsest-level&#x00A0;solver
&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;ML&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;PAR_AGGR_ALG&#8217;,&#8217;COUPLED&#8217;,info)
call&#x00A0;P%set(&#8217;AGGR_TYPE&#8217;,&#8217;MATCHBOXP&#8217;,info)
call&#x00A0;P%set(&#8217;AGGR_SIZE&#8217;,8,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;ML_CYCLE&#8217;,&#8217;WCYCLE&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SMOOTHER_TYPE&#8217;,&#8217;FBGS&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SMOOTHER_SWEEPS&#8217;,2,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;COARSE_SOLVE&#8217;,&#8217;KRM&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%hierarchy_build(A,desc_A,info)
&#x00A0;&#x00A0;call&#x00A0;P%smoothers_build(A,desc_A,info)
...&#x00A0;...
</pre>
<!--l. 379--><p class="nopar" ></div></div>
<br /> <div class="caption"
><span class="id">Listing 3: </span><span
class="content">setup of a multilevel preconditioner based on the coupled coarsening using
weighted matching</span></div><!--tex4ht:label?: x15-14003r3 -->
</div><hr class="endfloat" />
<!--l. 386--><p class="indent" > <a
id="x15-14004r4"></a><hr class="float"><div class="float"
>
<div class="center"
>
<!--l. 398--><p class="noindent" >
<div class="minipage"><pre class="verbatim" id="verbatim-9">
...&#x00A0;...
!&#x00A0;set&#x00A0;RAS&#x00A0;with&#x00A0;overlap&#x00A0;2&#x00A0;and&#x00A0;ILU(0)&#x00A0;on&#x00A0;the&#x00A0;local&#x00A0;blocks
&#x00A0;&#x00A0;call&#x00A0;P%init(&#8217;AS&#8217;,info)
&#x00A0;&#x00A0;call&#x00A0;P%set(&#8217;SUB_OVR&#8217;,2,info)
&#x00A0;&#x00A0;call&#x00A0;P%bld(A,desc_A,info)
...&#x00A0;...
!&#x00A0;solve&#x00A0;Ax=b&#x00A0;with&#x00A0;preconditioned&#x00A0;BiCGSTAB
&#x00A0;&#x00A0;call&#x00A0;psb_krylov(&#8217;BICGSTAB&#8217;,A,P,b,x,tol,desc_A,info)
</pre>
<!--l. 410--><p class="nopar" ></div></div>
<br /> <div class="caption"
><span class="id">Listing 4: </span><span
class="content">setup of a one-level Schwarz preconditioner.</span></div><!--tex4ht:label?: x15-14004r4 -->
</div><hr class="endfloat" />
<!--l. 420--><div class="crosslinks"><p class="noindent"><span
class="cmr-12">[</span><a
href="userhtmlsu7.html" ><span
class="cmr-12">next</span></a><span
class="cmr-12">] [</span><a
href="userhtmlsu6.html" ><span
class="cmr-12">front</span></a><span
class="cmr-12">] [</span><a
href="userhtmlse4.html#userhtmlsu6.html" ><span
class="cmr-12">up</span></a><span
class="cmr-12">] </span></p></div>
<!--l. 420--><p class="indent" > <a
id="tailuserhtmlsu6.html"></a>
</body></html>