[petsc-users] neumann failure in my version of snes ex12

Geoffrey Irving irving at naml.us
Fri Nov 22 17:36:05 CST 2013


I have a duplicate of snes ex12 (FEM Poisson) which works with
Dirichlet boundary conditions, but it's breaking for me with Neumann
conditions.  In particular, with Neumann conditions I get results
which explode even though I believe I am setting a constant nullspace.

For example, if I use two first order elements (the unit square
divided into two triangles), the resulting solution has

    L2 error = 1.75514e+08
    u = [-175513825.75680602, -175513825.66302037,
-175513825.48390722, -175513824.84436429]

This looks rather a lot like the null space isn't getting through.  I
am creating the constant nullspace with

      MatNullSpace null;
      CHECK(MatNullSpaceCreate(comm(),PETSC_TRUE,0,0,&null));
      CHECK(MatSetNullSpace(m,null));
      CHECK(MatNullSpaceDestroy(&null));

If I pass "-ksp_view -mat_view", I get the following.  The matrix
entries seem right (they do indeed have the constant nullspace), and
ksp_view shows that a nullspace is attached.  Is attaching the
nullspace to the matrix with MatSetNullSpace enough, or do I need to
additionally attach it to the KSP object?

Thanks,
Geoffrey

---------------------------------------------

Mat Object: 1 MPI processes
  type: seqaij
row 0: (0, 1)  (1, -0.5)  (2, -0.5)
row 1: (0, -0.5)  (1, 1)  (2, 0)  (3, -0.5)
row 2: (0, -0.5)  (1, 0)  (2, 1)  (3, -0.5)
row 3: (1, -0.5)  (2, -0.5)  (3, 1)
KSP Object: 1 MPI processes
  type: gmres
    GMRES: restart=30, using Classical (unmodified) Gram-Schmidt
Orthogonalization with no iterative refinement
    GMRES: happy breakdown tolerance 1e-30
  maximum iterations=10000, initial guess is zero
  tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
  left preconditioning
  has attached null space
  using PRECONDITIONED norm type for convergence test
PC Object: 1 MPI processes
  type: ilu
    ILU: out-of-place factorization
    0 levels of fill
    tolerance for zero pivot 2.22045e-14
    using diagonal shift on blocks to prevent zero pivot [INBLOCKS]
    matrix ordering: natural
    factor fill ratio given 1, needed 1
      Factored matrix follows:
        Mat Object:         1 MPI processes
          type: seqaij
          rows=4, cols=4
          package used to perform factorization: petsc
          total: nonzeros=14, allocated nonzeros=14
          total number of mallocs used during MatSetValues calls =0
            using I-node routines: found 3 nodes, limit used is 5
  linear system matrix = precond matrix:
  Mat Object:   1 MPI processes
    type: seqaij
    rows=4, cols=4
    total: nonzeros=14, allocated nonzeros=14
    total number of mallocs used during MatSetValues calls =0
      has attached null space
      using I-node routines: found 3 nodes, limit used is 5


More information about the petsc-users mailing list