flexible block matrix

Satish Balay balay at mcs.anl.gov
Mon Apr 21 10:33:39 CDT 2008


On Mon, 21 Apr 2008, Jed Brown wrote:

> I am solving a Stokes problem with nonlinear slip boundary conditions.  I don't
> think I can take advantage of block structure since the normal component of
> velocity has a Dirichlet constraint and this must be built into the velocity
> space in order to preserve conditioning.  An alternative formulation involves a
> Lagrange multiplier for the constraint, but even with clever preconditioning,
> this system is still more expensive to solve according to [1].
> 
> In solving the (velocity-pressure) saddle point problem, many approximate solves
> with the velocity system is needed in the preconditioner, hence I need a strong
> preconditioner for the velocity system.  Currently, I am using algebraic
> multigrid on a low-order discretization which works fairly well.  Since Hypre
> and ML only take AIJ matrices, perhaps I shouldn't worry about blocking after
> all.  Is there a way to use MATBAIJ when some nodes have fewer degrees of
> freedom?  Should I bother?

I'll say - don't bother. BAIJ can't support varing block size. The
code that supports it is INODE code - which is already part of AIJ
type - and is the default for AIJ.

You can run your code with -mat_no_inode to see the performance
difference between basic AIJ and INODE-AIJ. [The primary thing to look
for in -log_summary is MatMult()]

Inode code looks for consequitive *rows* with same column indices, and
marks them as a single inode. For each inode - [i.e say 5 rows] the
column indices are loaded only once, and used for all 5 rows - thus
improving the performance. A matrix can have an inode structure of
[2,2,3,3,1,3] etc.. i.e 14x14 matrix.

Satish

> Note that my method (currently just a single element) uses a high order
> discretization on some elements and low order on others.  The global matrix for
> the low order elements is assembled, but it is applied locally for the high order
> elements taking advantage of the tensor product basis.  For the preconditioner,
> a low order discretization on the nodes of the high order elements is globally
> assembled and added to the global matrix from the low-order elements.
> Experiments with a single element (spectral rather than spectral/hp element)
> show this to be effective, converging in a constant number of iterations
> independent of polynomial order when using a V-cycle of AMG as a preconditioner.
> 
> Thanks.
> 
> Jed
> 
> 
> [1] Bänsch, Höhn 2000, `Numerical treatment of the Navier-Stokes equations with
> slip boundary conditions', SIAM J. Sci. Comput.
> 
> 


More information about the petsc-users mailing list