<div dir="ltr"><div class="gmail_extra"><div class="gmail_quote">On Fri, Nov 7, 2014 at 11:04 AM, Lawrence Mitchell <span dir="ltr"><<a href="mailto:lawrence.mitchell@imperial.ac.uk" target="_blank">lawrence.mitchell@imperial.ac.uk</a>></span> wrote:<br><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">Hi petsc-dev,<br>
<br>
I'm solving a pure Neumann mixed Poisson-like problem, preconditioning with a schur complement.  The pressure space has a nullspace of the constant functions and so I attach the appropriate nullspace to the krylov solver, and compose the constant nullspace with the IS defining the pressure space.  My RHS is consistent.<br></blockquote><div><br></div><div>That is supposed to work, and I think it does in my tests. The code is here</div><div><br></div><div>  <a href="https://bitbucket.org/petsc/petsc/src/1f0d623c8336219eb98f7ded6f95c151ca603fe7/src/ksp/pc/impls/fieldsplit/fieldsplit.c?at=master#cl-562">https://bitbucket.org/petsc/petsc/src/1f0d623c8336219eb98f7ded6f95c151ca603fe7/src/ksp/pc/impls/fieldsplit/fieldsplit.c?at=master#cl-562</a></div><div><br></div><div>so maybe we can track this down in the debugger.</div><div><br></div><div>  Thanks,</div><div><br></div><div>    Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0px 0px 0px 0.8ex;border-left-width:1px;border-left-color:rgb(204,204,204);border-left-style:solid;padding-left:1ex">
When I precondition with:<br>
<br>
-pc_type fieldsplit -pc_fieldsplit_type schur -pc_fieldsplit_schur_fact_type full \<br>
-pc_fieldsplit_schur_precondition selfp<br>
<br>
I notice that the nullspace is not transferred over to the preconditioning matrix for S.  Is this a deliberate choice?<br>
<br>
ksp_view output below, note that the schurcomplement mat has an attached nullspace, but the generated pmat does not.<br>
<br>
Cheers,<br>
<br>
Lawrence<br>
<br>
KSP Object: 1 MPI processes<br>
  type: gmres<br>
    GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br>
    GMRES: happy breakdown tolerance 1e-30<br>
  maximum iterations=30, initial guess is zero<br>
  tolerances:  relative=1e-07, absolute=1e-50, divergence=10000<br>
  left preconditioning<br>
  has attached null space<br>
  using PRECONDITIONED norm type for convergence test<br>
