[petsc-users] Why use MATMPIBAIJ?

Hoang Giang Bui hgbk2008 at gmail.com
Tue Jan 26 02:58:20 CST 2016


Hi

I assert this line to the hypre.c to see what block size it set to

/* special case for BoomerAMG */
  if (jac->setup == HYPRE_BoomerAMGSetup) {
    ierr = MatGetBlockSize(pc->pmat,&bs);CHKERRQ(ierr);

    // check block size passed to HYPRE
    PetscPrintf(PetscObjectComm((PetscObject)pc),"the block size passed to
HYPRE is %d\n",bs);

    if (bs > 1)
PetscStackCallStandard(HYPRE_BoomerAMGSetNumFunctions,(jac->hsolver,bs));
  }

It shows that the passing block size is 1. So my hypothesis is correct.

In the manual of MatSetBlockSize (
http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/Mat/MatSetBlockSize.html),
it has to be called before MatSetUp. Hence I guess the matrix passed to
HYPRE is created before I set the block size. Given that, I set the block
size after the call to PCFieldSplitSetIS

        ierr = PCFieldSplitSetIS(pc, "u", IS_u); CHKERRQ(ierr);
        ierr = PCFieldSplitSetIS(pc, "p", IS_p); CHKERRQ(ierr);

        /*
            Set block size for sub-matrix,
        */
        ierr = PCFieldSplitGetSubKSP(pc, &nsplits, &sub_ksp); CHKERRQ(ierr);
        ksp_U = sub_ksp[0];
        ierr = KSPGetOperators(ksp_U, &A_U, &P_U); CHKERRQ(ierr);
        ierr = MatSetBlockSize(A_U, 3); CHKERRQ(ierr);
        ierr = MatSetBlockSize(P_U, 3); CHKERRQ(ierr);

I guess the sub-matrices is created at PCFieldSplitSetIS. If that's correct
then it's not possible to set the block size this way.


Giang

On Mon, Jan 25, 2016 at 7:43 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

>
> > On Jan 25, 2016, at 11:13 AM, Hoang Giang Bui <hgbk2008 at gmail.com>
> wrote:
> >
> > OK, let's come back to my problem. I got your point about the
> interaction between components in one block. In my case, the interaction is
> strong.
> >
> > As you said, I try this:
> >
> >         ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
> >         ierr = PCFieldSplitGetSubKSP(pc, &nsplits, &sub_ksp);
> CHKERRQ(ierr);
> >         ksp_U = sub_ksp[0];
> >         ierr = KSPGetOperators(ksp_U, &A_U, &P_U); CHKERRQ(ierr);
> >         ierr = MatSetBlockSize(A_U, 3); CHKERRQ(ierr);
> >         ierr = MatSetBlockSize(P_U, 3); CHKERRQ(ierr);
> >         ierr = PetscFree(sub_ksp); CHKERRQ(ierr);
> >
> > But it seems doesn't work. The output from -ksp_view shows that matrix
> passed to Hypre still has bs=1
>
>    Hmm, this is strange. MatSetBlockSize() should have either set the
> block size to 3 or generated an error.  Can you run in the debugger on one
> process and put a break point in MatSetBlockSize() and see what it is
> setting the block size to. Then in PCSetUp_hypre() you can see what it is
> passing to hypre as the block size and maybe figure out how it becomes 1.
>
>   Barry
>
>
> >
> > KSP Object:    (fieldsplit_u_)     8 MPI processes
> >       type: preonly
> >       maximum iterations=10000, initial guess is zero
> >       tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
> >       left preconditioning
> >       using NONE norm type for convergence test
> >     PC Object:    (fieldsplit_u_)     8 MPI processes
> >       type: hypre
> >         HYPRE BoomerAMG preconditioning
> >         HYPRE BoomerAMG: Cycle type V
> >         HYPRE BoomerAMG: Maximum number of levels 25
> >         HYPRE BoomerAMG: Maximum number of iterations PER hypre call 1
> >         HYPRE BoomerAMG: Convergence tolerance PER hypre call 0
> >         HYPRE BoomerAMG: Threshold for strong coupling 0.25
> >         HYPRE BoomerAMG: Interpolation truncation factor 0
> >         HYPRE BoomerAMG: Interpolation: max elements per row 0
> >         HYPRE BoomerAMG: Number of levels of aggressive coarsening 0
> >         HYPRE BoomerAMG: Number of paths for aggressive coarsening 1
> >         HYPRE BoomerAMG: Maximum row sums 0.9
> >         HYPRE BoomerAMG: Sweeps down         1
> >         HYPRE BoomerAMG: Sweeps up           1
> >         HYPRE BoomerAMG: Sweeps on coarse    1
> >         HYPRE BoomerAMG: Relax down          symmetric-SOR/Jacobi
> >         HYPRE BoomerAMG: Relax up            symmetric-SOR/Jacobi
> >         HYPRE BoomerAMG: Relax on coarse     Gaussian-elimination
> >         HYPRE BoomerAMG: Relax weight  (all)      1
> >         HYPRE BoomerAMG: Outer relax weight (all) 1
> >         HYPRE BoomerAMG: Using CF-relaxation
> >         HYPRE BoomerAMG: Measure type        local
> >         HYPRE BoomerAMG: Coarsen type        PMIS
> >         HYPRE BoomerAMG: Interpolation type  classical
> >       linear system matrix = precond matrix:
> >       Mat Object:      (fieldsplit_u_)       8 MPI processes
> >         type: mpiaij
> >         rows=792333, cols=792333
> >         total: nonzeros=1.39004e+08, allocated nonzeros=1.39004e+08
> >         total number of mallocs used during MatSetValues calls =0
> >           using I-node (on process 0) routines: found 30057 nodes, limit
> used is 5
> >
> > In other test, I can see the block size bs=3 in the section of Mat Object
> >
> > Regardless the setup cost of Hypre AMG, I saw it gives quite a radical
> performance, providing that the material parameters does not vary strongly,
> and the geometry is regular enough.
> >
> >
> > Giang
> >
> > On Fri, Jan 22, 2016 at 2:57 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
> > On Fri, Jan 22, 2016 at 7:27 AM, Hoang Giang Bui <hgbk2008 at gmail.com>
> wrote:
> > DO you mean the option pc_fieldsplit_block_size? In this thread:
> >
> > http://petsc-users.mcs.anl.narkive.com/qSHIOFhh/fieldsplit-error
> >
> > No. "Block Size" is confusing on PETSc since it is used to do several
> things. Here block size
> > is being used to split the matrix. You do not need this since you are
> prescribing your splits. The
> > matrix block size is used two ways:
> >
> >   1) To indicate that matrix values come in logically dense blocks
> >
> >   2) To change the storage to match this logical arrangement
> >
> > After everything works, we can just indicate to the submatrix which is
> extracted that it has a
> > certain block size. However, for the Laplacian I expect it not to matter.
> >
> > It assumes you have a constant number of fields at each grid point, am I
> right? However, my field split is not constant, like
> > [u1_x   u1_y    u1_z    p_1    u2_x    u2_y    u2_z    u3_x    u3_y
> u3_z    p_3    u4_x    u4_y    u4_z]
> >
> > Subsequently the fieldsplit is
> > [u1_x   u1_y    u1_z    u2_x    u2_y    u2_z    u3_x    u3_y    u3_z
>  u4_x    u4_y    u4_z]
> > [p_1    p_3]
> >
> > Then what is the option to set block size 3 for split 0?
> >
> > Sorry, I search several forum threads but cannot figure out the options
> as you said.
> >
> >
> >
> > You can still do that. It can be done with options once the
> decomposition is working. Its true that these solvers
> > work better with the block size set. However, if its the P2 Laplacian it
> does not really matter since its uncoupled.
> >
> > Yes, I agree it's uncoupled with the other field, but the crucial factor
> defining the quality of the block preconditioner is the approximate
> inversion of individual block. I would merely try block Jacobi first,
> because it's quite simple. Nevertheless, fieldsplit implements other nice
> things, like Schur complement, etc.
> >
> > I think concepts are getting confused here. I was talking about the
> interaction of components in one block (the P2 block). You
> > are talking about interaction between blocks.
> >
> >   Thanks,
> >
> >      Matt
> >
> > Giang
> >
> >
> >
> > On Fri, Jan 22, 2016 at 11:15 AM, Matthew Knepley <knepley at gmail.com>
> wrote:
> > On Fri, Jan 22, 2016 at 3:40 AM, Hoang Giang Bui <hgbk2008 at gmail.com>
> wrote:
> > Hi Matt
> > I would rather like to set the block size for block P2 too. Why?
> >
> > Because in one of my test (for problem involves only [u_x u_y u_z]), the
> gmres + Hypre AMG converges in 50 steps with block size 3, whereby it
> increases to 140 if block size is 1 (see attached files).
> >
> > You can still do that. It can be done with options once the
> decomposition is working. Its true that these solvers
> > work better with the block size set. However, if its the P2 Laplacian it
> does not really matter since its uncoupled.
> >
> > This gives me the impression that AMG will give better inversion for
> "P2" block if I can set its block size to 3. Of course it's still an
> hypothesis but worth to try.
> >
> > Another question: In one of the Petsc presentation, you said the Hypre
> AMG does not scale well, because set up cost amortize the iterations. How
> is it quantified? and what is the memory overhead?
> >
> > I said the Hypre setup cost is not scalable, but it can be amortized
> over the iterations. You can quantify this
> > just by looking at the PCSetUp time as your increase the number of
> processes. I don't think they have a good
> > model for the memory usage, and if they do, I do not know what it is.
> However, generally Hypre takes more
> > memory than the agglomeration MG like ML or GAMG.
> >
> >   Thanks,
> >
> >     Matt
> >
> >
> > Giang
> >
> > On Mon, Jan 18, 2016 at 5:25 PM, Jed Brown <jed at jedbrown.org> wrote:
> > Hoang Giang Bui <hgbk2008 at gmail.com> writes:
> >
> > > Why P2/P2 is not for co-located discretization?
> >
> > Matt typed "P2/P2" when me meant "P2/P1".
> >
> >
> >
> >
> > --
> > 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
> >
> >
> >
> >
> > --
> > 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
> >
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160126/27999d32/attachment-0001.html>


More information about the petsc-users mailing list