[petsc-users] Computing residual norm in KSPFGMRESCycle()

Smith, Barry F. bsmith at mcs.anl.gov
Thu Jul 4 21:59:38 CDT 2019



> On Jul 4, 2019, at 7:29 PM, Dave Lee via petsc-users <petsc-users at mcs.anl.gov> wrote:
> 
> Hi PETSc,
> 
> I have a problem for which I need to exclude certain degrees of freedom from the evaluation of the norm as used to monitor the residual in FGMRES (don't ask!). Therefore I would like to replace the existing norm evaluation with my own version.
> 
> Currently this is done within:
> KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm)
> as:
> *res = PetscAbsScalar(*RS(it+1));
> 
> Firstly, I would like to make sure I can replicate the existing residual norm. I tried to do this as:
> KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);
> KSPFGMRESBuildSoln(fgmres->nrs,ksp >vec_sol,tmp_sol,ksp,loc_it);

   I don't understand the lines of code below. I don't see how they could compute the norm of the residual. 

>From KSPMonitorTrueResidualNorm() you can compute the unpreconditioned residual with 

  ierr = KSPBuildResidual(ksp,NULL,NULL,&resid);CHKERRQ(ierr);
  ierr = VecNorm(resid,NORM_2,&truenorm);CHKERRQ(ierr);
  ierr = VecDestroy(&resid);CHKERRQ(ierr);

since FGMRES always uses right preconditioning (and hence the unpreconditioned residual) this should match the value computed internally by FGMRES (subject to roundoff differences)

To exclude certain values in the residual from the norm calculation you could write a custom "VecNorm" norm that skipped those entries.

  Barry





> VecNorm(tmp_sol,NORM_2,&res_tmp);
> VecNorm(ksp->vec_rhs,NORM_2,&res_rhs);
> res_norm = fabs(res_tmp-res_rhs);
> 
> But this doesn't match the norm that comes out of UpdateHessenberg.
> 
> Any ideas on how I can re-create the norm derived from UpdateHessenberg as an expansion over the Kyrlov vectors in FGMRES?
> 
> Cheers, Dave.



More information about the petsc-users mailing list