[petsc-users] Residual at the fieldsplit level

Barry Smith bsmith at petsc.dev
Tue Sep 17 12:03:54 CDT 2024


static PetscErrorCode KSPFGMRESCycle(PetscInt *itcount, KSP ksp)
{
...
  PetscCall(VecNorm(VEC_VV(0), NORM_2, &res_norm));


  /* The first entry in the right-hand side of the Hessenberg system is just
     the initial residual norm */
  *RS(0) = res_norm;

  ksp->rnorm = res_norm;

...

  /* scale VEC_VV (the initial residual) */
  PetscCall(VecScale(VEC_VV(0), 1.0 / res_norm));


So, it is expected behavior. GMRES scales the initial residual so it has norm 1, chugs along doing its thing and then when it builds the solution at the end, it knows it has solved the scaled problem so unscales it in KSPFGMRESBuildSoln()  with PetscCall(VecMAXPBY(VEC_TEMP, it + 1, nrs, 0, &PREVEC(0))); the first entry of nrs[0] contains what was set in *RS(0) = res_norm;

The KSPMonitor on the outer GMRES knows that the GMRES is solving a scaled right-hand side, so it still prints the correct original unscaled residual at each iteration.

It is a peculiar feature of the GMES family, the subproblems that are solved are scaled from the original problem. Presumably, this is also true for block Jacobi etc.




> On Sep 17, 2024, at 11:48 AM, miguel.salazar <miguel.salazar at corintis.com> wrote:
> 
> Is this the expected behavior though? I am not sure if it might be an issue with the Firedrake implementation. Getting a minimum example on PETSc like this one is not trivial
>  
> Miguel
>  
> 
> 
> 
> MIGUEL ANGEL SALAZAR DE TROYA
> Head of Software Engineering
> miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>
> Corintis SA
> EPFL Innovation Park Building C
> 1015 Lausanne
> <2024-08-1609_02_11-re_mailmigrationfromgoogletooffice365-sebastien.gobel at corintis.com-co_4456a1bf-dd3e-46b7-85f9-21f278e66a79.png>
> Here at Corintis we care for your privacy. That is why we have taken appropriate measures to ensure that the data you have provided to us is always secure.
> From: Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>>
> Date: Tuesday, 17 September 2024 at 16:58
> To: miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>>
> Cc: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov> <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>>
> Subject: Re: [petsc-users] Residual at the fieldsplit level
> 
>  
>    Ok. Please ignore my last message about fieldsplit_0_fieldsplit_0 I misunderstood 
>  
> 
> 
> On Sep 17, 2024, at 10:52 AM, miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>> wrote:
>  
> That is correct
>  
>  
> MIGUEL ANGEL SALAZAR DE TROYA
> Head of Software Engineering
> miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>
> Corintis SA
> EPFL Innovation Park Building C
> 1015 Lausanne
> 
> <2024-08-1609_02_11-re_mailmigrationfromgoogletooffice365-sebastien.gobel at corintis.com-co_4456a1bf-dd3e-46b7-85f9-21f278e66a79.png>
> Here at Corintis we care for your privacy. That is why we have taken appropriate measures to ensure that the data you have provided to us is always secure.
> From: Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>>
> Date: Tuesday, 17 September 2024 at 16:49
> To: miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>>
> Cc: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov> <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>>
> Subject: Re: [petsc-users] Residual at the fieldsplit level
> 
>  
>    So firedrake_0 is the object prefix of your SNES object?
>  
>  
> 
> On Sep 17, 2024, at 10:37 AM, miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>> wrote:
>  
>  >In my reading, it says 
> >Residual norms for firedrake_0_ solve.
> >    0 KSP Residual norm 4.621368953027e+02
>  
> >     so the residual in the first (0th) split  does correspond to the outer-most KSP residual norm (not printed but equal to the SNES at the first iteration of Newton.
>  
> >     It is the residual norm of the first field of the first field that is 1
>  
> >Residual norms for firedrake_0_fieldsplit_0_ solve.
>  >     0 KSP Residual norm 1.000000000000e+00
>  
> Residual norms for firedrake_0_ solve.
>    0 KSP Residual norm 4.621368953027e+02
>  
> Is the residual of the outer solve (it is being printed).
>  
> Residual norms for firedrake_0_fieldsplit_0_ solve.
>    0 KSP Residual norm 1.000000000000e+00
>  
> Is the residual of the first field. I am not sure what you are referring to with “the residual norm of the first field of the first field”
>  
> Thanks,
> Miguel
>  
>  
>  
>  
> MIGUEL ANGEL SALAZAR DE TROYA
> Head of Software Engineering
> miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>
> Corintis SA
> EPFL Innovation Park Building C
> 1015 Lausanne
> 
> <2024-08-1609_02_11-re_mailmigrationfromgoogletooffice365-sebastien.gobel at corintis.com-co_4456a1bf-dd3e-46b7-85f9-21f278e66a79.png>
> Here at Corintis we care for your privacy. That is why we have taken appropriate measures to ensure that the data you have provided to us is always secure.
> From: Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>>
> Date: Tuesday, 17 September 2024 at 16:24
> To: miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>>
> Cc: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov> <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>>
> Subject: Re: [petsc-users] Residual at the fieldsplit level
> 
>  
> Initial residual: 462.13689530272404
>   0 SNES Function norm 4.621368953027e+02
>     Residual norms for firedrake_0_ solve.
>     0 KSP Residual norm 4.621368953027e+02
>       Residual norms for firedrake_0_fieldsplit_0_ solve.
>       0 KSP Residual norm 1.000000000000e+00
>       1 KSP Residual norm 3.501082228626e-15
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>       0 KSP Residual norm 0.000000000000e+00
>       1 KSP Residual norm 0.000000000000e+00
>     1 KSP Residual norm 1.612167203819e-12
>   1 SNES Function norm 1.599350481360e-12
> ```
>  
> Using the fieldsplit additive preconditioner, the problem converges in a single KSP iteration, as expected. However, I do not understand why the residual of fieldsplit_0 (1e+0) does not coincide with the outer residual (462.13689530272404). It should be the case given that only fieldsplit_0 has a non-zero residual contribution. The fact that it is just 1 is suspicious. Is there something about how the fieldsplit works that I am missing?
>  
>      In my reading, it says 
> Residual norms for firedrake_0_ solve.
>     0 KSP Residual norm 4.621368953027e+02
>  
>      so the residual in the first (0th) split  does correspond to the outer-most KSP residual norm (not printed but equal to the SNES at the first iteration of Newton.
>  
>      It is the residual norm of the first field of the first field that is 1
>  
> Residual norms for firedrake_0_fieldsplit_0_ solve.
>       0 KSP Residual norm 1.000000000000e+00
>  
>    Or am I missing something
>  
>    Barry
>  
>  
>  
>  
> 
>  
>  
> Initial residual: 462.13689530272404
>   0 SNES Function norm 4.621368953027e+02
>     Residual norms for firedrake_0_ solve.
>     0 KSP Residual norm 4.621368953027e+02
>       Residual norms for firedrake_0_ solve.
>       0 KSP unpreconditioned resid norm 4.621368953027e+02 true resid norm 4.621368953027e+02 ||r(i)||/||b|| 1.000000000000e+00
>         Residual norms for firedrake_0_fieldsplit_0_ solve.
>         0 KSP Residual norm 1.000000000000e+00
>         Residual norms for firedrake_0_fieldsplit_0_ solve.
>         0 KSP none resid norm 1.000000000000e+00 true resid norm 3.501082228626e-15 ||r(i)||/||b|| 3.501082228626e-15
>         1 KSP Residual norm 3.501082228626e-15
>       1 KSP none resid norm 3.501082228626e-15 true resid norm 3.501082228626e-15 ||r(i)||/||b|| 3.501082228626e-15
>         Residual norms for firedrake_0_fieldsplit_1_ solve.
>         0 KSP Residual norm 0.000000000000e+00
>         Residual norms for firedrake_0_fieldsplit_1_ solve.
>         0 KSP none resid norm 0.000000000000e+00 true resid norm 0.000000000000e+00 ||r(i)||/||b|| inf
>         1 KSP Residual norm 0.000000000000e+00
>       1 KSP none resid norm 0.000000000000e+00 true resid norm 0.000000000000e+00 ||r(i)||/||b|| inf
>       1 KSP Residual norm 1.612167203819e-12
>     1 KSP unpreconditioned resid norm 1.612167203819e-12 true resid norm 1.589286585800e-12 ||r(i)||/||b|| 3.438995245681e-15
>   1 SNES Function norm 1.599350481360e-12
>  
>  
> 
> On Sep 17, 2024, at 1:04 AM, miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>> wrote:
>  
> These are the new options (with ksp_monitor_true_residual)
> 
> solver_mumps_assembled = {
>     "ksp_type": "preonly",
>     "ksp_monitor": None,
>     "ksp_monitor_true_residual": None,
>     "pc_type": "python",
>     "pc_python_type": "firedrake.AssembledPC",
>     "assembled_pc_type": "lu",
>     "assembled_pc_factor_mat_solver_type": "mumps",
>     "assembled_mat_mumps_icntl_14": 200,
>     "assembled_mat_mumps_icntl_24": 1,
> }
> solver_fieldsplit = {
>     "mat_type": "matfree",
>     "snes_monitor": None,
>     "snes_type": "newtonls",
>     "ksp_type": "fgmres",
>     "ksp_monitor_true_residual": None,
>     "ksp_rtol": 1e-1,
>     "ksp_monitor": None,
>     "pc_type": "fieldsplit",
>     "pc_fieldsplit_type": "additive",
>     "fieldsplit_0": solver_mumps_assembled,
>     "fieldsplit_1": solver_mumps_assembled,
> }
>  
> And this is the output
>  
> Initial residual: 462.13689530272404
>   0 SNES Function norm 4.621368953027e+02
>     Residual norms for firedrake_0_ solve.
>     0 KSP Residual norm 4.621368953027e+02
>     Residual norms for firedrake_0_ solve.
>     0 KSP unpreconditioned resid norm 4.621368953027e+02 true resid norm 4.621368953027e+02 ||r(i)||/||b|| 1.000000000000e+00
>       Residual norms for firedrake_0_fieldsplit_0_ solve.
>       0 KSP Residual norm 1.000000000000e+00
>       Residual norms for firedrake_0_fieldsplit_0_ solve.
>       0 KSP none resid norm 1.000000000000e+00 true resid norm 3.501082228626e-15 ||r(i)||/||b|| 3.501082228626e-15
>       1 KSP Residual norm 3.501082228626e-15
>       1 KSP none resid norm 3.501082228626e-15 true resid norm 3.501082228626e-15 ||r(i)||/||b|| 3.501082228626e-15
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>       0 KSP Residual norm 0.000000000000e+00
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>       0 KSP none resid norm 0.000000000000e+00 true resid norm 0.000000000000e+00 ||r(i)||/||b|| inf
>       1 KSP Residual norm 0.000000000000e+00
>       1 KSP none resid norm 0.000000000000e+00 true resid norm 0.000000000000e+00 ||r(i)||/||b|| inf
>     1 KSP Residual norm 1.612167203819e-12
>     1 KSP unpreconditioned resid norm 1.612167203819e-12 true resid norm 1.589286585800e-12 ||r(i)||/||b|| 3.438995245681e-15
>   1 SNES Function norm 1.599350481360e-12
>  
> The true residual is very low. (3.501082228626e-15)
>  
>  
> MIGUEL ANGEL SALAZAR DE TROYA
> Head of Software Engineering
> miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>
> Corintis SA
> EPFL Innovation Park Building C
> 1015 Lausanne
> 
> <2024-08-1609_02_11-re_mailmigrationfromgoogletooffice365-sebastien.gobel at corintis.com-co_4456a1bf-dd3e-46b7-85f9-21f278e66a79.png>
> Here at Corintis we care for your privacy. That is why we have taken appropriate measures to ensure that the data you have provided to us is always secure.
> From: Barry Smith <bsmith at petsc.dev <mailto:bsmith at petsc.dev>>
> Date: Monday, 16 September 2024 at 19:34
> To: miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>>
> Cc: petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov> <petsc-users at mcs.anl.gov <mailto:petsc-users at mcs.anl.gov>>
> Subject: Re: [petsc-users] Residual at the fieldsplit level
> 
>  
>   Try ksp_monitor_true_residual. Let's see if it is an issue of preconditioned vs unpreconditioned residual.
>  
>  
>  
> 
> On Sep 16, 2024, at 1:05 PM, miguel.salazar <miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>> wrote:
>  
> Hello,
>  
> I have this simple example on Firedrake to illustrate my point.  I am solving for a two-component poisson equation (uncoupled). Only the first component has a non-zero residual.
> 
> ```
> import firedrake as fd
> from firedrake import inner, grad, dx, sin, pi
>  
> N = 10
> mesh = fd.UnitSquareMesh(N, N)
> V = fd.FunctionSpace(mesh, "CG", 1)
> W = V * V
> u = fd.Function(W)
> v = fd.TestFunction(W)
> a = inner(grad(u[0]), grad(v[0])) * dx + inner(grad(u[1]), grad(v[1])) * dx
> x = fd.SpatialCoordinate(mesh)
> f = fd.Function(V)
> f.interpolate(fd.Constant(1e4) * sin(x[0] * pi) * sin(2 * x[1] * pi))
> L = f * v[0] * dx
> F = a - L
> bcs = [fd.DirichletBC(W.sub(0), fd.Constant(2.0), (1,))]
>  
>  
> def snes_firedrake_residual(F, u, bcs):
>     for bcs_ in bcs:
>         bcs_.apply(u)
>     residual = fd.assemble(F, bcs=bcs, zero_bc_nodes=True)
>     with residual.dat.vec_ro as r:
>         print("Initial residual:", r.norm())
>  
>  
> snes_firedrake_residual(F, u, bcs)
>  
>  
> problem = fd.NonlinearVariationalProblem(F, u, bcs=bcs)
> solver_mumps_assembled = {
>     "ksp_type": "preonly",
>     "ksp_monitor": None,
>     "pc_type": "python",
>     "pc_python_type": "firedrake.AssembledPC",
>     "assembled_pc_type": "lu",
>     "assembled_pc_factor_mat_solver_type": "mumps",
>     "assembled_mat_mumps_icntl_14": 200,
>     "assembled_mat_mumps_icntl_24": 1,
> }
> solver_fieldsplit = {
>     "mat_type": "matfree",
>     "snes_type": "newtonls",
>     "ksp_type": "fgmres",
>     "ksp_rtol": 1e-1,
>     "ksp_monitor": None,
>     "pc_type": "fieldsplit",
>     "pc_fieldsplit_type": "additive",
>     "fieldsplit_0": solver_mumps_assembled,
>     "fieldsplit_1": solver_mumps_assembled,
> }
>  
> solver = fd.NonlinearVariationalSolver(problem, solver_parameters=solver_fieldsplit)
> solver.solve()
> ```
> 
> The PETSc output is as follow
> 
> ```
> Initial residual: 462.13689530272404
>   0 SNES Function norm 4.621368953027e+02
>     Residual norms for firedrake_0_ solve.
>     0 KSP Residual norm 4.621368953027e+02
>       Residual norms for firedrake_0_fieldsplit_0_ solve.
>       0 KSP Residual norm 1.000000000000e+00
>       1 KSP Residual norm 3.501082228626e-15
>       Residual norms for firedrake_0_fieldsplit_1_ solve.
>       0 KSP Residual norm 0.000000000000e+00
>       1 KSP Residual norm 0.000000000000e+00
>     1 KSP Residual norm 1.612167203819e-12
>   1 SNES Function norm 1.599350481360e-12
> ```
>  
> Using the fieldsplit additive preconditioner, the problem converges in a single KSP iteration, as expected. However, I do not understand why the residual of fieldsplit_0 (1e+0) does not coincide with the outer residual (462.13689530272404). It should be the case given that only fieldsplit_0 has a non-zero residual contribution. The fact that it is just 1 is suspicious. Is there something about how the fieldsplit works that I am missing?
> 
> Thanks,
> Miguel
>  
> MIGUEL ANGEL SALAZAR DE TROYA
> Head of Software Engineering
> miguel.salazar at corintis.com <mailto:miguel.salazar at corintis.com>
> Corintis SA
> EPFL Innovation Park Building C
> 1015 Lausanne
> 
> <2024-08-1609_02_11-re_mailmigrationfromgoogletooffice365-sebastien.gobel at corintis.com-co_4456a1bf-dd3e-46b7-85f9-21f278e66a79.png>
> Here at Corintis we care for your privacy. That is why we have taken appropriate measures to ensure that the data you have provided to us is always secure.

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20240917/50fcc374/attachment-0001.html>


More information about the petsc-users mailing list