[petsc-users] PCFIELDSPLIT with MATSBAIJ

Carl-Johan Thore carljohanthore at gmail.com
Wed Sep 6 01:13:53 CDT 2023


Thanks for this!
I now get convergence with sbaij and fieldsplit. However, as can be seen in
the attached file, it requires around three times as many outer iterations
as with aij to reach the same point (to within rtol). Switching from
-pc_fieldsplit_schur_fact_type LOWER to FULL (which is the "correct"
choice?) helps quite a bit, but still sbaij takes many more iterations. So
it seems that your sbaij-implementation works but that I'm doing something
wrong when setting up the fieldsplit pc, or that some of the subsolvers
doesn't work properly with sbaij. I've tried bjacobi as you had in your
logs, but not hypre yet. But anyway, one doesn't expect the matrix format
to have this much impact on convergence? Any other suggestions?
/Carl-Johan

On Mon, Sep 4, 2023 at 9:31 PM Pierre Jolivet <Pierre.Jolivet at lip6.fr>
wrote:

> The branch should now be good to go (
> https://gitlab.com/petsc/petsc/-/merge_requests/6841).
> Sorry, I made a mistake before, hence the error on PetscObjectQuery().
> I’m not sure the code will be covered by the pipeline, but I have tested
> this on a Raviart—Thomas discretization with PCFIELDSPLIT.
> You’ll see in the attached logs that:
> 1) the numerics match
> 2) in the SBAIJ case, PCFIELDSPLIT extract the (non-symmetric) A_{01}
> block from the global (symmetric) A and we get the A_{10} block cheaply by
> just using MatCreateHermitianTranspose(), instead of calling another time
> MatCreateSubMatrix()
> Please let me know if you have some time to test the branch and whether it
> fails or succeeds on your test cases.
>
> Also, I do not agree with what Hong said.
> Sometimes, the assembly of a coefficient can be more expensive than the
> communication of the said coefficient.
> So they are instances where SBAIJ would be more efficient than AIJ even if
> it would require more communication, it is not a black and white picture.
>
> Thanks,
> Pierre
>
>
> On 28 Aug 2023, at 12:12 PM, Pierre Jolivet <Pierre.Jolivet at lip6.fr>
> wrote:
>
>
> On 28 Aug 2023, at 6:50 PM, Carl-Johan Thore <carljohanthore at gmail.com>
> wrote:
>
> I've tried the new files, and with them, PCFIELDSPLIT now gets set up
> without crashes (but the setup is significantly slower than for MATAIJ)
>
>
> I’ll be back from Japan at the end of this week, my schedule is too packed
> to get anything done in the meantime.
> But I’ll let you know when things are working properly (last I checked, I
> think it was working, but I may have forgotten about a corner case or two).
> But yes, though one would except things to be faster and less memory
> intensive with SBAIJ, it’s unfortunately not always the case.
>
> Thanks,
> Pierre
>
> Unfortunately I still get errors later in the process:
>
> [0]PETSC ERROR: Null argument, when expecting valid pointer
> [0]PETSC ERROR: Null Pointer: Parameter # 1
> [0]PETSC ERROR: Petsc Development GIT revision: v3.19.4-1023-ga6d78fcba1d
>  GIT Date: 2023-08-22 20:32:33 -0400
> [0]PETSC ERROR: Configure options -f --with-fortran-bindings=0 --with-cuda
> --with-cusp --download-scalapack --download-hdf5 --download-zlib
> --download-mumps --download-parmetis --download-metis --download-ptscotch
> --download-hypre --download-spai
> [0]PETSC ERROR: #1 PetscObjectQuery() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/sys/objects/inherit.c:742
> [0]PETSC ERROR: #2 MatCreateSubMatrix_MPISBAIJ() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/mat/impls/sbaij/mpi/mpisbaij.c:1414
> [0]PETSC ERROR: #3 MatCreateSubMatrix() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/mat/interface/matrix.c:8476
> [0]PETSC ERROR: #4 PCSetUp_FieldSplit() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/ksp/pc/impls/fieldsplit/fieldsplit.c:826
> [0]PETSC ERROR: #5 PCSetUp() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/ksp/pc/interface/precon.c:1069
> [0]PETSC ERROR: #6 KSPSetUp() at
> /mnt/c/mathware/petsc/petsc-v3-19-4/src/ksp/ksp/interface/itfunc.c:415
>
> The code I'm running here works without any problems for MATAIJ. To run it
> with MATSBAIJ I've simply used the command-line option
> -dm_mat_type sbaij
>
>
> Kind regards,
> Carl-Johan
>
>
> On Sat, Aug 26, 2023 at 5:21 PM Pierre Jolivet via petsc-users <
> petsc-users at mcs.anl.gov> wrote:
>
>>
>>
>> On 27 Aug 2023, at 12:14 AM, Carl-Johan Thore <carl-johan.thore at liu.se>
>> wrote:
>>
>> “Well, your A00 and A11 will possibly be SBAIJ also, so you’ll end up
>> with the same issue.”
>> I’m not sure I follow. Does PCFIELDSPLIT extract further submatrices from
>> these blocks, or is there
>> somewhere else in the code that things will go wrong?
>>
>>
>> Ah, no, you are right, in that case it should work.
>>
>> For the MATNEST I was thinking to get some savings from the
>> block-symmetry at least
>> even if symmetry in A00 and A11 cannot be exploited; using SBAIJ for them
>> would just be a
>> (pretty big) bonus.
>>
>> “I’ll rebase on top of main and try to get it integrated if it could be
>> useful to you (but I’m traveling
>> right now so everything gets done more slowly, sorry).”
>> Sound great, Thanks again!
>>
>>
>> The MR is there https://gitlab.com/petsc/petsc/-/merge_requests/6841.
>> I need to add a new code path in MatCreateRedundantMatrix() to make sure
>> the resulting Mat is indeed SBAIJ, but that is orthogonal to the
>> PCFIELDSPLIT issue.
>> The branch should be usable in its current state.
>>
>> Thanks,
>> Pierre
>>
>>
>> *From:* Pierre Jolivet <pierre.jolivet at lip6.fr>
>> *Sent:* Saturday, August 26, 2023 4:36 PM
>> *To:* Carl-Johan Thore <carl-johan.thore at liu.se>
>> *Cc:* Carl-Johan Thore <carljohanthore at gmail.com>;
>> petsc-users at mcs.anl.gov
>> *Subject:* Re: [petsc-users] PCFIELDSPLIT with MATSBAIJ
>>
>>
>>
>>
>> On 26 Aug 2023, at 11:16 PM, Carl-Johan Thore <carl-johan.thore at liu.se>
>> wrote:
>>
>> "(Sadly) MATSBAIJ is extremely broken, in particular, it cannot be used
>> to retrieve rectangular blocks in MatCreateSubMatrices, thus you cannot get
>> the A01 and A10 blocks in PCFIELDSPLIT.
>> I have a branch that fixes this, but I haven’t rebased in a while (and
>> I’m AFK right now), would you want me to rebase and give it a go, or must
>> you stick to a release tarball?"
>>
>> Ok, would be great if you could look at this! I don't need to stick to
>> any particular branch.
>>
>> Do you think MATNEST could be an alternative here?
>>
>>
>> Well, your A00 and A11 will possibly be SBAIJ also, so you’ll end up with
>> the same issue.
>> I’m using both approaches (monolithic SBAIJ or Nest + SBAIJ), it was
>> crashing but I think it was thoroughly fixed in
>> https://gitlab.com/petsc/petsc/-/commits/jolivet/feature-matcreatesubmatrices-rectangular-sbaij/
>> It is ugly code on top of ugly code, so I didn’t try to get it integrated
>> and just used the branch locally, and then moved to some other stuff.
>> I’ll rebase on top of main and try to get it integrated if it could be
>> useful to you (but I’m traveling right now so everything gets done more
>> slowly, sorry).
>>
>> Thanks,
>> Pierre
>>
>>
>> My matrix is
>> [A00 A01;
>> A01^t A11]
>> so perhaps with MATNEST I can make use of the block-symmetry at least,
>> and then use MATSBAIJ for
>> A00 and A11 if it's possible to combine matrix types which the manual
>> seems to imply.
>>
>> Kind regards
>> Carl-Johan
>>
>>
>>
>> On 26 Aug 2023, at 10:09 PM, Carl-Johan Thore <carljohanthore at gmail.com>
>> wrote:
>> 
>> Hi,
>>
>> I'm trying to use PCFIELDSPLIT with MATSBAIJ in PETSc 3.19.4.
>> According to the manual "[t]he fieldsplit preconditioner cannot
>> currently be used with the MATBAIJ or MATSBAIJ data formats if the
>> blocksize is larger than 1". Since my blocksize is exactly 1 it would
>> seem that I can use PCFIELDSPLIT. But this fails with "PETSC ERROR: For
>> symmetric format, iscol must equal isrow"
>> from MatCreateSubMatrix_MPISBAIJ. Tracing backwards one ends up in
>> fieldsplit.c at
>>
>> /* extract the A01 and A10 matrices */ ilink = jac->head;
>> PetscCall(ISComplement(ilink->is_col, rstart, rend, &ccis)); if
>> (jac->offdiag_use_amat) { PetscCall(MatCreateSubMatrix(pc->mat,
>> ilink->is, ccis, MAT_INITIAL_MATRIX, &jac->B)); } else {
>>        PetscCall(MatCreateSubMatrix(pc->pmat, ilink->is, ccis,
>> MAT_INITIAL_MATRIX, &jac->B)); }
>>
>> This, since my A01 and A10 are not square, seems to explain why iscol is
>> not equal to isrow.
>> From this I gather that it is in general NOT possible to use
>> PCFIELDSPLIT with MATSBAIJ even with block size 1?
>>
>> Kind regards,
>> Carl-Johan
>>
>>
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230906/f09a8a92/attachment-0001.html>
-------------- next part --------------
%%% With AIJ:
...
307 KSP Residual norm 1.854926995470e-09
308 KSP Residual norm 1.751893310256e-09
309 KSP Residual norm 1.646593329093e-09
310 KSP Residual norm 1.543009494478e-09
311 KSP Residual norm 1.452641639917e-09
312 KSP Residual norm 1.369620969855e-09
313 KSP Residual norm 1.302155745365e-09
314 KSP Residual norm 1.252174189783e-09
315 KSP Residual norm 1.192872051154e-09
316 KSP Residual norm 1.138612609555e-09
Linear solve converged due to CONVERGED_RTOL iterations 316


