[petsc-users] pcfieldsplitting not working for VI

Justin Chang jychang48 at gmail.com
Mon Nov 28 20:02:58 CST 2016


Hi all,

We are working on enforcing maximum principles for the Darcy equation using
the variational inequality (with Firedrake). We are employing the VMS
finite element formulation for the Darcy equation which results in a
saddle-point system, so we are using the schur complement approach to solve
the saddle-point system. These are the petsc4py options that we've set:

petsc_options = PETSc.Options()
petsc_options.prefixPush("rs_")
petsc_options.setValue("ksp_type", "gmres")
petsc_options.setValue("pc_type", "fieldsplit")
petsc_options.setValue("pc_fieldsplit_type", "schur")
petsc_options.setValue("pc_fieldsplit_schur_fact_type","full")
petsc_options.setValue("pc_fieldsplit_schur_precondition", "selfp")
petsc_options.setValue("fieldsplit_0_ksp_type", "preonly")
petsc_options.setValue("fieldsplit_0_pc_type", "bjacobi")
petsc_options.setValue("fieldsplit_1_ksp_type", "preonly")
petsc_options.setValue("fieldsplit_1_pc_type", "hypre")
petsc_options.setValue("fieldsplit_1_pc_hypre_type", "boomeramg")
petsc_options.setValue("fieldsplit_1_pc_hypre_boomeramg_strong_threshold",
0.75)
petsc_options.setValue("fieldsplit_1_pc_hypre_boomeramg_agg_nl", 2)
petsc_options.prefixPop()

Attached is the full Firedrake project file (as well as the corresponding
GMSH file in case any Firedrakers out there want to reproduce this) so that
you guys have a basic understanding of what we're trying to do.

If we do not apply any bounds and do a regular Newton method (with
"-snes_type ksponly"), our boundary value problem works. However, if we
apply VINEWTONRSLS we get the following error:

$ python sphere.py -rs_snes_monitor -rs_ksp_monitor
  0 SNES Function norm 3.862953179603e-01
Traceback (most recent call last):
  File "sphere.py", line 179, in <module>
    virs()
  File "sphere.py", line 162, in virs
    rs_solver.solve(None,sol_vec)
  File "PETSc/SNES.pyx", line 537, in petsc4py.PETSc.SNES.solve
(src/petsc4py.PETSc.c:169159)
petsc4py.PETSc.Error: error code 77
[0] SNESSolve() line 4061 in /tmp/pip-plvgeV-build/src/snes/interface/snes.c
[0] SNESSolve_VINEWTONRSLS() line 506 in
/tmp/pip-plvgeV-build/src/snes/impls/vi/rs/virs.c
[0] KSPSetUp() line 393 in
/tmp/pip-plvgeV-build/src/ksp/ksp/interface/itfunc.c
[0] PCSetUp() line 968 in
/tmp/pip-plvgeV-build/src/ksp/pc/interface/precon.c
[0] PCSetUp_FieldSplit() line 487 in
/tmp/pip-plvgeV-build/src/ksp/pc/impls/fieldsplit/fieldsplit.c
[0] PCFieldSplitSetDefaults() line 470 in
/tmp/pip-plvgeV-build/src/ksp/pc/impls/fieldsplit/fieldsplit.c
[0] Petsc has generated inconsistent data
[0] Unhandled case, must have at least two fields, not 1

Why is this happening? It seems to me VINEWTONRSLS is not understanding
that I have a two-field formulation.

Thanks,
Justin
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20161128/cc829c24/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: sphere.py
Type: text/x-python
Size: 5644 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20161128/cc829c24/attachment-0001.py>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: try.msh
Type: model/mesh
Size: 632781 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20161128/cc829c24/attachment-0001.msh>


More information about the petsc-users mailing list