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.
384 lines
18 KiB
HTML
384 lines
18 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. 109--><div class="crosslinks"><p class="noindent"><span
|
|
class="cmr-12">[</span><a
|
|
href="#tailuserhtmlsu9.html"><span
|
|
class="cmr-12">tail</span></a><span
|
|
class="cmr-12">] [</span><a
|
|
href="userhtmlse5.html#userhtmlsu9.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">5.1 </span></span> <a
|
|
id="x18-170005.1"></a><span
|
|
class="cmr-12">Examples</span></h4>
|
|
<!--l. 111--><p class="noindent" ><span
|
|
class="cmr-12">The code reported in Figure</span><span
|
|
class="cmr-12"> </span><a
|
|
href="#x18-170012"><span
|
|
class="cmr-12">2</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"> </span><a
|
|
href="userhtmlse5.html#x17-160151"><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="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">’ML’</span></span></span> <span
|
|
class="cmr-12">as the</span>
|
|
<span
|
|
class="cmr-12">second argument of </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">P%init</span></span></span> <span
|
|
class="cmr-12">(a call to </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">P%set</span></span></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="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">psb_base_mod</span></span></span><span
|
|
class="cmr-12">, </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">amg_prec_mod</span></span></span> <span
|
|
class="cmr-12">and </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">psb_krylov_mod</span></span></span> <span
|
|
class="cmr-12">must be used by the example</span>
|
|
<span
|
|
class="cmr-12">program.</span>
|
|
<!--l. 121--><p class="indent" > <span
|
|
class="cmr-12">The part of the code concerning the reading and assembling of the sparse matrix</span>
|
|
<span
|
|
class="cmr-12">and the right-hand side vector, performed through the PSBLAS routines for sparse</span>
|
|
<span
|
|
class="cmr-12">matrix and vector management, is not reported here for brevity; the statements</span>
|
|
<span
|
|
class="cmr-12">concerning the deallocation of the PSBLAS data structure are neglected too. The</span>
|
|
<span
|
|
class="cmr-12">complete code can be 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</span>
|
|
<span
|
|
class="cmr-12">directory </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">examples/fileread</span></span></span> <span
|
|
class="cmr-12">of the AMG4PSBLAS implementation (see</span>
|
|
<span
|
|
class="cmr-12">Section</span><span
|
|
class="cmr-12"> </span><a
|
|
href="userhtmlsu5.html#x12-110003.5"><span
|
|
class="cmr-12">3.5</span><!--tex4ht:ref: sec:ex_and_test --></a><span
|
|
class="cmr-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="cmr-12">. For details on the use of the PSBLAS routines, see the</span>
|
|
<span
|
|
class="cmr-12">PSBLAS User’s Guide</span><span
|
|
class="cmr-12"> </span><span class="cite"><span
|
|
class="cmr-12">[</span><a
|
|
href="userhtmlli4.html#XPSBLASGUIDE"><span
|
|
class="cmr-12">15</span></a><span
|
|
class="cmr-12">]</span></span><span
|
|
class="cmr-12">.</span>
|
|
|
|
|
|
|
|
<!--l. 133--><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"> </span><a
|
|
href="userhtmlse6.html#x19-180006"><span
|
|
class="cmr-12">6</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">examples/fileread/</span></span></span><span
|
|
class="cmr-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"><pre class="verbatim" id="verbatim-6">
|
|
  use psb_base_mod
|
|
  use amg_prec_mod
|
|
  use psb_krylov_mod
|
|
... ...
|
|
!
|
|
! sparse matrix
|
|
  type(psb_dspmat_type) :: A
|
|
! sparse matrix descriptor
|
|
  type(psb_desc_type)   :: desc_A
|
|
! preconditioner
|
|
  type(amg_dprec_type)  :: P
|
|
! right-hand side and solution vectors
|
|
  type(psb_d_vect_type) :: b, x
|
|
... ...
|
|
!
|
|
! initialize the parallel environment
|
|
  call psb_init(ictxt)
|
|
  call psb_info(ictxt,iam,np)
|
|
... ...
|
|
!
|
|
! read and assemble the spd matrix A and the right-hand side b
|
|
! using PSBLAS routines for sparse matrix / vector management
|
|
... ...
|
|
!
|
|
! initialize the default multilevel preconditioner, i.e. V-cycle
|
|
! with basic smoothed aggregation, 1 hybrid forward/backward
|
|
! GS sweep as pre/post-smoother and UMFPACK as coarsest-level
|
|
! solver
|
|
  call P%init(’ML’,info)
|
|
!
|
|
! build the preconditioner
|
|
  call P%hierarchy_build(A,desc_A,info)
|
|
  call P%smoothers_build(A,desc_A,info)
|
|
!
|
|
! set the solver parameters and the initial guess
|
|
  ... ...
|
|
!
|
|
! solve Ax=b with preconditioned CG
|
|
  call psb_krylov(’CG’,A,P,b,x,tol,desc_A,info)
|
|
  ... ...
|
|
!
|
|
! deallocate the preconditioner
|
|
  call P%free(info)
|
|
!
|
|
! deallocate other data structures
|
|
  ... ...
|
|
!
|
|
! exit the parallel environment
|
|
  call psb_exit(ictxt)
|
|
  stop
|
|
</pre>
|
|
<!--l. 195--><p class="nopar" >
|
|
|
|
|
|
</div>
|
|
<br /> <div class="caption"
|
|
><span class="id">Figure 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="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"> </span><a
|
|
href="#x18-170023"><span
|
|
class="cmr-12">3</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- and</span>
|
|
<span
|
|
class="cmr-12">post-smoother, and solves the coarsest-level system with 8 block-Jacobi sweeps. Note</span>
|
|
<span
|
|
class="cmr-12">that the ILU(0) factorization (plus triangular solve) is used as local solver for the</span>
|
|
<span
|
|
class="cmr-12">block-Jacobi sweeps, since this is the default associated with block-Jacobi and set</span>
|
|
<span
|
|
class="cmr-12">by</span><span
|
|
class="cmr-12"> </span><span class="obeylines-h"><span class="verb"><span
|
|
class="cmtt-12">P%init</span></span></span><span
|
|
class="cmr-12">. Furthermore, specifying block-Jacobi as coarsest-level solver implies that</span>
|
|
<span
|
|
class="cmr-12">the coarsest-level matrix is distributed among the processes. Figure</span><span
|
|
class="cmr-12"> </span><a
|
|
href="#x18-170034"><span
|
|
class="cmr-12">4</span><!--tex4ht:ref: fig:ex3 --></a> <span
|
|
class="cmr-12">shows how to set</span>
|
|
<span
|
|
class="cmr-12">a W-cycle preconditioner which applies 2 hybrid Gauss-Seidel sweeps as pre- and</span>
|
|
<span
|
|
class="cmr-12">post-smoother, and solves the coarsest-level system with the multifrontal LU</span>
|
|
<span
|
|
class="cmr-12">factorization implemented in MUMPS. It is specified that the coarsest-level matrix is</span>
|
|
<span
|
|
class="cmr-12">distributed, since MUMPS can be used on both replicated and distributed matrices,</span>
|
|
<span
|
|
class="cmr-12">and by default it is used on replicated ones. The code fragments shown in</span>
|
|
<span
|
|
class="cmr-12">Figures</span><span
|
|
class="cmr-12"> </span><a
|
|
href="#x18-170023"><span
|
|
class="cmr-12">3</span><!--tex4ht:ref: fig:ex2 --></a> <span
|
|
class="cmr-12">and </span><a
|
|
href="#x18-170034"><span
|
|
class="cmr-12">4</span><!--tex4ht:ref: fig:ex3 --></a> <span
|
|
class="cmr-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="cmr-12">too.</span>
|
|
<!--l. 227--><p class="indent" > <span
|
|
class="cmr-12">Finally, Figure</span><span
|
|
class="cmr-12"> </span><a
|
|
href="#x18-170045"><span
|
|
class="cmr-12">5</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. 234--><p class="indent" > <span
|
|
class="cmr-12">For all the previous preconditioners, example programs where the sparse matrix and</span>
|
|
<span
|
|
class="cmr-12">the right-hand side are generated by discretizing a PDE with Dirichlet boundary</span>
|
|
<span
|
|
class="cmr-12">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="cmr-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"><pre class="verbatim" id="verbatim-7">
|
|
... ...
|
|
! build a V-cycle preconditioner with 1 block-Jacobi sweep (with
|
|
! ILU(0) on the blocks) as pre- and post-smoother, and 8  block-Jacobi
|
|
! sweeps (with ILU(0) on the blocks) as coarsest-level solver
|
|
  call P%init(’ML’,info)
|
|
  call_P%set(’SMOOTHER_TYPE’,’BJAC’,info)
|
|
  call P%set(’COARSE_SOLVE’,’BJAC’,info)
|
|
  call P%set(’COARSE_SWEEPS’,8,info)
|
|
  call P%hierarchy_build(A,desc_A,info)
|
|
  call P%smoothers_build(A,desc_A,info)
|
|
... ...
|
|
</pre>
|
|
<!--l. 254--><p class="nopar" ></div>
|
|
<br /> <div class="caption"
|
|
><span class="id">Figure 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"><pre class="verbatim" id="verbatim-8">
|
|
... ...
|
|
! build a W-cycle preconditioner with 2 hybrid Gauss-Seidel sweeps
|
|
! as pre- and post-smoother, a distributed coarsest
|
|
! matrix, and MUMPS as coarsest-level solver
|
|
  call P%init(’ML’,info)
|
|
  call P%set(’ML_CYCLE’,’WCYCLE’,info)
|
|
  call P%set(’SMOOTHER_TYPE’,’FBGS’,info)
|
|
  call P%set(’SMOOTHER_SWEEPS’,2,info)
|
|
  call P%set(’COARSE_SOLVE’,’MUMPS’,info)
|
|
  call P%set(’COARSE_MAT’,’DIST’,info)
|
|
  call P%hierarchy_build(A,desc_A,info)
|
|
  call P%smoothers_build(A,desc_A,info)
|
|
... ...
|
|
</pre>
|
|
<!--l. 280--><p class="nopar" ></div>
|
|
<br /> <div class="caption"
|
|
><span class="id">Figure 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"><pre class="verbatim" id="verbatim-9">
|
|
... ...
|
|
! set RAS with overlap 2 and ILU(0) on the local blocks
|
|
  call P%init(’AS’,info)
|
|
  call P%set(’SUB_OVR’,2,info)
|
|
  call P%bld(A,desc_A,info)
|
|
... ...
|
|
! solve Ax=b with preconditioned BiCGSTAB
|
|
  call psb_krylov(’BICGSTAB’,A,P,b,x,tol,desc_A,info)
|
|
</pre>
|
|
<!--l. 300--><p class="nopar" ></div>
|
|
<br /> <div class="caption"
|
|
><span class="id">Figure 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="cmr-12">[</span><a
|
|
href="userhtmlsu9.html" ><span
|
|
class="cmr-12">front</span></a><span
|
|
class="cmr-12">] [</span><a
|
|
href="userhtmlse5.html#userhtmlsu9.html" ><span
|
|
class="cmr-12">up</span></a><span
|
|
class="cmr-12">] </span></p></div>
|
|
<!--l. 1--><p class="indent" > <a
|
|
id="tailuserhtmlsu9.html"></a>
|
|
</body></html>
|