%%% With SBAIJ (More detailed output at the end of this file, also output from factorization FULL)
...
965 KSP Residual norm 1.346789114008e-09
966 KSP Residual norm 1.306538742039e-09
967 KSP Residual norm 1.247901915750e-09
968 KSP Residual norm 1.186932599681e-09
969 KSP Residual norm 1.141186989549e-09
Linear solve converged due to CONVERGED_RTOL iterations 969

%%% KSPVIEW with SBAIJ :
KSP Object: 4 MPI processes
  type: fgmres
    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
    happy breakdown tolerance 1e-30
  maximum iterations=1000, nonzero initial guess
  minimum iterations=3
  tolerances:  relative=1e-10, absolute=1e-50, divergence=1e+09
  right preconditioning
  using UNPRECONDITIONED norm type for convergence test
PC Object: 4 MPI processes
  type: fieldsplit
    FieldSplit with Schur preconditioner, factorization LOWER
    Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's diagonal's inverse
    Split info:
    Split number 0 Defined by IS
    Split number 1 Defined by IS
    KSP solver for A00 block
      KSP Object: (fieldsplit_0_) 4 MPI processes
        type: preonly
        maximum iterations=10000, initial guess is zero
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
        left preconditioning
        using NONE norm type for convergence test
      PC Object: (fieldsplit_0_) 4 MPI processes
        type: asm
          total subdomain blocks = 4, amount of overlap = 1
          restriction/interpolation type - RESTRICT
          Local solver information for first block is in the following KSP and PC objects on rank 0:
          Use -fieldsplit_0_ksp_view ::ascii_info_detail to display information for all blocks
        KSP Object: (fieldsplit_0_sub_) 1 MPI process
          type: preonly
          maximum iterations=10000, initial guess is zero
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
          left preconditioning
          using NONE norm type for convergence test
        PC Object: (fieldsplit_0_sub_) 1 MPI process
          type: icc
            out-of-place factorization
            0 levels of fill
            tolerance for zero pivot 2.22045e-14
            using Manteuffel shift [POSITIVE_DEFINITE]
            matrix ordering: natural
            factor fill ratio given 1., needed 1.
              Factored matrix follows:
                Mat Object: (fieldsplit_0_sub_) 1 MPI process
                  type: seqsbaij
                  rows=97788, cols=97788
                  package used to perform factorization: petsc
                  total: nonzeros=1153455, allocated nonzeros=1153455
                      block size is 1
          linear system matrix = precond matrix:
          Mat Object: (fieldsplit_0_sub_) 1 MPI process
            type: seqsbaij
            rows=97788, cols=97788
            total: nonzeros=1153455, allocated nonzeros=1153455
            total number of mallocs used during MatSetValues calls=0
                block size is 1
        linear system matrix = precond matrix:
        Mat Object: (fieldsplit_0_) 4 MPI processes
          type: mpisbaij
          rows=392931, cols=392931
          total: nonzeros=5739651, allocated nonzeros=5739651
          total number of mallocs used during MatSetValues calls=0
              block size is 1
    KSP solver for S = A11 - A10 inv(A00) A01
      KSP Object: (fieldsplit_1_) 4 MPI processes
        type: preonly
        maximum iterations=10000, initial guess is zero
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
        left preconditioning
        using NONE norm type for convergence test
      PC Object: (fieldsplit_1_) 4 MPI processes
        type: asm
          total subdomain blocks = 4, amount of overlap = 1
          restriction/interpolation type - RESTRICT
          Local solver information for first block is in the following KSP and PC objects on rank 0:
          Use -fieldsplit_1_ksp_view ::ascii_info_detail to display information for all blocks
        KSP Object: (fieldsplit_1_sub_) 1 MPI process
          type: preonly
          maximum iterations=10000, initial guess is zero
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
          left preconditioning
          using NONE norm type for convergence test
        PC Object: (fieldsplit_1_sub_) 1 MPI process
          type: ilu
            out-of-place factorization
            0 levels of fill
            tolerance for zero pivot 2.22045e-14
            matrix ordering: natural
            factor fill ratio given 1., needed 1.
              Factored matrix follows:
                Mat Object: (fieldsplit_1_sub_) 1 MPI process
                  type: seqaij
                  rows=31038, cols=31038
                  package used to perform factorization: petsc
                  total: nonzeros=223972, allocated nonzeros=223972
                    not using I-node routines
          linear system matrix = precond matrix:
          Mat Object: (fieldsplit_1_sub_) 1 MPI process
            type: seqaij
            rows=31038, cols=31038
            total: nonzeros=223972, allocated nonzeros=223972
            total number of mallocs used during MatSetValues calls=0
              not using I-node routines
        linear system matrix followed by preconditioner matrix:
        Mat Object: (fieldsplit_1_) 4 MPI processes
          type: schurcomplement
          rows=121600, cols=121600
            Schur complement A11 - A10 inv(A00) A01
            A11
              Mat Object: (fieldsplit_1_) 4 MPI processes
                type: mpisbaij
                rows=121600, cols=121600
                total: nonzeros=242540, allocated nonzeros=242540
                total number of mallocs used during MatSetValues calls=0
                    block size is 1
            A10
              Mat Object: 4 MPI processes
                type: hermitiantranspose
                rows=121600, cols=392931
            KSP solver for A00 block viewable with the additional option -fieldsplit_0_ksp_view
            A01
              Mat Object: 4 MPI processes
                type: mpiaij
                rows=392931, cols=121600
                total: nonzeros=1018872, allocated nonzeros=1018872
                total number of mallocs used during MatSetValues calls=0
                  using I-node (on process 0) routines: found 23081 nodes, limit used is 5
        Mat Object: 4 MPI processes
          type: mpiaij
          rows=121600, cols=121600
          total: nonzeros=1120336, allocated nonzeros=1120336
          total number of mallocs used during MatSetValues calls=0
            not using I-node (on process 0) routines
  linear system matrix = precond matrix:
  Mat Object: 4 MPI processes
    type: mpisbaij
    rows=514531, cols=514531
    total: nonzeros=7001063, allocated nonzeros=26569279
    total number of mallocs used during MatSetValues calls=0
        block size is 1


