[petsc-users] Error message for DMShell using MG preconditioner to solve S

Luc Berger-Vergiat lb2653 at columbia.edu
Tue May 20 16:33:18 CDT 2014


Hi all,
I am running an FEM simulation that uses Petsc as a linear solver.
I am setting up ISs and pass them to a DMShell in order to use the 
FieldSplit capabilities of Petsc.

When I pass the following options to Petsc:

    " -ksp_type gmres -pc_type fieldsplit -pc_fieldsplit_type schur
    -pc_fieldsplit_schur_factorization_type full
    -pc_fieldsplit_schur_precondition selfp -pc_fieldsplit_0_fields 1,2
    -pc_fieldsplit_1_fields 0 -fieldsplit_0_ksp_type preonly
    -fieldsplit_0_pc_type ilu -fieldsplit_Field_0_ksp_type gmres
    -fieldsplit_Field_0_pc_type mg -malloc_log mlog -log_summary time.log"

I get an error message:

[0]PETSC ERROR: --------------------- Error Message 
--------------------------------------------------------------
[0]PETSC ERROR:
[0]PETSC ERROR: Must call DMShellSetGlobalVector() or 
DMShellSetCreateGlobalVector()
[0]PETSC ERROR: See 
http://http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble 
shooting.
[0]PETSC ERROR: Petsc Development GIT revision: v3.4.4-5071-g1163a46  
GIT Date: 2014-03-26 22:20:51 -0500
[0]PETSC ERROR: /home/luc/research/feap_repo/ShearBands/parfeap-dev/feap 
on a arch-linux2-c-opt named euler by luc Tue May 20 11:31:11 2014
[0]PETSC ERROR: Configure options --download-cmake --download-hypre 
--download-metis --download-mpich --download-parmetis 
--with-debugging=no --with-share-libraries=no
[0]PETSC ERROR: #1 DMCreateGlobalVector_Shell() line 245 in 
/home/luc/research/petsc/src/dm/impls/shell/dmshell.c
[0]PETSC ERROR: #2 DMCreateGlobalVector() line 669 in 
/home/luc/research/petsc/src/dm/interface/dm.c
[0]PETSC ERROR: #3 DMGetGlobalVector() line 154 in 
/home/luc/research/petsc/src/dm/interface/dmget.c

I am not really sure why this happens but it only happens when 
-fieldsplit_Field_0_pc_type mg, with other preconditioners, I have no 
problems. I attached the ksp_view in case that's any use.

-- 
Best,
Luc

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140520/0c1e7f70/attachment-0001.html>
-------------- next part --------------
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-08, absolute=1e-16, divergence=1e+16
  left preconditioning
  using PRECONDITIONED norm type for convergence test
