[petsc-users] how to switch ILU to its block version?

Barry Smith bsmith at mcs.anl.gov
Sun Mar 1 17:53:16 CST 2015


  Probably you want a better preconditioner than ILU (point or point block) for your problem. Just blindly using ILU on a matrix arising from complex physics is not very efficient and is fragile.

   You might start with the PCFIELDSPLIT where you have one "field" for the fluid and one for the "solid" part. Then you can use whatever appropriate preconditioners on each part; for example on the "solid" part you want to use PCGAMG (make sure you set the block size and the near null space for that one). The fluid part you might again use PCFIELDSPLIT on (yes it can be used recursively).

   Barry

> On Mar 1, 2015, at 5:03 PM, Fande Kong <fd.kong at siat.ac.cn> wrote:
> 
> Thanks, Matt.
> 
> Could bs be different from rows. I am solving a fluid-structure problem, which has different bs in fluid and solid parts.
> 
> Fande,
> 
> On Sun, Mar 1, 2015 at 3:57 PM, Matthew Knepley <knepley at gmail.com> wrote:
> On Sun, Mar 1, 2015 at 4:55 PM, Fande Kong <fd.kong at siat.ac.cn> wrote:
> Thanks, Barry.
> 
> More questions. Is bs usually the number of unknowns associated with a grid point. Do I need to change MatSetValues if use a MATSEQBAIJ matrix?
> 
> No, but it does enable you to use MatSetValuesBlock() if you want.
> 
>    Matt
>  
> Fande,
> 
> On Sun, Mar 1, 2015 at 2:36 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> 
>   Fande,
> 
>     The MATSEQBAIJ matrix automatically uses a point block ILU; with the block size of the matrix. So if your matrix naturally has a block size bs if you use the MATSEQBAIJ matrix type it will be point block ILU of block size bs.
> 
>     We don't have a variable point block ILU in PETSc for SeqAIJ. It is possible to write but would be rather involved.
> 
>    Barry
> 
> > On Mar 1, 2015, at 1:47 PM, Fande Kong <fd.kong at siat.ac.cn> wrote:
> >
> > Hi all,
> >
> > I have a code with ILU as subdomain solver, and it works well. Now I want to try a block version of ILU, there are some easy ways to implement this?
> >
> > Fande,
> 
> 
> 
> 
> 
> 
> -- 
> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
> -- Norbert Wiener
> 



More information about the petsc-users mailing list