%% SBAIJ AGAIN WITH MORE DETAILED OUTPUT:
....
848 KSP Residual norm 7.104369530073e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997543072e-01
    1 KSP Residual norm 3.002328976968e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.004801268598e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.106461054312e-01
    1 KSP Residual norm 9.032329726351e-03
    1 KSP Residual norm 6.971835623451e-04
849 KSP Residual norm 6.950875510226e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997408161e-01
    1 KSP Residual norm 2.988085947538e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.193961257338e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.387629868742e-01
    1 KSP Residual norm 8.835757046845e-03
    1 KSP Residual norm 7.160201171420e-04
850 KSP Residual norm 6.792871940409e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997211093e-01
    1 KSP Residual norm 2.806859529571e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.463113208838e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.463497493063e-01
    1 KSP Residual norm 8.969428389311e-03
    1 KSP Residual norm 7.429684096268e-04
851 KSP Residual norm 6.633256020443e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997044515e-01
    1 KSP Residual norm 2.645025742398e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.681680406391e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.591871186086e-01
    1 KSP Residual norm 9.326676402818e-03
    1 KSP Residual norm 7.649711247555e-04
852 KSP Residual norm 6.488673644160e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997190302e-01
    1 KSP Residual norm 2.614714200275e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.489198714589e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.588925897545e-01
    1 KSP Residual norm 9.242833701551e-03
    1 KSP Residual norm 7.455683026848e-04
