[petsc-dev] asm / gasm

Fande Kong fdkong.jd at gmail.com
Sun Jul 3 11:41:20 CDT 2016


On Sun, Jul 3, 2016 at 10:33 AM, Mark Adams <mfadams at lbl.gov> wrote:

> I am using KSP ex56.
>
> I also found a bug and the code is running now.
>
> I do not understand the difference between the last two args here:
>
> PetscErrorCode  PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])
>
>
>
"is" is the overlapping subdomain, and "is_local" is the subdomain without
overlap.



> In looking at the web doc I would think I want to set the first, but I get
> an error when I do that.  I am using the second and am not getting any
> reduction in the number of iterations.  I'm still looking into this, but
> what should I use here?
>
>
> On Sun, Jul 3, 2016 at 6:17 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>
>>
>>    Mark,
>>
>>    When I try to run your example in your branch with one process I get
>> the error below; why is that if you do not get the error?
>>
>> ~/Src/petsc/src/snes/examples/tutorials (mark/fix-gamg-asm-aggs=)
>> arch-basic
>> $ ./ex56 -ne 3 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type agg
>> -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10
>> -pc_gamg_reuse_interpolation true -two_solves true -ksp_converged_reason
>> -use_mat_nearnullspace -mg_levels_esteig_ksp_type cg -pc_gamg_square_graph
>> 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig
>> 0,0.05,0,1.05 -gamg_est_ksp_type cg -gamg_est_ksp_max_it 20
>> -pc_gamg_asm_use_agg -mg_levels_sub_pc_type lu -ksp_monitor
>> -mat_block_size 3
>> [0]PETSC ERROR: --------------------- Error Message
>> --------------------------------------------------------------
>> [0]PETSC ERROR: Petsc has generated inconsistent data
>> [0]PETSC ERROR: !(matA_1 && !matA_1->compressedrow.use)
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
>> for trouble shooting.
>> [0]PETSC ERROR: Petsc Development GIT revision: v3.7.2-757-g1d88083  GIT
>> Date: 2016-06-28 21:49:49 -0500
>> [0]PETSC ERROR: ./ex56 on a arch-basic named Barrys-MacBook-Pro.local by
>> barrysmith Sun Jul  3 11:15:36 2016
>> [0]PETSC ERROR: Configure options --download-metis --download-parmetis
>> --with-mpi-dir=/Users/barrysmith/libraries PETSC_ARCH=arch-basic
>> [0]PETSC ERROR: #1 smoothAggs() line 333 in
>> /Users/barrysmith/Src/petsc/src/ksp/pc/impls/gamg/agg.c
>> [0]PETSC ERROR: #2 PCGAMGCoarsen_AGG() line 1001 in
>> /Users/barrysmith/Src/petsc/src/ksp/pc/impls/gamg/agg.c
>> [0]PETSC ERROR: #3 PCSetUp_GAMG() line 524 in
>> /Users/barrysmith/Src/petsc/src/ksp/pc/impls/gamg/gamg.c
>> [0]PETSC ERROR: #4 PCSetUp() line 968 in
>> /Users/barrysmith/Src/petsc/src/ksp/pc/interface/precon.c
>> [0]PETSC ERROR: #5 KSPSetUp() line 393 in
>> /Users/barrysmith/Src/petsc/src/ksp/ksp/interface/itfunc.c
>> [0]PETSC ERROR: #6 main() line 384 in
>> /Users/barrysmith/Src/petsc/src/snes/examples/tutorials/ex56.c
>> [0]PETSC ERROR: PETSc Option Table entries:
>> [0]PETSC ERROR: -alpha 1.e-3
>> [0]PETSC ERROR: -gamg_est_ksp_max_it 20
>> [0]PETSC ERROR: -gamg_est_ksp_type cg
>> [0]PETSC ERROR: -ksp_converged_reason
>> [0]PETSC ERROR: -ksp_monitor
>> [0]PETSC ERROR: -ksp_type cg
>> [0]PETSC ERROR: -malloc_test
>> [0]PETSC ERROR: -mat_block_size 3
>> [0]PETSC ERROR: -mg_levels_esteig_ksp_type cg
>> [0]PETSC ERROR: -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05
>> [0]PETSC ERROR: -mg_levels_ksp_type chebyshev
>> [0]PETSC ERROR: -mg_levels_sub_pc_type lu
>> [0]PETSC ERROR: -ne 3
>> [0]PETSC ERROR: -pc_gamg_agg_nsmooths 1
>> [0]PETSC ERROR: -pc_gamg_asm_use_agg
>> [0]PETSC ERROR: -pc_gamg_coarse_eq_limit 10
>> [0]PETSC ERROR: -pc_gamg_reuse_interpolation true
>> [0]PETSC ERROR: -pc_gamg_square_graph 1
>> [0]PETSC ERROR: -pc_gamg_type agg
>> [0]PETSC ERROR: -pc_type gamg
>> [0]PETSC ERROR: -two_solves true
>> [0]PETSC ERROR: -use_mat_nearnullspace
>> [0]PETSC ERROR: ----------------End of Error Message -------send entire
>> error message to petsc-maint at mcs.anl.gov----------
>> application called MPI_Abort(MPI_COMM_WORLD, 77) - process 0
>> [unset]: write_line error; fd=-1 buf=:cmd=abort exitcode=77
>> :
>> system msg for write_line failure : Bad file descriptor
>> ~/Src/petsc/src/snes/examples/tutorials (mark/fix-gamg-asm-aggs=)
>> arch-basic
>>
>> > On Jul 3, 2016, at 9:04 AM, Mark Adams <mfadams at lbl.gov> wrote:
>> >
>> > Barry, I have never seen an empty matrix for this problem.  I am in my
>> branch (mark/fix-gamg-asm-aggs), which basically just has the fix to add a
>> dummy block for BCs to make ASM work properly, added to your branch.
>> >
>> > This (appended) is what I get with 1 proc (it is fine, it is solving) .
>> With 8 procs I get this stack trace.  I will keep digging.
>> >
>> >
>> > #8 main (argc=36, args=0x7fffffff4628) at
>> /global/homes/m/madams/petsc/src/ksp/ksp/examples/tutorials/ex56.c:280 (at
>> 0x0000000000409f2e)
>> > #7 KSPSetUp (ksp=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/ksp/interface/itfunc.c:393 (at
>> 0x000000000075d99f)
>> > #6 PCSetUp (pc=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/pc/interface/precon.c:968 (at
>> 0x00000000007531da)
>> > #5 PCSetUp_GAMG (pc=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/pc/impls/gamg/gamg.c:672 (at
>> 0x0000000000c28a81)
>> > #4 PCSetUp_MG (pc=0x4760330) at
>> /global/u2/m/madams/petsc/src/ksp/pc/impls/mg/mg.c:781 (at
>> 0x0000000000c87f7e)
>> > #3 KSPSetUp (ksp=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/ksp/interface/itfunc.c:393 (at
>> 0x000000000075d99f)
>> > #2 PCSetUp (pc=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/pc/interface/precon.c:968 (at
>> 0x00000000007531da)
>> > #1 PCSetUp_ASM (pc=0x46e2f20) at
>> /global/u2/m/madams/petsc/src/ksp/pc/impls/asm/asm.c:383 (at
>> 0x0000000000c5c499)
>> > #0 MatGetSubMatrices (mat=0x46e2f20, n=63418000, irow=0x21, icol=0x4a8,
>> scall=MAT_INITIAL_MATRIX, submat=0x2) at
>> /global/u2/m/madams/petsc/src/mat/interface/matrix.c:6780 (at
>> 0x00000000004e8bcd)
>> >
>> >
>> >
>> > 06:32 134 nid00015  ~/petsc/src/ksp/ksp/examples/tutorials$ srun -n 1
>> ./ex56 -ne 3 -alpha 1.e-3 -ksp_type cg -pc_type gamg -pc_gamg_type agg
>> -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10
>> -pc_gamg_reuse_interpolation true -two_solves true -ksp_converged_reason
>> -use_mat_nearnullspace -mg_levels_esteig_ksp_type cg -pc_gamg_square_graph
>> 1 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig
>> 0,0.05,0,1.05 -gamg_est_ksp_type cg -gamg_est_ksp_max_it 20
>> -pc_gamg_asm_use_agg -mg_levels_sub_pc_type lu -ksp_monitor
>> >   0 KSP Residual norm 2.553528059796e+02
>> >   1 KSP Residual norm 1.285760737388e+01
>> >   2 KSP Residual norm 5.896306650166e-01
>> >   3 KSP Residual norm 3.206061310861e-02
>> >   4 KSP Residual norm 3.287356003543e-03
>> >   5 KSP Residual norm 1.695272457611e-04
>> > Linear solve converged due to CONVERGED_RTOL iterations 5
>> >   0 KSP Residual norm 2.553528059796e-03
>> >   1 KSP Residual norm 1.285760737388e-04
>> >   2 KSP Residual norm 5.896306650166e-06
>> >   3 KSP Residual norm 3.206061310861e-07
>> >   4 KSP Residual norm 3.287356003543e-08
>> >   5 KSP Residual norm 1.695272457611e-09
>> > Linear solve converged due to CONVERGED_RTOL iterations 5
>> >   0 KSP Residual norm 2.553528059796e-08
>> >   1 KSP Residual norm 1.285760737388e-09
>> >   2 KSP Residual norm 5.896306650166e-11
>> >   3 KSP Residual norm 3.206061310861e-12
>> >   4 KSP Residual norm 3.287356003543e-13
>> >   5 KSP Residual norm 1.695272457611e-14
>> > Linear solve converged due to CONVERGED_RTOL iterations 5
>> > [0]main |b-Ax|/|b|=9.327411e-06, |b|=7.453560e+00, emax=9.907752e-01
>> >
>> > On Wed, Jun 29, 2016 at 1:05 AM, Barry Smith <bsmith at mcs.anl.gov>
>> wrote:
>> >
>> >    Mark
>> >
>> >    When I run with exactly the options in this file on either 1 or 8
>> process and MatView() the matrix you pass to KSP the matrix has ALL zeros.
>> So of course the GAMG solver is not happy. Is the matrix suppose to have
>> all zeros? When I run the command line options in the make file for runex56
>> the matrix has nonzero values in it and GAMG runs fine.
>> >
>> >    Barry
>> >
>> > > On Jun 27, 2016, at 6:56 AM, Mark Adams <mfadams at lbl.gov> wrote:
>> > >
>> > > Barry,
>> > >
>> > > I am having problems with your branch on ksp/ex56. First I was
>> getting an error with the increase of the overlap that you added (I don't
>> want that to be the default).  So I set this back to 0 in gamg.c.
>> > >
>> > > Next, I get the same error with the vec out of range because ex56 has
>> BCs in the matrix.  I created a branch mark/fix-gamg-asm-aggs from your
>> branch and added the dummy ASM block. But I am getting all kinds of
>> errors.  I get errors in get MatSubMatrix and a funny error when it tries
>> to free a prefix string.
>> > >
>> > > I've attached the valgrind output of one process, of 8.  There is a
>> lot of noise with string operations but it looks like the first non-string
>> thing is in MatGetSubMatrices.  ex56 is elasticity. Maybe there is a block
>> size problem?
>> > >
>> > >
>> > > <out_val_1>
>> >
>> >
>>
>>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-dev/attachments/20160703/9dbd19ca/attachment.html>


More information about the petsc-dev mailing list