The code reported in Figure 2 shows how to set and apply the default
multi-level preconditioner available in the real double precision version
of MLD2P4 (see Table 1). This preconditioner is chosen
by simply specifying 'ML'
as second argument of mld_precinit
(a call to mld_precset
is not needed) and is applied with the BiCGSTAB
solver provided by PSBLAS. As previously observed, the modules psb_base_mod
,
mld_prec_mod
and psb_krylov_mod
must be used by the example program.
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 mld_dexample_ml.f90
,
in the directory examples/fileread
of the MLD2P4 tree (see
Section 3.5).
For details on the use of the PSBLAS routines, see the PSBLAS User's
Guide [15].
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 6 for details). If these versions are installed,
the corresponding Fortran 95 codes are available in examples/fileread/
.
|
Different versions of multi-level preconditioners can be obtained by changing
the default values of the preconditioner parameters. The code reported in
Figure 3 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 [9].
The number of levels is specified by using mld_precinit
; the other
preconditioner parameters are set by calling mld_precset
. 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 mld_precinit
.
Figure 4 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, mld_precset
is used only to set
non-default values of the parameters (see Tables 2-5).
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 3-4 are
included in the example program file mld_dexample_ml.f90
too.
Finally, Figure 5 shows the setup of a one-level
additive Schwarz preconditioner, i.e. RAS with overlap 2. The corresponding
example program is available in mld_dexample_
1lev.f90
.
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 examples/pdegen
.
Remark 3. 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.
|
|
|