853 KSP Residual norm 6.373470520032e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997676684e-01
    1 KSP Residual norm 2.700531175710e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.811087731693e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.402443837288e-01
    1 KSP Residual norm 9.505505737547e-03
    1 KSP Residual norm 6.774446464568e-04
854 KSP Residual norm 6.276004379688e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997984917e-01
    1 KSP Residual norm 2.974092924210e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.339611385305e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.252539481806e-01
    1 KSP Residual norm 9.734258685392e-03
    1 KSP Residual norm 6.302508777492e-04
855 KSP Residual norm 6.179752138080e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997971420e-01
    1 KSP Residual norm 2.993089069965e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.362463802669e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.329578047938e-01
    1 KSP Residual norm 9.380111691525e-03
    1 KSP Residual norm 6.324339399561e-04
856 KSP Residual norm 6.075772917666e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997755989e-01
    1 KSP Residual norm 2.969200049264e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.693003909781e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.418766186278e-01
    1 KSP Residual norm 9.399399487587e-03
    1 KSP Residual norm 6.656554238335e-04
857 KSP Residual norm 5.964043172989e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997513769e-01
    1 KSP Residual norm 2.957682018798e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.044315948253e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.312769868691e-01
    1 KSP Residual norm 9.609354497625e-03
    1 KSP Residual norm 7.011054335511e-04
858 KSP Residual norm 5.836067800148e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997035756e-01
    1 KSP Residual norm 2.872816376259e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.694755749377e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.541154918731e-01
    1 KSP Residual norm 9.433228788725e-03
    1 KSP Residual norm 7.662132564765e-04
