[petsc-users] How to do the point-block ILU in PETSc

Jed Brown jed at jedbrown.org
Sun May 4 08:56:44 CDT 2014


Please use the mailing list for questions like this.

Lulu Liu <lulu.liu at kaust.edu.sa> writes:

> Dear Jed,
>
> I saw in man-page of PCILU:
> For BAIJ matrices this implements a point block ILU
>
> Take /src/snes/examples/tutorials/ex19.c for examples, I add the following
> lines
>   ierr = MatCreate(PETSC_COMM_WORLD,&J);
>   ierr = MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mx*my,mx*my);

This example uses DMDA so creating your own layout won't generally be
the partition you want.  You should use DMCreateMatrix().  The matrix
type is set via DMSetMatType() and -dm_mat_type.

>   ierr = MatSetType(J,MATBAIJ);
>   ierr = MatSetFromOptions(J);
>   ierr = MatSetUp(J);
>   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
>   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);

This assembly should not exist.

>   ierr = SNESSetJacobian(snes,J,J,NULL,NULL);
>
> but I got errors, could you tell me how to do the point-block ILU in ex19.c
> ( the small block should be 4x4). Thanks!
>
> ./ex19 -da_grid_x 64 -da_grid_y 64 -contours -draw_pause 1  -snes_monitor
> -snes_rtol 1.e-6 -pc_type ilu

Don't modify the source at all.  Instead, run this:

$ mpiexec -n 4 ./ex19 -da_grid_x 64 -da_grid_y 64 -snes_monitor -snes_view -dm_mat_type baij                                                                                        
lid velocity = 0.000244141, prandtl # = 1, grashof # = 1
  0 SNES Function norm 1.573890417811e-02 
  1 SNES Function norm 1.602010905072e-06 
  2 SNES Function norm 1.580493963868e-11 
SNES Object: 4 MPI processes
  type: newtonls
  maximum iterations=50, maximum function evaluations=10000
  tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
  total number of linear solver iterations=368
  total number of function evaluations=3
  SNESLineSearch Object:   4 MPI processes
    type: bt
      interpolation: cubic
      alpha=1.000000e-04
    maxstep=1.000000e+08, minlambda=1.000000e-12
    tolerances: relative=1.000000e-08, absolute=1.000000e-15, lambda=1.000000e-08
    maximum iterations=40
  KSP Object:   4 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:   4 MPI processes
    type: bjacobi
      block Jacobi: number of blocks = 4
      Local solve is same for all blocks, in the following KSP and PC objects:
    KSP Object:    (sub_)     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:    (sub_)     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: seqbaij
              rows=4096, cols=4096, bs=4
              package used to perform factorization: petsc
              total: nonzeros=79872, allocated nonzeros=79872
              total number of mallocs used during MatSetValues calls =0
                  block size is 4
      linear system matrix = precond matrix:
      Mat Object:       1 MPI processes
        type: seqbaij
        rows=4096, cols=4096, bs=4
        total: nonzeros=79872, allocated nonzeros=79872
        total number of mallocs used during MatSetValues calls =0
            block size is 4
    linear system matrix = precond matrix:
    Mat Object:     4 MPI processes
      type: mpibaij
      rows=16384, cols=16384, bs=4
      total: nonzeros=323584, allocated nonzeros=323584
      total number of mallocs used during MatSetValues calls =0
Number of SNES iterations = 2


> lid velocity = 0.000244141, prandtl # = 1, grashof # = 1
>   0 SNES Function norm 1.573890417811e-02
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Caught signal number 11 SEGV: Segmentation Violation,
> probably memory access out of range
> [0]PETSC ERROR: Try option -start_in_debugger or -on_error_attach_debugger
> [0]PETSC ERROR: or see
> http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind[0]PETSC ERROR:
> or try http://valgrind.org on GNU/linux and Apple Mac OS X to find memory
> corruption errors
> [0]PETSC ERROR: likely location of problem given in stack below
> [0]PETSC ERROR: ---------------------  Stack Frames
> ------------------------------------
> [0]PETSC ERROR: Note: The EXACT line numbers in the stack are not available,
> [0]PETSC ERROR:       INSTEAD the line number of the start of the function
> [0]PETSC ERROR:       is given.
> [0]PETSC ERROR: [0] MatFDColoringCreate_SeqAIJ line 20
> src/mat/impls/aij/seq/fdaij.c
> [0]PETSC ERROR: [0] MatFDColoringCreate line 367 src/mat/matfd/fdmatrix.c
> [0]PETSC ERROR: [0] SNESComputeJacobian_DMDA line 165
> src/snes/utils/dmdasnes.c
> [0]PETSC ERROR: [0] SNES user Jacobian function line 2151
> src/snes/interface/snes.c
> [0]PETSC ERROR: [0] SNESComputeJacobian line 2106 src/snes/interface/snes.c
> [0]PETSC ERROR: [0] SNESSolve_NEWTONLS line 144 src/snes/impls/ls/ls.c
> [0]PETSC ERROR: [0] SNESSolve line 3589 src/snes/interface/snes.c
> [0]PETSC ERROR: [0] main line 106 src/snes/examples/tutorials/ex19.c
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Signal received!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.4.3, unknown
> [0]PETSC ERROR: See docs/changes/index.html for recent updates.
> [0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> [0]PETSC ERROR: See docs/index.html for manual pages.
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: ./ex19 on a arch-darwin-c-debug named kl-12681.local by
> liul Sun May  4 15:43:58 2014
> [0]PETSC ERROR: Libraries linked from
> /Users/liul/soft/petsc-3.4.3/petsc/arch-darwin-c-debug/lib
> [0]PETSC ERROR: Configure run at Sun Mar  9 17:02:57 2014
> [0]PETSC ERROR: Configure options --with-cc=gcc --with-fc=gfortran
> --download-f-blas-lapack --download-mpich
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: User provided function() line 0 in unknown directory
> unknown file
> application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
> [unset]: aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 59) - process 0
>
> -- 
>
> ------------------------------
> This message and its contents, including attachments are intended solely 
> for the original recipient. If you are not the intended recipient or have 
> received this message in error, please notify me immediately and delete 
> this message from your computer system. Any unauthorized use or 
> distribution is prohibited. Please consider the environment before printing 
> this email.
-------------- next part --------------
A non-text attachment was scrubbed...
Name: not available
Type: application/pgp-signature
Size: 835 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140504/8b65089d/attachment.pgp>


More information about the petsc-users mailing list