[petsc-users] Issue of mg_coarse_ksp not converge

Matthew Knepley knepley at gmail.com
Mon Oct 2 10:29:56 CDT 2017


On Mon, Oct 2, 2017 at 11:21 AM, Mark Adams <mfadams at lbl.gov> wrote:

> Yea, it fails in the eigen estimator, but the Cheby eigen estimator works
> in the solve that works:
>
>         eigenvalue estimates used:  min = 0.100004, max = 1.10004
>         eigenvalues estimate via gmres min 0.0118548, max 1.00004
>
> Why would it just give "KSPSolve has not converged". It is not supposed to
> converge ...
>

This sounds like a mistake with


http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/KSP/KSPSetErrorIfNotConverged.html

somewhere.

   Matt


> On Mon, Oct 2, 2017 at 11:11 AM, Matthew Knepley <knepley at gmail.com>
> wrote:
>
>> On Mon, Oct 2, 2017 at 11:06 AM, Wenbo Zhao <zhaowenbo.npic at gmail.com>
>> wrote:
>>
>>> Matt,
>>>
>>
>> Thanks Wenbo.
>>
>>
>>> Test 1 nonsmooth
>>> zhaowenbo at ubuntu:~/test_slepc/SPARK/spark$ make NCORE=1 runkr_nonsmooth
>>> mpirun -n 1 ./step-41 \
>>>    -st_ksp_type gmres  -st_ksp_view -st_ksp_monitor  \
>>>    -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 0 \
>>>    -mata AMAT.dat -matb BMAT.dat \
>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor  -log_view > log_nonsmooth 2>&1
>>>
>>> Test 2 smooth
>>> zhaowenbo at ubuntu:~/test_slepc/SPARK/spark$ make NCORE=1 runkr_smooth
>>> mpirun -n 1 ./step-41 \
>>>    -st_ksp_type gmres -st_ksp_view -st_ksp_monitor  \
>>>    -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 1 \
>>>    -mata AMAT.dat -matb BMAT.dat \
>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor -log_view > log_smooth 2>&1
>>> makefile:43: recipe for target 'runkr_smooth' failed
>>> make: *** [runkr_smooth] Error 91
>>>
>>>
>> Mark, the solve is not failing, its the construction of the interpolator
>> I think. Check out this stack
>>
>> [0]PETSC ERROR: --------------------- Error Message
>> --------------------------------------------------------------
>> [0]PETSC ERROR:
>> [0]PETSC ERROR: KSPSolve has not converged
>> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
>> for trouble shooting.
>> [0]PETSC ERROR: Petsc Release Version 3.8.0, unknown
>> [0]PETSC ERROR: ./step-41 on a arch-linux2-c-debug named ubuntu by
>> zhaowenbo Mon Oct  2 08:00:58 2017
>> [0]PETSC ERROR: Configure options --with-mpi=1 --with-shared-libraries=1
>> --with-64-bit-indices=1 --with-debugging=1
>> [0]PETSC ERROR: #1 KSPSolve() line 855 in /home/zhaowenbo/research/petsc
>> /petsc_git/src/ksp/ksp/interface/itfunc.c
>> [0]PETSC ERROR: #2 PCGAMGOptProlongator_AGG() line 1186 in
>> /home/zhaowenbo/research/petsc/petsc_git/src/ksp/pc/impls/gamg/agg.c
>> [0]PETSC ERROR: #3 PCSetUp_GAMG() line 528 in
>> /home/zhaowenbo/research/petsc/petsc_git/src/ksp/pc/impls/gamg/gamg.c
>> [0]PETSC ERROR: #4 PCSetUp() line 924 in /home/zhaowenbo/research/petsc
>> /petsc_git/src/ksp/pc/interface/precon.c
>> [0]PETSC ERROR: #5 KSPSetUp() line 378 in /home/zhaowenbo/research/petsc
>> /petsc_git/src/ksp/ksp/interface/itfunc.c
>> [0]PETSC ERROR: #6 STSetUp_Shift() line 129 in
>> /home/zhaowenbo/research/slepc/slepc_git/src/sys/classes/st/
>> impls/shift/shift.c
>> [0]PETSC ERROR: #7 STSetUp() line 281 in /home/zhaowenbo/research/slepc
>> /slepc_git/src/sys/classes/st/interface/stsolve.c
>> [0]PETSC ERROR: #8 EPSSetUp() line 273 in /home/zhaowenbo/research/slepc
>> /slepc_git/src/eps/interface/epssetup.c
>> [0]PETSC ERROR: #9 solve_diffusion_3d() line 1029 in src/diffu.c
>> [0]PETSC ERROR: #10 main() line 25 in src/main.c
>> [0]PETSC ERROR: PETSc Option Table entries:
>> [0]PETSC ERROR: -eps_monitor
>> [0]PETSC ERROR: -eps_ncv 10
>> [0]PETSC ERROR: -eps_nev 1
>> [0]PETSC ERROR: -log_view
>> [0]PETSC ERROR: -mata AMAT.dat
>> [0]PETSC ERROR: -matb BMAT.dat
>> [0]PETSC ERROR: -st_ksp_monitor
>> [0]PETSC ERROR: -st_ksp_type gmres
>> [0]PETSC ERROR: -st_ksp_view
>> [0]PETSC ERROR: -st_pc_gamg_agg_nsmooths 1
>> [0]PETSC ERROR: -st_pc_gamg_type agg
>> [0]PETSC ERROR: -st_pc_type gamg
>> [0]PETSC ERROR: ----------------End of Error Message -------send entire
>> error message to petsc-maint at mcs.anl.gov----------
>> ------------------------------------------------------------
>> --------------
>> MPI_ABORT was invoked on rank 0 in communicator MPI_COMM_WORLD
>> with errorcode 91.
>>
>>   Thanks,
>>
>>      Matt
>>
>>
>>> Thanks,
>>>
>>> Wenbo
>>>
>>> On Mon, Oct 2, 2017 at 10:48 PM, Matthew Knepley <knepley at gmail.com>
>>> wrote:
>>>
>>>> On Mon, Oct 2, 2017 at 10:43 AM, Wenbo Zhao <zhaowenbo.npic at gmail.com>
>>>> wrote:
>>>>
>>>>> Mark,
>>>>>
>>>>> Thanks for your reply.
>>>>>
>>>>> On Mon, Oct 2, 2017 at 9:51 PM, Mark Adams <mfadams at lbl.gov> wrote:
>>>>>
>>>>>> Please send the output with -st_ksp_view and -st_ksp_monitor and we
>>>>>> can start to debug it.
>>>>>>
>>>>>> Test 1 with nonsmooth and preonly is OK
>>>>> zhaowenbo at ubuntu:~/test_slepc/SPARK/spark$ make NCORE=1
>>>>> runkr_nonsmooth
>>>>> mpirun -n 1 ./step-41 \
>>>>>    -st_ksp_type gmres  -st_ksp_view -st_ksp_monitor  \
>>>>>    -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 0 \
>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>    -st_mg_coarse_ksp_type preonly   -st_mg_coarse_ksp_monitor  \
>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor  -log_view > log_nonsmooth 2>&1
>>>>>
>>>>> Test 2 smooth and preonly is not OK
>>>>> zhaowenbo at ubuntu:~/test_slepc/SPARK/spark$ make NCORE=1 runkr_smooth
>>>>> mpirun -n 1 ./step-41 \
>>>>>    -st_ksp_type gmres -st_ksp_view -st_ksp_monitor  \
>>>>>    -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 1 \
>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>    -st_mg_coarse_ksp_type preonly   -st_mg_coarse_ksp_monitor  \
>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor -log_view > log_smooth 2>&1
>>>>> makefile:43: recipe for target 'runkr_smooth' failed
>>>>> make: *** [runkr_smooth] Error 91
>>>>>
>>>>> Test 3 nonsmooth and gmres is not OK
>>>>> zhaowenbo at ubuntu:~/test_slepc/SPARK/spark$ make NCORE=1 runkr_gmres
>>>>> mpirun -n 1 ./step-41 \
>>>>>    -st_ksp_type gmres  -st_ksp_view -st_ksp_monitor  \
>>>>>    -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 0 \
>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>    -st_mg_coarse_ksp_type gmres  -st_mg_coarse_ksp_monitor
>>>>> -st_mg_coarse_ksp_rtol 1.0e-6 \
>>>>>
>>>>
>>>> DO NOT DO THIS. Please send the output where you do NOTHING to the
>>>> coarse solver.
>>>>
>>>>   Thanks,
>>>>
>>>>      Matt
>>>>
>>>>
>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor  -log_view > log_gmres 2>&1
>>>>> makefile:59: recipe for target 'runkr_gmres' failed
>>>>> make: *** [runkr_gmres] Error 91
>>>>>
>>>>> log-files is attached.
>>>>>
>>>>>
>>>>> You mentioned that B is not symmetric. I assume it is elliptic
>>>>>> (diffusion). Where does the asymmetry come from?
>>>>>>
>>>>>>
>>>>> It is a two-group diffusion equations, where group denotes neutron
>>>>> enegry discretisation.
>>>>> Matrix B consists of neutron diffusion/leakage term, removal term and
>>>>> minus neutron scatter source term between different energies, when matrix A
>>>>> denotes neutron fission source.
>>>>>
>>>>> Diffusion term(Laplace operator) is elliptic and symmetric. Removal
>>>>> term is diagonal only. However scatter term is asymmetry since scatter term
>>>>> from high energy to low energy is far greater than the term from low to
>>>>> high.
>>>>>
>>>>>
>>>>> Wenbo
>>>>>
>>>>>
>>>>>> On Mon, Oct 2, 2017 at 9:39 AM, Wenbo Zhao <zhaowenbo.npic at gmail.com>
>>>>>> wrote:
>>>>>>
>>>>>>> Matt,
>>>>>>> Thanks for your reply.
>>>>>>> For the defalt option doesnt work firstly( -st_ksp_type gmres
>>>>>>> -st_pc_type gamg -st_pc_gamg_type agg -st_pc_gamg_agg_nsmooths 1), I tried
>>>>>>> to test those options.
>>>>>>>
>>>>>>> Wenbo
>>>>>>>
>>>>>>> On Mon, Oct 2, 2017 at 9:08 PM, Matthew Knepley <knepley at gmail.com>
>>>>>>> wrote:
>>>>>>>
>>>>>>>> On Mon, Oct 2, 2017 at 8:30 AM, Wenbo Zhao <
>>>>>>>> zhaowenbo.npic at gmail.com> wrote:
>>>>>>>>
>>>>>>>>> Matt
>>>>>>>>>
>>>>>>>>> Because I am not clear about what will happen using 'preonly' for
>>>>>>>>> large scale problem.
>>>>>>>>>
>>>>>>>>
>>>>>>>> The size of the problem has nothing to do with 'preonly'. All it
>>>>>>>> means is to apply a preconditioner without a Krylov solver.
>>>>>>>>
>>>>>>>>
>>>>>>>>> It seems to use a direct solver from below,
>>>>>>>>> http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/
>>>>>>>>> KSP/KSPPREONLY.html
>>>>>>>>>
>>>>>>>>
>>>>>>>> However, I still cannot understand why you would change the default?
>>>>>>>>
>>>>>>>>   Matt
>>>>>>>>
>>>>>>>>
>>>>>>>>>
>>>>>>>>> Thanks!
>>>>>>>>> Wenbo
>>>>>>>>>
>>>>>>>>> On Mon, Oct 2, 2017 at 5:09 PM, Matthew Knepley <knepley at gmail.com
>>>>>>>>> > wrote:
>>>>>>>>>
>>>>>>>>>> On Sun, Oct 1, 2017 at 9:53 PM, Wenbo Zhao <
>>>>>>>>>> zhaowenbo.npic at gmail.com> wrote:
>>>>>>>>>>
>>>>>>>>>>> Matt,
>>>>>>>>>>> Thanks for your reply.
>>>>>>>>>>> It DOES make no sense for this problem.
>>>>>>>>>>> But I am not clear about the 'preonly' option. Which solver is
>>>>>>>>>>> used in preonly? I wonder if 'preonly' is suitable for large scale problem
>>>>>>>>>>> such as 400,000,000 unknowns.
>>>>>>>>>>> So I tried 'gmres' option and found these error messages.
>>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>> I mean, why are you setting this at all. Just do not set the
>>>>>>>>>> coarse solver. The default should work fine.
>>>>>>>>>>
>>>>>>>>>>   Thanks,
>>>>>>>>>>
>>>>>>>>>>     Matt
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>> Could you give me some suggestions?
>>>>>>>>>>>
>>>>>>>>>>> Thanks.
>>>>>>>>>>>
>>>>>>>>>>> Wenbo
>>>>>>>>>>>
>>>>>>>>>>>
>>>>>>>>>>> On Mon, Oct 2, 2017 at 12:34 AM, Matthew Knepley <
>>>>>>>>>>> knepley at gmail.com> wrote:
>>>>>>>>>>>
>>>>>>>>>>>> On Sun, Oct 1, 2017 at 6:49 AM, Wenbo Zhao <
>>>>>>>>>>>> zhaowenbo.npic at gmail.com> wrote:
>>>>>>>>>>>>
>>>>>>>>>>>>> Hi,
>>>>>>>>>>>>>
>>>>>>>>>>>>> I met some questions when I use PETSC/SLEPC to solve two-group
>>>>>>>>>>>>> neutron diffusion equations with finite difference method. The grid is
>>>>>>>>>>>>> 3*3*3, when DOF on each points is 2. So the matrix size is 54*54.
>>>>>>>>>>>>> It is generalized eigenvalue problem Ax=\lamda Bx, where B is
>>>>>>>>>>>>> diagonally dominant matrix but not symmetry.
>>>>>>>>>>>>> EPS is set as below,
>>>>>>>>>>>>>  ierr = EPSSetProblemType(eps,EPS_GNHEP);CHKERRQ(ierr);¬
>>>>>>>>>>>>>  ierr = EPSSetWhichEigenpairs(eps,EPS_
>>>>>>>>>>>>> LARGEST_REAL);CHKERRQ(ierr);¬
>>>>>>>>>>>>>
>>>>>>>>>>>>> Krylovschur is used as eps sovler. GAMG is used as PC.
>>>>>>>>>>>>> I tried agg_nsmooths and mg_coarse_ksp_type. Only non-smooths
>>>>>>>>>>>>> and preonly is OK.
>>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>> Why are you setting the coarse solver. This makes no sense.
>>>>>>>>>>>>
>>>>>>>>>>>>    Thanks,
>>>>>>>>>>>>
>>>>>>>>>>>>     Matt
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>>>
>>>>>>>>>>>>> Test 1
>>>>>>>>>>>>> $ make NCORE=1 runkr_nonsmooth
>>>>>>>>>>>>> mpirun -n 1 ./step-41 \
>>>>>>>>>>>>>    -st_ksp_type gmres  \
>>>>>>>>>>>>>    -st_pc_type gamg -st_pc_gamg_type agg
>>>>>>>>>>>>> -st_pc_gamg_agg_nsmooths 0 \
>>>>>>>>>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>>>>>>>>>    -st_mg_coarse_ksp_type preonly   -st_mg_coarse_ksp_monitor
>>>>>>>>>>>>> \
>>>>>>>>>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor  -log_view >
>>>>>>>>>>>>> log_nonsmooth 2>&1
>>>>>>>>>>>>>
>>>>>>>>>>>>> Test 2
>>>>>>>>>>>>> $ make NCORE=1 runkr_smooth
>>>>>>>>>>>>> mpirun -n 1 ./step-41 \
>>>>>>>>>>>>>    -st_ksp_type gmres  \
>>>>>>>>>>>>>    -st_pc_type gamg -st_pc_gamg_type agg
>>>>>>>>>>>>> -st_pc_gamg_agg_nsmooths 1 \
>>>>>>>>>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>>>>>>>>>    -st_mg_coarse_ksp_type preonly   -st_mg_coarse_ksp_monitor
>>>>>>>>>>>>> \
>>>>>>>>>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor -log_view > log_smooth
>>>>>>>>>>>>> 2>&1
>>>>>>>>>>>>> makefile:43: recipe for target 'runkr_smooth' failed
>>>>>>>>>>>>> make: *** [runkr_smooth] Error 91
>>>>>>>>>>>>>
>>>>>>>>>>>>> Test 3
>>>>>>>>>>>>> $ make NCORE=1 runkr_gmres
>>>>>>>>>>>>> mpirun -n 1 ./step-41 \
>>>>>>>>>>>>>    -st_ksp_type gmres  \
>>>>>>>>>>>>>    -st_pc_type gamg -st_pc_gamg_type agg
>>>>>>>>>>>>> -st_pc_gamg_agg_nsmooths 0 \
>>>>>>>>>>>>>    -st_ksp_view  -mata AMAT.dat -matb BMAT.dat \
>>>>>>>>>>>>>    -st_mg_coarse_ksp_type gmres  -st_mg_coarse_ksp_monitor
>>>>>>>>>>>>> -st_mg_coarse_ksp_rtol 1.0e-6 \
>>>>>>>>>>>>>    -eps_nev 1 -eps_ncv 10  -eps_monitor  -log_view > log_gmres
>>>>>>>>>>>>> 2>&1
>>>>>>>>>>>>> makefile:59: recipe for target 'runkr_gmres' failed
>>>>>>>>>>>>> make: *** [runkr_gmres] Error 91
>>>>>>>>>>>>>
>>>>>>>>>>>>> Log files were attched.
>>>>>>>>>>>>> The matrix file were also attched as AMAT.dat and BMAT.dat.
>>>>>>>>>>>>>
>>>>>>>>>>>>> Is it correct? Or something wrong with my code or commad-line?
>>>>>>>>>>>>>
>>>>>>>>>>>>> Thanks!
>>>>>>>>>>>>>
>>>>>>>>>>>>> Wenbo
>>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>>
>>>>>>>>>>>> --
>>>>>>>>>>>> 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
>>>>>>>>>>>>
>>>>>>>>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>>>>>>>>> <http://www.caam.rice.edu/%7Emk51/>
>>>>>>>>>>>>
>>>>>>>>>>>
>>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>>
>>>>>>>>>> --
>>>>>>>>>> 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
>>>>>>>>>>
>>>>>>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>>>>>>> <http://www.caam.rice.edu/%7Emk51/>
>>>>>>>>>>
>>>>>>>>>
>>>>>>>>>
>>>>>>>>
>>>>>>>>
>>>>>>>> --
>>>>>>>> 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
>>>>>>>>
>>>>>>>> https://www.cse.buffalo.edu/~knepley/
>>>>>>>> <http://www.caam.rice.edu/%7Emk51/>
>>>>>>>>
>>>>>>>
>>>>>>>
>>>>>>
>>>>>
>>>>
>>>>
>>>> --
>>>> 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
>>>>
>>>> https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
>>>>
>>>
>>>
>>
>>
>> --
>> 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
>>
>> https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
>>
>
>


-- 
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

https://www.cse.buffalo.edu/~knepley/ <http://www.caam.rice.edu/~mk51/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171002/9d7337de/attachment-0001.html>


More information about the petsc-users mailing list