859 KSP Residual norm 5.706276444297e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999996843365e-01
    1 KSP Residual norm 2.880361159842e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.938909858461e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.830052055871e-01
    1 KSP Residual norm 9.499862082307e-03
    1 KSP Residual norm 7.906577788230e-04
860 KSP Residual norm 5.617583613817e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997765805e-01
    1 KSP Residual norm 2.777252773447e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.678302112735e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.424164602895e-01
    1 KSP Residual norm 8.849911913621e-03
    1 KSP Residual norm 6.642132850462e-04
861 KSP Residual norm 5.561524925175e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999998577848e-01
    1 KSP Residual norm 2.988351670810e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 5.323587005378e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.912319135071e-01
    1 KSP Residual norm 8.263259855887e-03
    1 KSP Residual norm 5.281070035682e-04
862 KSP Residual norm 5.506766166665e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999998617170e-01
    1 KSP Residual norm 2.845994126117e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 5.250474194834e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.260796113127e-01
    1 KSP Residual norm 9.811742550036e-03
    1 KSP Residual norm 5.203946455451e-04
863 KSP Residual norm 5.439915701058e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999998293776e-01
    1 KSP Residual norm 2.727333121556e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 5.832172367101e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 6.441137952198e-01
    1 KSP Residual norm 1.095694768745e-02
    1 KSP Residual norm 5.783308321627e-04
864 KSP Residual norm 5.403748531336e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999076505e-01
    1 KSP Residual norm 2.946574403519e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 4.284720832619e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.468942171326e-01
    1 KSP Residual norm 1.008737195428e-02
    1 KSP Residual norm 4.226715398658e-04
865 KSP Residual norm 5.390145255217e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999680097e-01
    1 KSP Residual norm 2.682345341857e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.513585515632e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.917467950759e-01
    1 KSP Residual norm 8.647135665273e-03
    1 KSP Residual norm 2.417216343103e-04
866 KSP Residual norm 5.380625495115e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999810804e-01
    1 KSP Residual norm 2.716840513330e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.924843435164e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.107512032702e-01
    1 KSP Residual norm 9.149302558168e-03
    1 KSP Residual norm 1.792387424209e-04
867 KSP Residual norm 5.372335239833e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999876245e-01
    1 KSP Residual norm 2.568578816344e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.558982763978e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.905650299337e-01
    1 KSP Residual norm 8.903714362140e-03
    1 KSP Residual norm 1.400939970262e-04
868 KSP Residual norm 5.360438558152e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999864784e-01
    1 KSP Residual norm 2.699420168709e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.614606040709e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.810945846549e-01
    1 KSP Residual norm 9.180051719635e-03
    1 KSP Residual norm 1.443011750042e-04
869 KSP Residual norm 5.348297140501e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999894059e-01
    1 KSP Residual norm 2.650620998028e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.427955851890e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.386304184444e-01
    1 KSP Residual norm 9.455476649684e-03
    1 KSP Residual norm 1.243623777490e-04
870 KSP Residual norm 5.332375540684e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999907753587e-01
    1 KSP Residual norm 1.397049625888e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 4.294988477916e-03
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.943871296139e-01
    1 KSP Residual norm 1.398563480253e-02
    1 KSP Residual norm 4.293847237178e-03
871 KSP Residual norm 5.321236219875e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999651286e-01
    1 KSP Residual norm 3.081023329909e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.645953802241e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 3.642004434000e-01
    1 KSP Residual norm 7.705417280598e-03
    1 KSP Residual norm 2.566454735225e-04
872 KSP Residual norm 5.311640476046e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999744831e-01
    1 KSP Residual norm 2.490733310825e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.235044817215e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.548255612027e-01
    1 KSP Residual norm 8.856512457730e-03
    1 KSP Residual norm 2.122746642681e-04
873 KSP Residual norm 5.304028133887e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999828602e-01
    1 KSP Residual norm 2.549044979155e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.828871960928e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.622312578666e-01
    1 KSP Residual norm 9.072031411109e-03
    1 KSP Residual norm 1.683717190982e-04
874 KSP Residual norm 5.297958730745e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999897053e-01
    1 KSP Residual norm 2.530617076025e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.405564820886e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.907914926001e-01
    1 KSP Residual norm 8.583515677868e-03
    1 KSP Residual norm 1.227759379058e-04
875 KSP Residual norm 5.289797950943e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999910240e-01
    1 KSP Residual norm 2.589040425820e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.313546481415e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.516770543897e-01
    1 KSP Residual norm 8.446815143534e-03
    1 KSP Residual norm 1.131801128714e-04
876 KSP Residual norm 5.265823936593e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999816999e-01
    1 KSP Residual norm 2.626792233723e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.886975236830e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.923140627679e-01
    1 KSP Residual norm 8.298469925643e-03
    1 KSP Residual norm 1.758810885662e-04
877 KSP Residual norm 5.206582915262e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999613590e-01
    1 KSP Residual norm 2.645003250652e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.762824921229e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.371189304468e-01
    1 KSP Residual norm 9.648412215549e-03
    1 KSP Residual norm 2.670491947781e-04