PC Object: 1 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_)       1 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_)       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=2000, cols=2000
                package used to perform factorization: petsc
                total: nonzeros=40000, allocated nonzeros=40000
                total number of mallocs used during MatSetValues calls =0
                  using I-node routines: found 400 nodes, limit used is 5
        linear system matrix = precond matrix:
        Mat Object:        (fieldsplit_0_)         1 MPI processes
          type: seqaij
          rows=2000, cols=2000
          total: nonzeros=40000, allocated nonzeros=40000
          total number of mallocs used during MatSetValues calls =0
            using I-node routines: found 400 nodes, limit used is 5
    KSP solver for S = A11 - A10 inv(A00) A01 
      KSP Object:      (fieldsplit_Field_0_)       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
        using PRECONDITIONED norm type for convergence test
      PC Object:      (fieldsplit_Field_0_)       1 MPI processes
        type: mg
          MG: type is MULTIPLICATIVE, levels=1 cycles=v
            Cycles per PCApply=1
            Not using Galerkin computed coarse grid matrices
        Coarse grid solver -- level -------------------------------
          KSP Object:          (fieldsplit_Field_0_mg_levels_0_)           1 MPI processes
            type: chebyshev
              Chebyshev: eigenvalue estimates:  min = 0.937483, max = 10.3123
              Chebyshev: estimated using:  [0 0.1; 0 1.1]
              KSP Object:              (fieldsplit_Field_0_mg_levels_0_est_)               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=10, 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_Field_0_mg_levels_0_)               1 MPI processes
                type: sor
                  SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
                linear system matrix followed by preconditioner matrix:
                Mat Object:                (fieldsplit_Field_0_)                 1 MPI processes
                  type: schurcomplement
                  rows=209, cols=209
                    Schur complement A11 - A10 inv(A00) A01
                    A11
                      Mat Object:                      (fieldsplit_Field_0_)                       1 MPI processes
                        type: seqaij
                        rows=209, cols=209
                        total: nonzeros=3209, allocated nonzeros=3209
                        total number of mallocs used during MatSetValues calls =0
                          using I-node routines: found 119 nodes, limit used is 5
                    A10
                      Mat Object:                       1 MPI processes
                        type: seqaij
                        rows=209, cols=2000
                        total: nonzeros=14800, allocated nonzeros=14800
                        total number of mallocs used during MatSetValues calls =0
                          using I-node routines: found 119 nodes, limit used is 5
                    KSP of A00
                      KSP Object:                      (fieldsplit_0_)                       1 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_)                       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=2000, cols=2000
                                package used to perform factorization: petsc
                                total: nonzeros=40000, allocated nonzeros=40000
                                total number of mallocs used during MatSetValues calls =0
                                  using I-node routines: found 400 nodes, limit used is 5
                        linear system matrix = precond matrix:
                        Mat Object:                        (fieldsplit_0_)                         1 MPI processes
                          type: seqaij
                          rows=2000, cols=2000
                          total: nonzeros=40000, allocated nonzeros=40000
                          total number of mallocs used during MatSetValues calls =0
                            using I-node routines: found 400 nodes, limit used is 5
                    A01
                      Mat Object:                       1 MPI processes
                        type: seqaij
                        rows=2000, cols=209
                        total: nonzeros=14800, allocated nonzeros=14800
                        total number of mallocs used during MatSetValues calls =0
                          using I-node routines: found 400 nodes, limit used is 5
                Mat Object:                 1 MPI processes
                  type: seqaij
                  rows=209, cols=209
                  total: nonzeros=3209, allocated nonzeros=3209
                  total number of mallocs used during MatSetValues calls =0
                    using I-node routines: found 119 nodes, limit used is 5
            maximum iterations=1, 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_Field_0_mg_levels_0_)           1 MPI processes
            type: sor
              SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
            linear system matrix followed by preconditioner matrix:
            Mat Object:            (fieldsplit_Field_0_)             1 MPI processes
              type: schurcomplement
              rows=209, cols=209
                Schur complement A11 - A10 inv(A00) A01
                A11
                  Mat Object:                  (fieldsplit_Field_0_)                   1 MPI processes
                    type: seqaij
                    rows=209, cols=209
                    total: nonzeros=3209, allocated nonzeros=3209
                    total number of mallocs used during MatSetValues calls =0
                      using I-node routines: found 119 nodes, limit used is 5
                A10
                  Mat Object:                   1 MPI processes
                    type: seqaij
                    rows=209, cols=2000
                    total: nonzeros=14800, allocated nonzeros=14800
                    total number of mallocs used during MatSetValues calls =0
                      using I-node routines: found 119 nodes, limit used is 5
                KSP of A00
                  KSP Object:                  (fieldsplit_0_)                   1 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_)                   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=2000, cols=2000
                            package used to perform factorization: petsc
                            total: nonzeros=40000, allocated nonzeros=40000
                            total number of mallocs used during MatSetValues calls =0
                              using I-node routines: found 400 nodes, limit used is 5
                    linear system matrix = precond matrix:
                    Mat Object:                    (fieldsplit_0_)                     1 MPI processes
                      type: seqaij
                      rows=2000, cols=2000
                      total: nonzeros=40000, allocated nonzeros=40000
                      total number of mallocs used during MatSetValues calls =0
                        using I-node routines: found 400 nodes, limit used is 5
                A01
                  Mat Object:                   1 MPI processes
                    type: seqaij
                    rows=2000, cols=209
                    total: nonzeros=14800, allocated nonzeros=14800
                    total number of mallocs used during MatSetValues calls =0
                      using I-node routines: found 400 nodes, limit used is 5
            Mat Object:             1 MPI processes
              type: seqaij
              rows=209, cols=209
              total: nonzeros=3209, allocated nonzeros=3209
              total number of mallocs used during MatSetValues calls =0
                using I-node routines: found 119 nodes, limit used is 5
        linear system matrix followed by preconditioner matrix:
        Mat Object:        (fieldsplit_Field_0_)         1 MPI processes
          type: schurcomplement
          rows=209, cols=209
            Schur complement A11 - A10 inv(A00) A01
            A11
              Mat Object:              (fieldsplit_Field_0_)               1 MPI processes
                type: seqaij
                rows=209, cols=209
                total: nonzeros=3209, allocated nonzeros=3209
                total number of mallocs used during MatSetValues calls =0
                  using I-node routines: found 119 nodes, limit used is 5
            A10
              Mat Object:               1 MPI processes
                type: seqaij
                rows=209, cols=2000
                total: nonzeros=14800, allocated nonzeros=14800
                total number of mallocs used during MatSetValues calls =0
                  using I-node routines: found 119 nodes, limit used is 5
            KSP of A00
              KSP Object:              (fieldsplit_0_)               1 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_)               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=2000, cols=2000
                        package used to perform factorization: petsc
                        total: nonzeros=40000, allocated nonzeros=40000
                        total number of mallocs used during MatSetValues calls =0
                          using I-node routines: found 400 nodes, limit used is 5
                linear system matrix = precond matrix:
                Mat Object:                (fieldsplit_0_)                 1 MPI processes
                  type: seqaij
                  rows=2000, cols=2000
                  total: nonzeros=40000, allocated nonzeros=40000
                  total number of mallocs used during MatSetValues calls =0
                    using I-node routines: found 400 nodes, limit used is 5
            A01
              Mat Object:               1 MPI processes
                type: seqaij
                rows=2000, cols=209
                total: nonzeros=14800, allocated nonzeros=14800
                total number of mallocs used during MatSetValues calls =0
                  using I-node routines: found 400 nodes, limit used is 5
        Mat Object:         1 MPI processes
          type: seqaij
          rows=209, cols=209
          total: nonzeros=3209, allocated nonzeros=3209
          total number of mallocs used during MatSetValues calls =0
            using I-node routines: found 119 nodes, limit used is 5
  linear system matrix = precond matrix:
  Mat Object:   1 MPI processes
    type: seqaij
    rows=2209, cols=2209
    total: nonzeros=72809, allocated nonzeros=72809
    total number of mallocs used during MatSetValues calls =0
      using I-node routines: found 519 nodes, limit used is 5


More information about the petsc-users mailing list