<!DOCTYPE HTML PUBLIC "-//W3C//DTD HTML 3.2 Final//EN">

<!--Converted with LaTeX2HTML 2008 (1.71)
original version by:  Nikos Drakos, CBLU, University of Leeds
* revised and updated by:  Marcus Hennecke, Ross Moore, Herb Swan
* with significant contributions from:
  Jens Lippmann, Marek Rouchal, Martin Wilck and others -->
<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 v2008">
<META HTTP-EQUIV="Content-Style-Type" CONTENT="text/css">

<LINK REL="STYLESHEET" HREF="userhtml.css">

<LINK REL="previous" HREF="node14.html">
<LINK REL="up" HREF="node14.html">
<LINK REL="next" HREF="node16.html">
</HEAD>

<BODY >
<!--Navigation Panel-->
<A NAME="tex2html242"
  HREF="node16.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html238"
  HREF="node14.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html234"
  HREF="node14.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html240"
  HREF="node2.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>  
<BR>
<B> Next:</B> <A NAME="tex2html243"
  HREF="node16.html">User Interface</A>
<B> Up:</B> <A NAME="tex2html239"
  HREF="node14.html">Getting Started</A>
<B> Previous:</B> <A NAME="tex2html235"
  HREF="node14.html">Getting Started</A>
 &nbsp; <B>  <A NAME="tex2html241"
  HREF="node2.html">Contents</A></B> 
<BR>
<BR>
<!--End of Navigation Panel-->

<H2><A NAME="SECTION00071000000000000000"></A><A NAME="sec:examples"></A>
<BR>
Examples
</H2>

<P>
The code reported in Figure&nbsp;<A HREF="#fig:ex_default">2</A> shows how to set and apply the default
multi-level preconditioner available in the real double precision version
of MLD2P4 (see Table&nbsp;<A HREF="#tab:precinit">1</A>). This preconditioner is chosen
by simply specifying <code>'ML'</code> as second argument of <code>mld_precinit</code>
(a call to <code>mld_precset</code> is not needed) and is applied with the BiCGSTAB
solver provided by PSBLAS. 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.

<P>
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 tree (see
Section&nbsp;<A HREF="node10.html#sec:ex_and_test">3.5</A>).
For details on the use of the PSBLAS routines, see the PSBLAS User's
Guide [<A
 HREF="node25.html#PSBLASGUIDE">15</A>].

<P>
The setup and application of the default multi-level
preconditioners for the real single precision and the complex, single and double
precision, versions are obtained with straightforward modifications of the previous
example (see Section&nbsp;<A HREF="node16.html#sec:userinterface">6</A> for details). If these versions are installed,
the corresponding Fortran 95 codes are available in <code>examples/fileread/</code>.

<P>

<DIV ALIGN="CENTER"><A NAME="fig:ex_default"></A><A NAME="932"></A>
<TABLE>
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 2:</STRONG>
Setup and application of the default multi-level Schwarz preconditioner.
</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
  real(kind(1.d0))      :: b(:), x(:)
... ...
!
! initialize the parallel environment
  call psb_init(ictxt)
  call psb_info(ictxt,iam,np)
... ...
!
! read and assemble the matrix A and the right-hand side b 
! using PSBLAS routines for sparse matrix / vector management 
... ...
!
! initialize the default multi-level preconditioner, i.e. hybrid
! Schwarz, using RAS (with overlap 1 and ILU(0) on the blocks) 
! as post-smoother and 4 block-Jacobi sweeps (with UMFPACK LU
! on the blocks) as distributed coarse-level solver
  call mld_precinit(P,'ML',info)
!
! build the preconditioner
  call mld_precbld(A,desc_A,P,info)
!
! set the solver parameters and the initial guess
  ... ...
!
! solve Ax=b with preconditioned BiCGSTAB
  call psb_krylov('BICGSTAB',A,P,b,x,tol,desc_A,info)
  ... ...