878 KSP Residual norm 5.098934275848e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999323304e-01
    1 KSP Residual norm 2.764456965888e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.668724930676e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.648246014540e-01
    1 KSP Residual norm 1.033089944388e-02
    1 KSP Residual norm 3.596241951589e-04
879 KSP Residual norm 4.962827438525e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999142388e-01
    1 KSP Residual norm 2.744339887142e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 4.130484934034e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.653261835339e-01
    1 KSP Residual norm 1.047142410600e-02
    1 KSP Residual norm 4.066851512318e-04
880 KSP Residual norm 4.826880938808e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999143318e-01
    1 KSP Residual norm 2.727921880673e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 4.124153963081e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.438052544487e-01
    1 KSP Residual norm 9.594505730934e-03
    1 KSP Residual norm 4.063984102255e-04
881 KSP Residual norm 4.702015308248e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999240689e-01
    1 KSP Residual norm 2.503756452900e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.883465220153e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.474951146516e-01
    1 KSP Residual norm 8.828569486257e-03
    1 KSP Residual norm 3.815594530755e-04
882 KSP Residual norm 4.589486600871e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999343353e-01
    1 KSP Residual norm 2.594473594854e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.607979673769e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.225293861142e-01
    1 KSP Residual norm 1.005442966904e-02
    1 KSP Residual norm 3.538440193504e-04
883 KSP Residual norm 4.470395006935e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999316322e-01
    1 KSP Residual norm 2.787682281979e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.681898764376e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.600575563823e-01
    1 KSP Residual norm 9.176573483881e-03
    1 KSP Residual norm 3.609997684272e-04
884 KSP Residual norm 4.362037158339e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999361837e-01
    1 KSP Residual norm 2.590346560085e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.556683123661e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.439810481135e-01
    1 KSP Residual norm 8.344777001285e-03
    1 KSP Residual norm 3.486444125855e-04
885 KSP Residual norm 4.264818647925e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999425552e-01
    1 KSP Residual norm 2.624793333346e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.378604015147e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.375951095518e-01
    1 KSP Residual norm 8.743941451579e-03
    1 KSP Residual norm 3.304139236914e-04
886 KSP Residual norm 4.162077613245e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999407085e-01
    1 KSP Residual norm 2.421854578387e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 3.427145476426e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 6.105235783101e-01
    1 KSP Residual norm 9.394677258658e-03
    1 KSP Residual norm 3.347473805044e-04
887 KSP Residual norm 4.104041459006e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999706516e-01
    1 KSP Residual norm 2.552848945740e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.408542694160e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.961473288921e-01
    1 KSP Residual norm 8.475489690816e-03
    1 KSP Residual norm 2.310730104914e-04
888 KSP Residual norm 4.051594385824e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999797607e-01
    1 KSP Residual norm 2.640672731384e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.987378501134e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.282294805377e-01
    1 KSP Residual norm 8.446521503361e-03
    1 KSP Residual norm 1.867003266355e-04
889 KSP Residual norm 3.994747088446e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999824683e-01
    1 KSP Residual norm 2.578572337607e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.847695200729e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.319020713275e-01
    1 KSP Residual norm 8.655303071130e-03
    1 KSP Residual norm 1.716565401001e-04
890 KSP Residual norm 3.927017960813e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999849653e-01
    1 KSP Residual norm 2.893843148694e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.714913923829e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 4.627146296301e-01
    1 KSP Residual norm 8.613991808792e-03
    1 KSP Residual norm 1.593386879312e-04
891 KSP Residual norm 3.781726178385e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999742928e-01
    1 KSP Residual norm 2.543822657309e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.242105260143e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.670862528763e-01
    1 KSP Residual norm 9.237998555973e-03
    1 KSP Residual norm 2.131571070562e-04
892 KSP Residual norm 3.597384146255e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999734849e-01
    1 KSP Residual norm 2.654311462169e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.293579633553e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.712717598573e-01
    1 KSP Residual norm 9.357065900445e-03
    1 KSP Residual norm 2.182717164772e-04
893 KSP Residual norm 3.407287854869e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999762235e-01
    1 KSP Residual norm 2.510163216278e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.155987237242e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.690915456306e-01
    1 KSP Residual norm 9.282701730130e-03
    1 KSP Residual norm 2.042695827417e-04
894 KSP Residual norm 3.243283927127e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999829098e-01
    1 KSP Residual norm 2.698253063950e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.820080242004e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.223887505680e-01
    1 KSP Residual norm 8.483960541761e-03
    1 KSP Residual norm 1.688953333530e-04
895 KSP Residual norm 3.066143682941e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999851791e-01
    1 KSP Residual norm 2.470643280024e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.696519487503e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.328976774892e-01
    1 KSP Residual norm 8.719873280040e-03
    1 KSP Residual norm 1.550364098784e-04
