[petsc-users] Question with preconditioned resid norm and true resid norm
Song Gao
song.gao2 at mail.mcgill.ca
Fri Apr 18 08:50:22 CDT 2014
Yeah. Thank you all. I think there are bugs in my function evaluation. I
just wondering if I could guess where the bug might be from the residual
history. I would check my function evaluation. Thank you very much.
I tried add -ksp_norm_type unpreconditioned, but got an error.
run with
-ksp_max_it 10 -ksp_gmres_restart 30 -ksp_monitor_true_residual
-ksp_converged_reason -ksp_view -ksp_norm_type unpreconditioned
-ksp_pc_side right -pc_type none
The error is
[0]PETSC ERROR: --------------------- Error Message
------------------------------------
[0]PETSC ERROR: No support for this operation for this object type!
[0]PETSC ERROR: KSP gmres does not support UNPRECONDITIONED with LEFT!
[0]PETSC ERROR:
------------------------------------------------------------------------
[0]PETSC ERROR: Petsc Release Version 3.3.0, Patch 6, Mon Feb 11 12:26:34
CST 2013
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.
[0]PETSC ERROR:
------------------------------------------------------------------------
On Thu, Apr 17, 2014 at 5:28 PM, Jed Brown <jed at jedbrown.org> wrote:
> Song Gao <song.gao2 at mail.mcgill.ca> writes:
>
> > Hello,
> >
> > I am using KSP framework to solve the problem A \Delta x = b in the
> > matrix free fashion. where A is the matrix free matrix. I have another
> > assembled matrix for preconditioning, but I'm NOT using it. My code is
> not
> > working so I'm debugging it.
> > I run the code with options -ksp_pc_side left -ksp_max_it 10
> > -ksp_gmres_restart 30 -ksp_monitor_true_residual -pc_type none -ksp_view
> >
> > I think if pc_type is none, the precondiitoned resid norm should equal to
> > true resid norm (Am I correct?). But this doesn't happen. So maybe it
> would
> > be helpful if I know how preconditioned resid norm and true resid norm
> are
> > computed.
> >
> > Website says true resid norm is just b - A \Delta x. But how is
> > preconditioned resid norm computed?
>
> Just apply the preconditioner P^{-1} to the residual above. In
> practice, it is usually computed indirectly via a recurrence in the
> Krylov method, but they should agree up to rounding error.
>
> > 0 KSP preconditioned resid norm 9.619278462343e-03 true resid norm
> > 9.619278462343e-03 ||r(i)||/||b|| 1.000000000000e+00
> > 1 KSP preconditioned resid norm 9.619210849854e-03 true resid norm
> > 2.552369536916e+06 ||r(i)||/||b|| 2.653389801437e+08
> > 2 KSP preconditioned resid norm 9.619210847390e-03 true resid norm
> > 2.552458142544e+06 ||r(i)||/||b|| 2.653481913988e+08
> > 3 KSP preconditioned resid norm 9.619210847385e-03 true resid norm
> > 2.552458343191e+06 ||r(i)||/||b|| 2.653482122576e+08
> > 4 KSP preconditioned resid norm 9.619210847385e-03 true resid norm
> > 2.552458344014e+06 ||r(i)||/||b|| 2.653482123432e+08
> > 5 KSP preconditioned resid norm 9.619210847385e-03 true resid norm
> > 2.552458344015e+06 ||r(i)||/||b|| 2.653482123433e+08
>
> Try the above with -ksp_norm_type unpreconditioned. The output below
> claims there is no preconditioner, so the most likely cause is that your
> operator is nonlinear. With MFFD, I would speculate that it is caused
> by your nonlinear function being discontinuous, or by reusing some
> memory without clearing it (thus computing nonsense after the first
> iteration).
>
> > Linear solve did not converge due to DIVERGED_ITS iterations 10
> > KSP Object: 1 MPI processes
> > type: gmres
> > GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
> > Orthogonalization with no iterative refinement
> > GMRES: happy breakdown tolerance 1e-30
> > maximum iterations=10, initial guess is zero
> > tolerances: relative=1e-06, absolute=1e-50, divergence=100000
> > left preconditioning
> > using PRECONDITIONED norm type for convergence test
> > PC Object: 1 MPI processes
> > type: none
> > linear system matrix followed by preconditioner matrix:
> > Matrix Object: 1 MPI processes
> > type: mffd
> > rows=22905, cols=22905
> > Matrix-free approximation:
> > err=1.49012e-08 (relative error in function evaluation)
> > Using wp compute h routine
> > Does not compute normU
> > Matrix Object: 1 MPI processes
> > type: seqbaij
> > rows=22905, cols=22905, bs=5
> > total: nonzeros=785525, allocated nonzeros=785525
> > total number of mallocs used during MatSetValues calls =0
> > block size is 5
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140418/4d0d6400/attachment.html>
More information about the petsc-users
mailing list