[petsc-users] Multigrid iterations for different versions of Petsc

Filippo Leonardi filippo.leonardi at sam.math.ethz.ch
Mon Oct 6 08:15:25 CDT 2014


On Monday 06 October 2014 07:36:47 Barry Smith wrote:
>   First make sure both runs think they are doing the exact same solver
> options so run both with -ksp_view and make sure the output is the same.

Ok, now I am really confused. I did a "diff" onto the ouput of kspview 
(attached) and the result is... weird. The size of the matrices is cubed for 
the version 3.4.3, although I am sure the mesh size and the initial matrix 
(BTW a DMDA get Jacobian) are the same in the two cases (I even output the 
result) (the code and the command line options are exactly the same except for 
API changes).

Is this a bug?

> 
>   I note that the 3.5 run has
> 
> KSPSolve            2200 1.0 8.3100e+00 1.0 7.28e+08 1.0 1.4e+05 2.5e+03
> 5.9e+03 55 67 45 44 79  83 83 47 80100   700 PCApply             5764 1.0
> 9.4324e-01 1.1 3.12e+08 1.0 2.5e+04 6.4e+01 0.0e+00  6 29  8  0  0   9 36 
> 8  0  0  2646
> 
> while the 3.4 has
> 
> KSPSolve            4125 1.0 7.7953e+00 1.0 1.39e+09 1.0 2.6e+05 2.5e+03
> 1.1e+04 60 65 45 53 82  82 83 47 80100  1424 PCApply            10769 1.0
> 4.4700e+00 1.0 5.80e+08 1.0 4.6e+04 6.4e+01 0.0e+00 34 27  8  0  0  47 35 
> 8  0  0  1038
> 
> why twice as many KSPSolves? Are you calling the KSPSolve more ?

I call KSPSolve 11 times here (I am looking at the stage 5), and then MG is 
calling KSPSolve by its own and the doubling of iterations is explained by the 
doubling of iterations.