896 KSP Residual norm 2.858172087964e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999855601e-01
    1 KSP Residual norm 2.568660438411e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.674644418657e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.161607679354e-01
    1 KSP Residual norm 8.669875299815e-03
    1 KSP Residual norm 1.528608785186e-04
897 KSP Residual norm 2.588954412031e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999825789e-01
    1 KSP Residual norm 2.753012309020e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.843315447895e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.212180337037e-01
    1 KSP Residual norm 8.684708427077e-03
    1 KSP Residual norm 1.710031196471e-04
898 KSP Residual norm 2.256718151201e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999760160e-01
    1 KSP Residual norm 2.682264127232e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.170701569099e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.702072364353e-01
    1 KSP Residual norm 9.698055434273e-03
    1 KSP Residual norm 2.049959694642e-04
899 KSP Residual norm 1.976711538311e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999999748076e-01
    1 KSP Residual norm 2.831721853798e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 2.222663050362e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.750514624447e-01
    1 KSP Residual norm 9.404384631283e-03
    1 KSP Residual norm 2.105552239259e-04
900 KSP Residual norm 1.786462293904e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999983468747e-01
    1 KSP Residual norm 1.364269546626e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 1.817756705757e-03
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.493098953106e-01
    1 KSP Residual norm 1.198945933762e-02
    1 KSP Residual norm 1.815174387380e-03
901 KSP Residual norm 1.690320893601e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999998115411e-01
    1 KSP Residual norm 3.269248862736e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.138700673278e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 3.758899308390e-01
    1 KSP Residual norm 7.951667253883e-03
    1 KSP Residual norm 6.108132284271e-04
902 KSP Residual norm 1.574339227111e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999997338107e-01
    1 KSP Residual norm 2.708769142179e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.288071587055e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.135596800446e-01
    1 KSP Residual norm 9.135297858300e-03
    1 KSP Residual norm 7.256674055521e-04
903 KSP Residual norm 1.430140155907e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999995817412e-01
    1 KSP Residual norm 2.780472261188e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 9.143042178830e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.763929125962e-01
    1 KSP Residual norm 1.002153942926e-02
    1 KSP Residual norm 9.114575342396e-04
904 KSP Residual norm 1.300000530638e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999995028410e-01
    1 KSP Residual norm 2.785877931845e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 9.968269364968e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.811489563167e-01
    1 KSP Residual norm 1.062384504647e-02
    1 KSP Residual norm 9.942193439674e-04
905 KSP Residual norm 1.211080970732e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.999995693878e-01
    1 KSP Residual norm 2.983280892357e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 9.275726050786e-04
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.275115786546e-01
    1 KSP Residual norm 9.761851016203e-03
    1 KSP Residual norm 9.249820914173e-04
906 KSP Residual norm 1.139759024126e-09



%%%% SBAIJ again with -pc_fieldsplit_schur_fact_type FULL. This helps quite a bit

503 KSP Residual norm 1.335198148782e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.972328162111e-01
    1 KSP Residual norm 4.518715383704e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 7.434142477591e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.659720275605e+01
    1 KSP Residual norm 6.189926845964e-01
    1 KSP Residual norm 7.403062361616e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.673081251569e+01
    1 KSP Residual norm 6.239752505178e-01
504 KSP Residual norm 1.279949951034e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.976716150443e-01
    1 KSP Residual norm 4.866857915634e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 6.820038403874e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.135988214250e+01
    1 KSP Residual norm 5.144348811046e-01
    1 KSP Residual norm 6.781410979594e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 5.144970605895e+01
    1 KSP Residual norm 5.148763638488e-01
505 KSP Residual norm 1.186275242447e-09
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 9.953591307188e-01
    1 KSP Residual norm 4.518765750227e-02
    Residual norms for fieldsplit_1_ solve.
    0 KSP Residual norm 9.622980462164e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 6.474321429841e+01
    1 KSP Residual norm 6.421404156501e-01
    1 KSP Residual norm 9.592264156459e-02
    Residual norms for fieldsplit_0_ solve.
    0 KSP Residual norm 6.482110598687e+01
    1 KSP Residual norm 6.359439449320e-01
506 KSP Residual norm 1.128137176971e-09
Linear solve converged due to CONVERGED_RTOL iterations 506
KSP Object: 4 MPI processes
  type: fgmres
    restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
    happy breakdown tolerance 1e-30
  maximum iterations=1000, nonzero initial guess
  minimum iterations=3
  tolerances:  relative=1e-10, absolute=1e-50, divergence=1e+09
  right preconditioning
  using UNPRECONDITIONED norm type for convergence test