PC Object: 1 MPI processes<br>
  type: fieldsplit<br>
    FieldSplit with Schur preconditioner, factorization FULL<br>
    Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses (the lumped) A00's diagonal's inverse<br>
    Split info:<br>
    Split number 0 Defined by IS<br>
    Split number 1 Defined by IS<br>
    KSP solver for A00 block<br>
      KSP Object:      (fieldsplit_0_)       1 MPI processes<br>
        type: gmres<br>
          GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br>
          GMRES: happy breakdown tolerance 1e-30<br>
        maximum iterations=10000, initial guess is zero<br>
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000<br>
        left preconditioning<br>
        using PRECONDITIONED norm type for convergence test<br>
      PC Object:      (fieldsplit_0_)       1 MPI processes<br>
        type: ilu<br>
          ILU: out-of-place factorization<br>
          0 levels of fill<br>
          tolerance for zero pivot 2.22045e-14<br>
          matrix ordering: natural<br>
          factor fill ratio given 1, needed 1<br>
            Factored matrix follows:<br>
              Mat Object:               1 MPI processes<br>
                type: seqaij<br>
                rows=72, cols=72<br>
                package used to perform factorization: petsc<br>
                total: nonzeros=1080, allocated nonzeros=1080<br>
                total number of mallocs used during MatSetValues calls =0<br>
                  using I-node routines: found 23 nodes, limit used is 5<br>
        linear system matrix = precond matrix:<br>
        Mat Object:        (fieldsplit_0_)         1 MPI processes<br>
          type: seqaij<br>
          rows=72, cols=72<br>
          total: nonzeros=1080, allocated nonzeros=0<br>
          total number of mallocs used during MatSetValues calls =0<br>
            using I-node routines: found 23 nodes, limit used is 5<br>
    KSP solver for S = A11 - A10 inv(A00) A01<br>
      KSP Object:      (fieldsplit_1_)       1 MPI processes<br>
        type: gmres<br>
          GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br>
          GMRES: happy breakdown tolerance 1e-30<br>
        maximum iterations=10000, initial guess is zero<br>
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000<br>
        left preconditioning<br>
        using PRECONDITIONED norm type for convergence test<br>
      PC Object:      (fieldsplit_1_)       1 MPI processes<br>
        type: ilu<br>
          ILU: out-of-place factorization<br>
          0 levels of fill<br>
          tolerance for zero pivot 2.22045e-14<br>
          matrix ordering: natural<br>
          factor fill ratio given 1, needed 1<br>
            Factored matrix follows:<br>
              Mat Object:               1 MPI processes<br>
                type: seqaij<br>
                rows=24, cols=24<br>
                package used to perform factorization: petsc<br>
                total: nonzeros=216, allocated nonzeros=216<br>
                total number of mallocs used during MatSetValues calls =0<br>
                  using I-node routines: found 8 nodes, limit used is 5<br>
        linear system matrix followed by preconditioner matrix:<br>
        Mat Object:        (fieldsplit_1_)         1 MPI processes<br>
          type: schurcomplement<br>
          rows=24, cols=24<br>
            has attached null space<br>
            Schur complement A11 - A10 inv(A00) A01<br>
            A11<br>
              Mat Object:              (fieldsplit_1_)               1 MPI processes<br>
                type: seqaij<br>
                rows=24, cols=24<br>
                total: nonzeros=72, allocated nonzeros=0<br>
                total number of mallocs used during MatSetValues calls =0<br>
                  has attached null space<br>
                  using I-node routines: found 8 nodes, limit used is 5<br>
            A10<br>
              Mat Object:               1 MPI processes<br>
                type: seqaij<br>
                rows=24, cols=72<br>
                total: nonzeros=288, allocated nonzeros=0<br>
                total number of mallocs used during MatSetValues calls =0<br>
                  using I-node routines: found 8 nodes, limit used is 5<br>
            KSP of A00<br>
              KSP Object:              (fieldsplit_0_)               1 MPI processes<br>
                type: gmres<br>
                  GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement<br>
                  GMRES: happy breakdown tolerance 1e-30<br>
                maximum iterations=10000, initial guess is zero<br>
                tolerances:  relative=1e-05, absolute=1e-50, divergence=10000<br>
                left preconditioning<br>
                using PRECONDITIONED norm type for convergence test<br>
              PC Object:              (fieldsplit_0_)               1 MPI processes<br>
                type: ilu<br>
                  ILU: out-of-place factorization<br>
                  0 levels of fill<br>
                  tolerance for zero pivot 2.22045e-14<br>
                  matrix ordering: natural<br>
                  factor fill ratio given 1, needed 1<br>
                    Factored matrix follows:<br>
                      Mat Object:                       1 MPI processes<br>
                        type: seqaij<br>
                        rows=72, cols=72<br>
                        package used to perform factorization: petsc<br>
                        total: nonzeros=1080, allocated nonzeros=1080<br>
                        total number of mallocs used during MatSetValues calls =0<br>
                          using I-node routines: found 23 nodes, limit used is 5<br>
                linear system matrix = precond matrix:<br>
                Mat Object:                (fieldsplit_0_)                 1 MPI processes<br>
                  type: seqaij<br>
                  rows=72, cols=72<br>
                  total: nonzeros=1080, allocated nonzeros=0<br>
                  total number of mallocs used during MatSetValues calls =0<br>
                    using I-node routines: found 23 nodes, limit used is 5<br>
            A01<br>
              Mat Object:               1 MPI processes<br>
                type: seqaij<br>
                rows=72, cols=24<br>
                total: nonzeros=288, allocated nonzeros=0<br>
                total number of mallocs used during MatSetValues calls =0<br>
                  using I-node routines: found 23 nodes, limit used is 5<br>
        Mat Object:         1 MPI processes<br>
          type: seqaij<br>
          rows=24, cols=24<br>
          total: nonzeros=216, allocated nonzeros=216<br>
          total number of mallocs used during MatSetValues calls =0<br>
            using I-node routines: found 8 nodes, limit used is 5<br>
  linear system matrix = precond matrix:<br>
  Mat Object:   1 MPI processes<br>
    type: nest<br>
    rows=96, cols=96<br>
      has attached null space<br>
      Matrix object:<br>
        type=nest, rows=2, cols=2<br>
        MatNest structure:<br>
        (0,0) : prefix="fieldsplit_0_", type=seqaij, rows=72, cols=72<br>
        (0,1) : type=seqaij, rows=72, cols=24<br>
        (1,0) : type=seqaij, rows=24, cols=72<br>
        (1,1) : prefix="fieldsplit_1_", type=seqaij, rows=24, cols=24<br>
<br>
</blockquote></div><br><br clear="all"><div><br></div>-- <br><div class="gmail_signature">What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div>
</div></div>