> 
>   Barry
> 
> On Oct 6, 2014, at 7:23 AM, Filippo Leonardi 
<filippo.leonardi at sam.math.ethz.ch> wrote:
> > Hi,
> > 
> > Again me with multigrid. Sorry for the annoyance, but I am a bit worried
> > about my multigrid solver: I still can't figure out where (if something)
> > my code is doing wrong.
> > 
> > I post the details of a run with my local machine and Petsc 3.5 compared
> > with the runs done with Petsc 3.4.3 on a cluster (the most recent version
> > available).
> > 
> > The problem is that the run in the cluster (3.4.3) does twice many
> > iteration per krylow solver (as the local 3.5) and goes even worse when
> > increasing the mesh. Also, while on my local machine I get 7 iterations
> > (which is quite nice), on the cluster I get 14 iteration, which is not
> > what I expect from multigrid.
> > 
> > Are there been any changes on the MG solver that can influence the
> > behaviour of the solver?
> > 
> > http://www.mcs.anl.gov/petsc/documentation/changes/35.html
> > seems not to include anything relevant.
> > 
> > Everything is the same: outer iteration with gmres, lu as coarse solver
> > etc
> > (all default), I just adapted the APIs. Thus I'd expect more or less the
> > same result.
> > 
> > Best,
> > Filippo
> > 
> > I attach the profiles.<3-4.txt><3-5.txt><ETHZ.vcf>
-------------- next part --------------
KSP Object: 8 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: 8 MPI processes
  type: mg
    MG: type is FULL, levels=5 cycles=v
      Not using Galerkin computed coarse grid matrices
  Coarse grid solver -- level -------------------------------
    KSP Object:    (mg_coarse_)     8 MPI processes
      type: preonly
      maximum iterations=1, initial guess is zero
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using NONE norm type for convergence test
    PC Object:    (mg_coarse_)     8 MPI processes
      type: redundant
        Redundant preconditioner: First (color=0) of 8 PCs follows
      KSP Object:      (mg_coarse_redundant_)       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:      (mg_coarse_redundant_)       1 MPI processes
        type: lu
          LU: out-of-place factorization
          tolerance for zero pivot 2.22045e-14
          using diagonal shift on blocks to prevent zero pivot [INBLOCKS]
          matrix ordering: nd
          factor fill ratio given 5, needed 4.42857
            Factored matrix follows:
              Matrix Object:               1 MPI processes
                type: seqaij
                rows=64, cols=64
                package used to perform factorization: petsc
                total: nonzeros=1984, allocated nonzeros=1984
                total number of mallocs used during MatSetValues calls =0
                  not using I-node routines
        linear system matrix = precond matrix:
        Matrix Object:         1 MPI processes
          type: seqaij
          rows=64, cols=64
          total: nonzeros=448, allocated nonzeros=448
          total number of mallocs used during MatSetValues calls =0
            not using I-node routines
      linear system matrix = precond matrix:
      Matrix Object:       8 MPI processes
        type: mpiaij
        rows=64, cols=64
        total: nonzeros=448, allocated nonzeros=448
        total number of mallocs used during MatSetValues calls =0
  Down solver (pre-smoother) on level 1 -------------------------------
    KSP Object:    (mg_levels_1_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.153663, max = 1.6903
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_1_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_1_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Matrix Object:           8 MPI processes
            type: mpiaij
            rows=512, cols=512
            total: nonzeros=3584, allocated nonzeros=3584
            total number of mallocs used during MatSetValues calls =0
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_1_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Matrix Object:       8 MPI processes
        type: mpiaij
        rows=512, cols=512
        total: nonzeros=3584, allocated nonzeros=3584
        total number of mallocs used during MatSetValues calls =0
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 2 -------------------------------
    KSP Object:    (mg_levels_2_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.153015, max = 1.68317
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_2_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_2_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Matrix Object:           8 MPI processes
            type: mpiaij
            rows=4096, cols=4096
            total: nonzeros=28672, allocated nonzeros=28672
            total number of mallocs used during MatSetValues calls =0
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_2_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Matrix Object:       8 MPI processes
        type: mpiaij
        rows=4096, cols=4096
        total: nonzeros=28672, allocated nonzeros=28672
        total number of mallocs used during MatSetValues calls =0
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 3 -------------------------------
    KSP Object:    (mg_levels_3_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.152678, max = 1.67946
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_3_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_3_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Matrix Object:           8 MPI processes
            type: mpiaij
            rows=32768, cols=32768
            total: nonzeros=229376, allocated nonzeros=229376
            total number of mallocs used during MatSetValues calls =0
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_3_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Matrix Object:       8 MPI processes
        type: mpiaij
        rows=32768, cols=32768
        total: nonzeros=229376, allocated nonzeros=229376
        total number of mallocs used during MatSetValues calls =0
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 4 -------------------------------
    KSP Object:    (mg_levels_4_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.147343, max = 1.62077
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_4_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_4_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Matrix Object:           8 MPI processes
            type: mpiaij
            rows=262144, cols=262144
            total: nonzeros=1835008, allocated nonzeros=1835008
            total number of mallocs used during MatSetValues calls =0
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_4_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Matrix Object:       8 MPI processes
        type: mpiaij
        rows=262144, cols=262144
        total: nonzeros=1835008, allocated nonzeros=1835008
        total number of mallocs used during MatSetValues calls =0
  Up solver (post-smoother) same as down solver (pre-smoother)
  linear system matrix = precond matrix:
  Matrix Object:   8 MPI processes
    type: mpiaij
    rows=262144, cols=262144
    total: nonzeros=1835008, allocated nonzeros=1835008
    total number of mallocs used during MatSetValues calls =0
-------------- next part --------------
KSP Object: 8 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: 8 MPI processes
  type: mg
    MG: type is FULL, levels=5 cycles=v
      Not using Galerkin computed coarse grid matrices
  Coarse grid solver -- level -------------------------------
    KSP Object:    (mg_coarse_)     8 MPI processes
      type: preonly
      maximum iterations=1, initial guess is zero
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using NONE norm type for convergence test
    PC Object:    (mg_coarse_)     8 MPI processes
      type: redundant
        Redundant preconditioner: First (color=0) of 8 PCs follows
      KSP Object:      (mg_coarse_redundant_)       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:      (mg_coarse_redundant_)       1 MPI processes
        type: lu
          LU: out-of-place factorization
          tolerance for zero pivot 2.22045e-14
          using diagonal shift on blocks to prevent zero pivot [INBLOCKS]
          matrix ordering: nd
          factor fill ratio given 5, needed 1.375
            Factored matrix follows:
              Mat Object:               1 MPI processes
                type: seqaij
                rows=8, cols=8
                package used to perform factorization: petsc
                total: nonzeros=44, allocated nonzeros=44
                total number of mallocs used during MatSetValues calls =0
                  not using I-node routines
        linear system matrix = precond matrix:
        Mat Object:         1 MPI processes
          type: seqaij
          rows=8, cols=8
          total: nonzeros=32, allocated nonzeros=32
          total number of mallocs used during MatSetValues calls =0
            not using I-node routines
      linear system matrix = precond matrix:
      Mat Object:       8 MPI processes
        type: mpiaij
        rows=8, cols=8
        total: nonzeros=32, allocated nonzeros=56
        total number of mallocs used during MatSetValues calls =0
          has attached null space
  Down solver (pre-smoother) on level 1 -------------------------------
    KSP Object:    (mg_levels_1_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.155469, max = 1.71016
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_1_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_1_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Mat Object:           8 MPI processes
            type: mpiaij
            rows=64, cols=64
            total: nonzeros=448, allocated nonzeros=448
            total number of mallocs used during MatSetValues calls =0
              has attached null space
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_1_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Mat Object:       8 MPI processes
        type: mpiaij
        rows=64, cols=64
        total: nonzeros=448, allocated nonzeros=448
        total number of mallocs used during MatSetValues calls =0
          has attached null space
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 2 -------------------------------
    KSP Object:    (mg_levels_2_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.15366, max = 1.69026
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_2_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_2_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Mat Object:           8 MPI processes
            type: mpiaij
            rows=512, cols=512
            total: nonzeros=3584, allocated nonzeros=3584
            total number of mallocs used during MatSetValues calls =0
              has attached null space
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_2_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Mat Object:       8 MPI processes
        type: mpiaij
        rows=512, cols=512
        total: nonzeros=3584, allocated nonzeros=3584
        total number of mallocs used during MatSetValues calls =0
          has attached null space
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 3 -------------------------------
    KSP Object:    (mg_levels_3_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.152956, max = 1.68252
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_3_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_3_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Mat Object:           8 MPI processes
            type: mpiaij
            rows=4096, cols=4096
            total: nonzeros=28672, allocated nonzeros=28672
            total number of mallocs used during MatSetValues calls =0
              has attached null space
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_3_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Mat Object:       8 MPI processes
        type: mpiaij
        rows=4096, cols=4096
        total: nonzeros=28672, allocated nonzeros=28672
        total number of mallocs used during MatSetValues calls =0
          has attached null space
  Up solver (post-smoother) same as down solver (pre-smoother)
  Down solver (pre-smoother) on level 4 -------------------------------
    KSP Object:    (mg_levels_4_)     8 MPI processes
      type: chebyshev
        Chebyshev: eigenvalue estimates:  min = 0.149148, max = 1.64063
        Chebyshev: estimated using:  [0 0.1; 0 1.1]
        KSP Object:        (mg_levels_4_est_)         8 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
          tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
          left preconditioning
          using nonzero initial guess
          using NONE norm type for convergence test
        PC Object:        (mg_levels_4_)         8 MPI processes
          type: sor
            SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
          linear system matrix = precond matrix:
          Mat Object:           8 MPI processes
            type: mpiaij
            rows=32768, cols=32768
            total: nonzeros=229376, allocated nonzeros=229376
            total number of mallocs used during MatSetValues calls =0
              has attached null space
      maximum iterations=2
      tolerances:  relative=1e-05, absolute=1e-50, divergence=10000
      left preconditioning
      has attached null space
      using nonzero initial guess
      using NONE norm type for convergence test
    PC Object:    (mg_levels_4_)     8 MPI processes
      type: sor
        SOR: type = local_symmetric, iterations = 1, local iterations = 1, omega = 1
      linear system matrix = precond matrix:
      Mat Object:       8 MPI processes
        type: mpiaij
        rows=32768, cols=32768
        total: nonzeros=229376, allocated nonzeros=229376
        total number of mallocs used during MatSetValues calls =0
          has attached null space
  Up solver (post-smoother) same as down solver (pre-smoother)
  linear system matrix = precond matrix:
  Mat Object:   8 MPI processes
    type: mpiaij
    rows=32768, cols=32768
    total: nonzeros=229376, allocated nonzeros=229376
    total number of mallocs used during MatSetValues calls =0
      has attached null space
-------------- next part --------------
A non-text attachment was scrubbed...
Name: ETHZ.vcf
Type: text/vcard
Size: 594 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20141006/453f60d9/attachment.bin>


More information about the petsc-users mailing list