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.
316 lines
11 KiB
HTML
316 lines
11 KiB
HTML
<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 4.0 Transitional//EN">
|
|
|
|
<!--Converted with LaTeX2HTML 2017.2 (Released Jan 23, 2017) -->
|
|
<HTML>
|
|
<HEAD>
|
|
<TITLE>Examples</TITLE>
|
|
<META NAME="description" CONTENT="Examples">
|
|
<META NAME="keywords" CONTENT="userhtml">
|
|
<META NAME="resource-type" CONTENT="document">
|
|
<META NAME="distribution" CONTENT="global">
|
|
|
|
<META NAME="Generator" CONTENT="LaTeX2HTML v2017.2">
|
|
<META HTTP-EQUIV="Content-Style-Type" CONTENT="text/css">
|
|
|
|
<LINK REL="STYLESHEET" HREF="userhtml.css">
|
|
|
|
<LINK REL="previous" HREF="node16.html">
|
|
<LINK REL="up" HREF="node16.html">
|
|
<LINK REL="next" HREF="node18.html">
|
|
</HEAD>
|
|
|
|
<BODY >
|
|
|
|
<DIV CLASS="navigation"><!--Navigation Panel-->
|
|
<A NAME="tex2html294"
|
|
HREF="node18.html">
|
|
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A>
|
|
<A NAME="tex2html290"
|
|
HREF="node16.html">
|
|
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A>
|
|
<A NAME="tex2html286"
|
|
HREF="node16.html">
|
|
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A>
|
|
<A NAME="tex2html292"
|
|
HREF="node2.html">
|
|
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>
|
|
<BR>
|
|
<B> Next:</B> <A NAME="tex2html295"
|
|
HREF="node18.html">User Interface</A>
|
|
<B> Up:</B> <A NAME="tex2html291"
|
|
HREF="node16.html">Getting Started</A>
|
|
<B> Previous:</B> <A NAME="tex2html287"
|
|
HREF="node16.html">Getting Started</A>
|
|
<B> <A NAME="tex2html293"
|
|
HREF="node2.html">Contents</A></B>
|
|
<BR>
|
|
<BR></DIV>
|
|
<!--End of Navigation Panel-->
|
|
|
|
<H2><A NAME="SECTION00071000000000000000"></A><A NAME="sec:examples"></A>
|
|
<BR>
|
|
Examples
|
|
</H2><BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">The code reported in Figure <A HREF="#fig:ex1">2</A> shows how to set and apply the default
|
|
multilevel preconditioner available in the real double precision version
|
|
of MLD2P4 (see Table <A HREF="#tab:precinit">1</A>). This preconditioner is chosen
|
|
by simply specifying <code>'ML'</code> as the second argument of <code>P%init</code>
|
|
(a call to <code>P%set</code> is not needed) and is applied with the CG
|
|
solver provided by PSBLAS (the matrix of the system to be solved is
|
|
assumed to be positive definite). As previously observed, the modules
|
|
<code>psb_base_mod</code>, <code>mld_prec_mod</code> and <code>psb_krylov_mod</code>
|
|
must be used by the example program.
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">The part of the code concerning the
|
|
reading and assembling of the sparse matrix and the right-hand side vector, performed
|
|
through the PSBLAS routines for sparse matrix and vector management, is not reported
|
|
here for brevity; the statements concerning the deallocation of the PSBLAS
|
|
data structure are neglected too.
|
|
The complete code can be found in the example program file <code>mld_dexample_ml.f90</code>,
|
|
in the directory <code>examples/fileread</code> of the MLD2P4 implementation (see
|
|
Section <A HREF="node11.html#sec:ex_and_test">3.5</A>). A sample test problem along with the relevant
|
|
input data is available in <code>examples/fileread/runs</code>.
|
|
For details on the use of the PSBLAS routines, see the PSBLAS User's
|
|
Guide [<A
|
|
HREF="node36.html#PSBLASGUIDE">13</A>].
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">The setup and application of the default multilevel preconditioner
|
|
for the real single precision and the complex, single and double
|
|
precision, versions are obtained with straightforward modifications of the previous
|
|
example (see Section <A HREF="node18.html#sec:userinterface">6</A> for details). If these versions are installed,
|
|
the corresponding codes are available in <code>examples/fileread/</code>.
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<DIV ALIGN="CENTER"><A NAME="fig:ex1"></A><A NAME="907"></A>
|
|
<TABLE>
|
|
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 2:</STRONG>
|
|
setup and application of the default multilevel preconditioner (example 1).
|
|
</CAPTION>
|
|
<TR><TD>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV><TABLE WIDTH="90%">
|
|
<TR><TD>
|
|
<PRE>
|
|
use psb_base_mod
|
|
use mld_prec_mod
|
|
use psb_krylov_mod
|
|
... ...
|
|
!
|
|
! sparse matrix
|
|
type(psb_dspmat_type) :: A
|
|
! sparse matrix descriptor
|
|
type(psb_desc_type) :: desc_A
|
|
! preconditioner
|
|
type(mld_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>
|
|
</TD></TR>
|
|
</TABLE>
|
|
<DIV ALIGN="CENTER">
|
|
|
|
</DIV></TD></TR>
|
|
</TABLE>
|
|
</DIV>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">Different versions of the multilevel preconditioner can be obtained by changing
|
|
the default values of the preconditioner parameters. The code reported in
|
|
Figure <A HREF="#fig:ex2">3</A> shows how to set a V-cycle preconditioner
|
|
which applies 1 block-Jacobi sweep as pre- and post-smoother,
|
|
and solves the coarsest-level system with 8 block-Jacobi sweeps.
|
|
Note that the ILU(0) factorization (plus triangular solve) is used as
|
|
local solver for the block-Jacobi sweeps, since this is the default associated
|
|
with block-Jacobi and set by <code>P%init</code>.
|
|
Furthermore, specifying block-Jacobi as coarsest-level
|
|
solver implies that the coarsest-level matrix is distributed
|
|
among the processes.
|
|
Figure <A HREF="#fig:ex3">4</A> shows how to set a W-cycle preconditioner which
|
|
applies 2 hybrid Gauss-Seidel sweeps as pre- and post-smoother,
|
|
and solves the coarsest-level system with the multifrontal LU factorization
|
|
implemented in MUMPS. It is specified that the coarsest-level
|
|
matrix is distributed, since MUMPS can be used on both
|
|
replicated and distributed matrices, and by default
|
|
it is used on replicated ones.
|
|
The code fragments shown in Figures <A HREF="#fig:ex2">3</A> and <A HREF="#fig:ex3">4</A> are
|
|
included in the example program file <code>mld_dexample_ml.f90</code> too.
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">Finally, Figure <A HREF="#fig:ex4">5</A> shows the setup of a one-level
|
|
additive Schwarz preconditioner, i.e., RAS with overlap 2.
|
|
Note also that a Krylov method different from CG must be used to solve
|
|
the preconditioned system, since the preconditione in nonsymmetric.
|
|
The corresponding example program is available in the file
|
|
<code>mld_dexample_1lev.f90</code>.
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE">For all the previous preconditioners, example programs where the sparse matrix and
|
|
the right-hand side are generated by discretizing a PDE with Dirichlet
|
|
boundary conditions are also available in the directory <code>examples/pdegen</code>.
|
|
</BIG></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<DIV ALIGN="CENTER"><A NAME="fig:ex2"></A><A NAME="909"></A>
|
|
<TABLE>
|
|
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 3:</STRONG>
|
|
setup of a multilevel preconditioner</CAPTION>
|
|
<TR><TD>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV><TABLE WIDTH="90%">
|
|
<TR><TD>
|
|
<PRE>
|
|
... ...
|
|
! 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>
|
|
</TD></TR>
|
|
</TABLE>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV>
|
|
<P>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV></TD></TR>
|
|
</TABLE>
|
|
</DIV>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<DIV ALIGN="CENTER"><A NAME="fig:ex3"></A><A NAME="911"></A>
|
|
<TABLE>
|
|
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 4:</STRONG>
|
|
setup of a multilevel preconditioner</CAPTION>
|
|
<TR><TD>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV><TABLE WIDTH="90%">
|
|
<TR><TD>
|
|
<PRE>
|
|
... ...
|
|
! 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>
|
|
</TD></TR>
|
|
</TABLE>
|
|
<DIV ALIGN="CENTER">
|
|
|
|
</DIV></TD></TR>
|
|
</TABLE>
|
|
</DIV>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<DIV ALIGN="CENTER"><A NAME="fig:ex4"></A><A NAME="913"></A>
|
|
<TABLE>
|
|
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 5:</STRONG>
|
|
setup of a one-level Schwarz preconditioner.</CAPTION>
|
|
<TR><TD>
|
|
<DIV ALIGN="CENTER">
|
|
</DIV><TABLE WIDTH="90%">
|
|
<TR><TD>
|
|
<PRE>
|
|
... ...
|
|
! 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>
|
|
</TD></TR>
|
|
</TABLE>
|
|
<DIV ALIGN="CENTER">
|
|
|
|
</DIV></TD></TR>
|
|
</TABLE>
|
|
</DIV>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG>
|
|
<P>
|
|
<BIG CLASS="LARGE"><BIG CLASS="LARGE"><BIG CLASS="LARGE"></BIG></BIG></BIG>
|
|
<DIV CLASS="navigation"><HR>
|
|
<!--Navigation Panel-->
|
|
<A NAME="tex2html294"
|
|
HREF="node18.html">
|
|
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A>
|
|
<A NAME="tex2html290"
|
|
HREF="node16.html">
|
|
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A>
|
|
<A NAME="tex2html286"
|
|
HREF="node16.html">
|
|
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A>
|
|
<A NAME="tex2html292"
|
|
HREF="node2.html">
|
|
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>
|
|
<BR>
|
|
<B> Next:</B> <A NAME="tex2html295"
|
|
HREF="node18.html">User Interface</A>
|
|
<B> Up:</B> <A NAME="tex2html291"
|
|
HREF="node16.html">Getting Started</A>
|
|
<B> Previous:</B> <A NAME="tex2html287"
|
|
HREF="node16.html">Getting Started</A>
|
|
<B> <A NAME="tex2html293"
|
|
HREF="node2.html">Contents</A></B> </DIV>
|
|
<!--End of Navigation Panel-->
|
|
|
|
</BODY>
|
|
</HTML>
|