<div dir="ltr">Hi Matt and Barry,<div><br></div><div>thanks for the good ideas.</div><div><br></div><div>I wasn't aware of the function KSPSetConvergenceTest(), thanks for alerting me to this Matt. Unfortunately I really do need to exclude a couple of specific degrees of freedom from the norm, and as such I want to make sure I can replicate the norm as computed within UpdateHessenberg() before I replace this.</div><div><br></div><div>Apologies Barry, there was a stupid bug in that code, it should have read:</div><div><font face="courier new, monospace">    KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);<br></font></div><div><font face="courier new, monospace">    KSPFGMRESBuildSoln(fgmres->nrs,ksp->vec_sol,tmp_sol,ksp,loc_it+1);<br>    VecAXPY(tmp_sol,-1.0,ksp->vec_rhs);<br>    VecNorm(tmp_sol,NORM_2,&res_tmp);</font><br></div><div><font face="arial, sans-serif"><br></font></div><div><font face="arial, sans-serif">Which does not give the same result as the norm for UpdateHessenberg(). Unfortunately I want to avoid using</font></div><div><font face="courier new, monospace">KSPBuildResidual(ksp,NULL,NULL,&resid);</font><font face="arial, sans-serif"><br></font></div><div>as this calls <font face="courier new, monospace">KSP_MatMult()</font>, which will in turn call my snes residual assembly routine. This routine in turn calls a Navier Stokes solver to run for a fixed (long) time period, and so is expensive.</div><div><br></div><div>Replacing <font face="courier new, monospace">VecNorm()</font> with my own version is indeed what I'm doing elsewhere, I'm just not sure of how to replace the UpdateHessenberg() norm right now.</div><div><br></div><div>Cheers, Dave.</div></div><br><div class="gmail_quote"><div dir="ltr" class="gmail_attr">On Fri, Jul 5, 2019 at 12:59 PM Smith, Barry F. <<a href="mailto:bsmith@mcs.anl.gov">bsmith@mcs.anl.gov</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left:1px solid rgb(204,204,204);padding-left:1ex"><br>
<br>
> On Jul 4, 2019, at 7:29 PM, Dave Lee via petsc-users <<a href="mailto:petsc-users@mcs.anl.gov" target="_blank">petsc-users@mcs.anl.gov</a>> wrote:<br>
> <br>
> Hi PETSc,<br>
> <br>
> 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.<br>
> <br>
> Currently this is done within:<br>
> KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm)<br>
> as:<br>
> *res = PetscAbsScalar(*RS(it+1));<br>
> <br>
> Firstly, I would like to make sure I can replicate the existing residual norm. I tried to do this as:<br>
> KSPFGMRESUpdateHessenberg(ksp,loc_it,hapend,&res_norm);<br>
> KSPFGMRESBuildSoln(fgmres->nrs,ksp >vec_sol,tmp_sol,ksp,loc_it);<br>
<br>
   I don't understand the lines of code below. I don't see how they could compute the norm of the residual. <br>
<br>
>From KSPMonitorTrueResidualNorm() you can compute the unpreconditioned residual with <br>
<br>
  ierr = KSPBuildResidual(ksp,NULL,NULL,&resid);CHKERRQ(ierr);<br>
  ierr = VecNorm(resid,NORM_2,&truenorm);CHKERRQ(ierr);<br>
  ierr = VecDestroy(&resid);CHKERRQ(ierr);<br>
<br>
since FGMRES always uses right preconditioning (and hence the unpreconditioned residual) this should match the value computed internally by FGMRES (subject to roundoff differences)<br>
<br>
To exclude certain values in the residual from the norm calculation you could write a custom "VecNorm" norm that skipped those entries.<br>
<br>
  Barry<br>
<br>
<br>
<br>
<br>
<br>
> VecNorm(tmp_sol,NORM_2,&res_tmp);<br>
> VecNorm(ksp->vec_rhs,NORM_2,&res_rhs);<br>
> res_norm = fabs(res_tmp-res_rhs);<br>
> <br>
> But this doesn't match the norm that comes out of UpdateHessenberg.<br>
> <br>
> Any ideas on how I can re-create the norm derived from UpdateHessenberg as an expansion over the Kyrlov vectors in FGMRES?<br>
> <br>
> Cheers, Dave.<br>
<br>
</blockquote></div>