PC Object: 4 MPI processes
  type: fieldsplit
    FieldSplit with Schur preconditioner, factorization FULL
    Preconditioner for the Schur complement formed from Sp, an assembled approximation to S, which uses A00's diagonal's inverse
    Split info:
    Split number 0 Defined by IS
    Split number 1 Defined by IS
    KSP solver for A00 block
      KSP Object: (fieldsplit_0_) 4 MPI processes
        type: preonly
        maximum iterations=10000, initial guess is zero
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
        left preconditioning
        using NONE norm type for convergence test
      PC Object: (fieldsplit_0_) 4 MPI processes
        type: asm
          total subdomain blocks = 4, amount of overlap = 1
          restriction/interpolation type - RESTRICT
          Local solver information for first block is in the following KSP and PC objects on rank 0:
          Use -fieldsplit_0_ksp_view ::ascii_info_detail to display information for all blocks
        KSP Object: (fieldsplit_0_sub_) 1 MPI process
          type: preonly
          maximum iterations=10000, initial guess is zero
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
          left preconditioning
          using NONE norm type for convergence test
        PC Object: (fieldsplit_0_sub_) 1 MPI process
          type: icc
            out-of-place factorization
            0 levels of fill
            tolerance for zero pivot 2.22045e-14
            using Manteuffel shift [POSITIVE_DEFINITE]
            matrix ordering: natural
            factor fill ratio given 1., needed 1.
              Factored matrix follows:
                Mat Object: (fieldsplit_0_sub_) 1 MPI process
                  type: seqsbaij
                  rows=97788, cols=97788
                  package used to perform factorization: petsc
                  total: nonzeros=1153455, allocated nonzeros=1153455
                      block size is 1
          linear system matrix = precond matrix:
          Mat Object: (fieldsplit_0_sub_) 1 MPI process
            type: seqsbaij
            rows=97788, cols=97788
            total: nonzeros=1153455, allocated nonzeros=1153455
            total number of mallocs used during MatSetValues calls=0
                block size is 1
        linear system matrix = precond matrix:
        Mat Object: (fieldsplit_0_) 4 MPI processes
          type: mpisbaij
          rows=392931, cols=392931
          total: nonzeros=5739651, allocated nonzeros=5739651
          total number of mallocs used during MatSetValues calls=0
              block size is 1
    KSP solver for S = A11 - A10 inv(A00) A01
      KSP Object: (fieldsplit_1_) 4 MPI processes
        type: preonly
        maximum iterations=10000, initial guess is zero
        tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
        left preconditioning
        using NONE norm type for convergence test
      PC Object: (fieldsplit_1_) 4 MPI processes
        type: asm
          total subdomain blocks = 4, amount of overlap = 1
          restriction/interpolation type - RESTRICT
          Local solver information for first block is in the following KSP and PC objects on rank 0:
          Use -fieldsplit_1_ksp_view ::ascii_info_detail to display information for all blocks
        KSP Object: (fieldsplit_1_sub_) 1 MPI process
          type: preonly
          maximum iterations=10000, initial guess is zero
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000.
          left preconditioning
          using NONE norm type for convergence test
        PC Object: (fieldsplit_1_sub_) 1 MPI process
          type: ilu
            out-of-place factorization
            0 levels of fill
            tolerance for zero pivot 2.22045e-14
            matrix ordering: natural
            factor fill ratio given 1., needed 1.
              Factored matrix follows:
                Mat Object: (fieldsplit_1_sub_) 1 MPI process
                  type: seqaij
                  rows=31038, cols=31038
                  package used to perform factorization: petsc
                  total: nonzeros=223972, allocated nonzeros=223972
                    not using I-node routines
          linear system matrix = precond matrix:
          Mat Object: (fieldsplit_1_sub_) 1 MPI process
            type: seqaij
            rows=31038, cols=31038
            total: nonzeros=223972, allocated nonzeros=223972
            total number of mallocs used during MatSetValues calls=0
              not using I-node routines
        linear system matrix followed by preconditioner matrix:
        Mat Object: (fieldsplit_1_) 4 MPI processes
          type: schurcomplement
          rows=121600, cols=121600
            Schur complement A11 - A10 inv(A00) A01
            A11
              Mat Object: (fieldsplit_1_) 4 MPI processes
                type: mpisbaij
                rows=121600, cols=121600
                total: nonzeros=242540, allocated nonzeros=242540
                total number of mallocs used during MatSetValues calls=0
                    block size is 1
            A10
              Mat Object: 4 MPI processes
                type: hermitiantranspose
                rows=121600, cols=392931
            KSP solver for A00 block viewable with the additional option -fieldsplit_0_ksp_view
            A01
              Mat Object: 4 MPI processes
                type: mpiaij
                rows=392931, cols=121600
                total: nonzeros=1018872, allocated nonzeros=1018872
                total number of mallocs used during MatSetValues calls=0
                  using I-node (on process 0) routines: found 23081 nodes, limit used is 5
        Mat Object: 4 MPI processes
          type: mpiaij
          rows=121600, cols=121600
          total: nonzeros=1120336, allocated nonzeros=1120336
          total number of mallocs used during MatSetValues calls=0
            not using I-node (on process 0) routines
  linear system matrix = precond matrix:
  Mat Object: 4 MPI processes
    type: mpisbaij
    rows=514531, cols=514531
    total: nonzeros=7001063, allocated nonzeros=26569279
    total number of mallocs used during MatSetValues calls=0
        block size is 1


More information about the petsc-users mailing list