!
! deallocate the preconditioner
  call mld_precfree(P,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>

<P>
Different versions of multi-level preconditioners can be obtained by changing
the default values of the preconditioner parameters. The code reported in
Figure&nbsp;<A HREF="#fig:ex_3lh">3</A> shows how to set a three-level hybrid Schwarz
preconditioner, which uses block Jacobi with ILU(0) on the
local blocks as post-smoother, has a coarsest matrix replicated on the processors,
and solves the coarsest-level system with the LU factorization from UMFPACK&nbsp;[<A
 HREF="node25.html#UMFPACK">9</A>].
The number of levels is specified by using <code>mld_precinit</code>; the other
preconditioner parameters are set by calling <code>mld_precset</code>. Note that
the type of multilevel framework (i.e. multiplicative among the levels
with post-smoothing only) is not specified since it is the default 
set by <code>mld_precinit</code>. 

<P>
Figure&nbsp;<A HREF="#fig:ex_3la">4</A> shows how to
set a three-level additive Schwarz preconditioner,
which uses RAS, with overlap 1 and ILU(0) on the blocks, 
as pre- and post-smoother, and applies five block-Jacobi sweeps, with
the UMFPACK LU factorization on the blocks, as distributed coarsest-level
solver. Again, <code>mld_precset</code> is used only to set
non-default values of the parameters (see Tables&nbsp;<A HREF="#tab:p_type">2</A>-<A HREF="#tab:p_coarse">5</A>).
In both cases, the construction and the application of the preconditioner
are carried out as for the default multi-level preconditioner.
The code fragments shown in in Figures&nbsp;<A HREF="#fig:ex_3lh">3</A>-<A HREF="#fig:ex_3la">4</A> are
included in the example program file <code>mld_dexample_ml.f90</code> too.

<P>
Finally, Figure&nbsp;<A HREF="#fig:ex_1l">5</A> shows the setup of a one-level
additive Schwarz preconditioner, i.e. RAS with overlap 2. The corresponding
example program is available in <code>mld_dexample_</code> <code>1lev.f90</code>.

<P>
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>.

<P>
 
<BR><B>Remark 3.</B> Any PSBLAS-based program using the basic preconditioners
implemented in PSBLAS 2.0, i.e. the diagonal and block-Jacobi ones,
can use the diagonal and block-Jacobi preconditioners
implemented in MLD2P4 without any change in the code.
The PSBLAS-based program must be only recompiled
and linked to the MLD2P4 library.

<BR>
<P>

<DIV ALIGN="CENTER"><A NAME="fig:ex_3lh"></A><A NAME="934"></A>
<TABLE>
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 3:</STRONG>
Setup of a hybrid three-level Schwarz preconditioner.</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
</DIV><TABLE  WIDTH="90%">
<TR><TD> 
<PRE>
... ...
! set a three-level hybrid Schwarz preconditioner, which uses 
! block Jacobi (with ILU(0) on the blocks) as post-smoother,
! a coarsest matrix replicated on the processors, and the 
! LU factorization from UMFPACK as coarse-level solver
  call mld_precinit(P,'ML',info,nlev=3)
  call_mld_precset(P,mld_smoother_type_,'BJAC',info)
  call mld_precset(P,mld_coarse_mat_,'REPL',info)
  call mld_precset(P,mld_coarse_solve_,'UMF',info)
... ...
</PRE>
</TD></TR>
</TABLE>
<DIV ALIGN="CENTER">
</DIV>
<P>
<DIV ALIGN="CENTER">
</DIV></TD></TR>
</TABLE>
</DIV>

<P>

<DIV ALIGN="CENTER"><A NAME="fig:ex_3la"></A><A NAME="936"></A>
<TABLE>
<CAPTION ALIGN="BOTTOM"><STRONG>Figure 4:</STRONG>
Setup of an additive three-level Schwarz preconditioner.</CAPTION>
<TR><TD>
<DIV ALIGN="CENTER">
</DIV><TABLE  WIDTH="90%">
<TR><TD> 
<PRE>
... ...
! set a three-level additive Schwarz preconditioner, which uses 
! RAS (with overlap 1 and ILU(0) on the blocks) as pre- and 
! post-smoother, and 5 block-Jacobi sweeps (with UMFPACK LU
! on the blocks) as distributed coarsest-level solver
  call mld_precinit(P,'ML',info,nlev=3)
  call mld_precset(P,mld_ml_type_,'ADD',info)
  call_mld_precset(P,mld_smoother_pos_,'TWOSIDE',info)
  call mld_precset(P,mld_coarse_sweeps_,5,info)
... ...
</PRE>
</TD></TR>
</TABLE>
<DIV ALIGN="CENTER">
</DIV>
<P>
<DIV ALIGN="CENTER">
</DIV></TD></TR>
</TABLE>
</DIV>

<P>

<DIV ALIGN="CENTER"><A NAME="fig:ex_1l"></A><A NAME="938"></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 mld_precinit(P,'AS',info)
  call mld_precset(P,mld_sub_ovr_,2,info)
... ...
</PRE>
</TD></TR>
</TABLE>
<DIV ALIGN="CENTER">

</DIV></TD></TR>
</TABLE>
</DIV>

<P>
<HR>
<!--Navigation Panel-->
<A NAME="tex2html242"
  HREF="node16.html">
<IMG WIDTH="37" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="next" SRC="next.png"></A> 
<A NAME="tex2html238"
  HREF="node14.html">
<IMG WIDTH="26" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="up" SRC="up.png"></A> 
<A NAME="tex2html234"
  HREF="node14.html">
<IMG WIDTH="63" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="previous" SRC="prev.png"></A> 
<A NAME="tex2html240"
  HREF="node2.html">
<IMG WIDTH="65" HEIGHT="24" ALIGN="BOTTOM" BORDER="0" ALT="contents" SRC="contents.png"></A>  
<BR>
<B> Next:</B> <A NAME="tex2html243"
  HREF="node16.html">User Interface</A>
<B> Up:</B> <A NAME="tex2html239"
  HREF="node14.html">Getting Started</A>
<B> Previous:</B> <A NAME="tex2html235"
  HREF="node14.html">Getting Started</A>
 &nbsp; <B>  <A NAME="tex2html241"
  HREF="node2.html">Contents</A></B> 
<!--End of Navigation Panel-->

</BODY>
</HTML>