[petsc-users] Issue of mg_coarse_ksp not converge
Mark Adams
mfadams at lbl.gov
Mon Oct 2 10:45:45 CDT 2017
This is normal:
Linear st_gamg_est_ solve did not converge due to DIVERGED_ITS iterations 10
It looks like ksp->errorifnotconverged got set somehow. If the default
changed in KSP then (SAGG) GAMG would not ever work.
I assume you don't have a .petscrc file with more (crazy) options in it ...
On Mon, Oct 2, 2017 at 11:39 AM, Wenbo Zhao <zhaowenbo.npic at gmail.com>
wrote:
>
>
> On Mon, Oct 2, 2017 at 11:30 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
>
>> On Mon, Oct 2, 2017 at 11:15 AM, Mark Adams <mfadams at lbl.gov> wrote:
>>
>>> non-smoothed aggregation is converging very fast. smoothed fails in the
>>> eigen estimator.
>>>
>>> Run this again with -st_gamg_est_ksp_view and -st_gamg_est_ksp_monitor,
>>> and see if you get more output (I'm not 100% sure about these args).
>>>
>>
>> I also want -st_gamg_est_ksp_converged_reason
>>
>> Thanks,
>>
>> Matt
>>
> $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 \
> -st_gamg_est_ksp_view -st_gamg_est_ksp_monitor \
> -st_gamg_est_ksp_converged_reason \
> -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
>
> Thanks
> Wenbo
>
>
>>
>>
>>>
>>> On Mon, Oct 2, 2017 at 11:06 AM, Wenbo Zhao <zhaowenbo.npic at gmail.com>
>>> wrote:
>>>
>>>> Matt,
>>>>
>>>> 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
>>>>
>>>>
>>>> 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/%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/>
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20171002/0261710a/attachment-0001.html>
More information about the petsc-users
mailing list