On Aug 28, 2012, at 8:34 AM, Jed Brown <jedbrown at mcs.anl.gov> wrote:
> On Tue, Aug 28, 2012 at 8:19 AM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> > How to interpret the results? For the first case, is the null space
> > information forwarded to the preconditioner? And even for using just gmres
> > without any preconditioner, the dot product is much closer to zero but
> > still far away.
>
>
> Are you sure the null space you provide is truly the null space of the
> matrix? Call MatNullSpaceTest().
>
> GMRES uses left preconditioning by default which forces the null space
> into the Krylov space at each iteration. FGMRES uses right preconditioning
> and doesn't do anything to the null space because it is naturally in the
> Krylov space (since the result with fgmres contains the constants this is why
> I suspect the matrix actually doesn't have the null space you think it has).
> You can try running -ksp_type gmres -ksp_pc_side right my guess is that
> you will not get near zero for the dot.
>
> Barry, it looks like null space removal is not implemented correctly for
> right preconditioning:
Is this possible to fix? Don't we need the null space of A B , which would
require solving with B to get (from the null space of A) and in general we
cannot solve with B, only apply it?
Barry
>
> diff --git a/src/ksp/ksp/examples/tutorials/ex29.c
> b/src/ksp/ksp/examples/tutorials/ex29.c
> --- a/src/ksp/ksp/examples/tutorials/ex29.c
> +++ b/src/ksp/ksp/examples/tutorials/ex29.c
> @@ -78,6 +78,12 @@
> ierr = KSPSolve(ksp,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
> ierr = KSPGetSolution(ksp,&x);CHKERRQ(ierr);
> ierr = KSPGetRhs(ksp,&b);CHKERRQ(ierr);
> + {
> + PetscScalar dot;
> + ierr = VecSet(b,1.0);CHKERRQ(ierr);
> + ierr = VecDot(b,x,&dot);CHKERRQ(ierr);
> + ierr = PetscPrintf(PETSC_COMM_WORLD,"dot %g\n",dot);CHKERRQ(ierr);
> + }
>
> ierr = DMDestroy(&da);CHKERRQ(ierr);
> ierr = KSPDestroy(&ksp);CHKERRQ(ierr);
>
> $ ./ex29 -bc_type neumann -ksp_monitor_true_residual -ksp_test_null_space
> -mat_null_space_test_view -pc_type jacobi -ksp_pc_side left
> Constants are likely null vector|| A * 1/N || = 0
> 0 KSP preconditioned resid norm 1.135199672180e-01 true resid norm
> 2.319967097991e-01 ||r(i)||/||b|| 1.000000000000e+00
> 1 KSP preconditioned resid norm 4.410876553096e-02 true resid norm
> 1.254974589275e-01 ||r(i)||/||b|| 5.409449945916e-01
> 2 KSP preconditioned resid norm 2.255575515098e-02 true resid norm
> 5.942041105213e-02 ||r(i)||/||b|| 2.561260937863e-01
> 3 KSP preconditioned resid norm 9.933310900089e-03 true resid norm
> 2.093328265899e-02 ||r(i)||/||b|| 9.023094627987e-02
> 4 KSP preconditioned resid norm 1.073259860881e-17 true resid norm
> 1.660357840786e-16 ||r(i)||/||b|| 7.156816328229e-16
> dot 1.73472e-16
> $ ./ex29 -bc_type neumann -ksp_monitor_true_residual -ksp_test_null_space
> -mat_null_space_test_view -pc_type jacobi -ksp_pc_side right
> Constants are likely null vector|| A * 1/N || = 0
> 0 KSP unpreconditioned resid norm 2.319967097991e-01 true resid norm
> 2.319967097991e-01 ||r(i)||/||b|| 1.000000000000e+00
> 1 KSP unpreconditioned resid norm 1.218421389101e-01 true resid norm
> 1.218421389101e-01 ||r(i)||/||b|| 5.251890814123e-01
> 2 KSP unpreconditioned resid norm 5.897341285066e-02 true resid norm
> 5.897341285066e-02 ||r(i)||/||b|| 2.541993500758e-01
> 3 KSP unpreconditioned resid norm 2.064441965393e-02 true resid norm
> 2.064441965393e-02 ||r(i)||/||b|| 8.898582946200e-02
> 4 KSP unpreconditioned resid norm 8.019849301680e-17 true resid norm
> 5.101974804388e-17 ||r(i)||/||b|| 2.199158259100e-16
> dot 0.0188443
>
> It seems to be caused by this:
>
> #define KSP_RemoveNullSpace(ksp,y) ((ksp->nullsp && ksp->pc_side == PC_LEFT)
> ? MatNullSpaceRemove(ksp->nullsp,y,PETSC_NULL) : 0)