[petsc-users] Problems about PCtype bjacobi

Mark Adams mfadams at lbl.gov
Wed Nov 7 10:01:54 CST 2018


On Wed, Nov 7, 2018 at 10:16 AM Yingjie Wu via petsc-users <
petsc-users at mcs.anl.gov> wrote:

> Dear Petsc developer:
> Hi,
> Recently, I'm solving the problems of nonlinear systems of PDEs, I
> encountered some problems about precondition and wanted to seek help.
>
> 1.I set the precondition matrix in SNES as MPIAIJ in the program, and then
> use Matrix Free method to solve my problem. The log information of the
> program is as follows:
>
> SNES Object: 1 MPI processes
>   type: newtonls
>   maximum iterations=50, maximum function evaluations=100000000
>   tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>   total number of linear solver iterations=177
>   total number of function evaluations=371
>   norm schedule ALWAYS
>   SNESLineSearch Object: 1 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: 1 MPI processes
>     type: gmres
>       restart=30, using Classical (unmodified) Gram-Schmidt
> Orthogonalization with no iterative refinement
>       happy breakdown tolerance 1e-30
>     maximum iterations=10000, initial guess is zero
>     tolerances:  relative=0.01, absolute=1e-50, divergence=10000.
>     left preconditioning
>     using PRECONDITIONED norm type for convergence test
>   PC Object: 1 MPI processes
>     type: bjacobi
>       number of blocks = 1
>       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: bjacobi
>           number of blocks = 1
>           Local solve is same for all blocks, in the following KSP and PC
> objects:
>           KSP Object: (sub_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_sub_) 1 MPI processes
>            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: 1 MPI processes
>                     type: seqaij
>                     rows=961, cols=961
>                    package used to perform factorization: petsc
>                     total: nonzeros=4129, allocated nonzeros=4129
>                     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=961, cols=961
>               total: nonzeros=4129, allocated nonzeros=4805
>               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: mpiaij
>           rows=961, cols=961
>           total: nonzeros=4129, allocated nonzeros=9610
>          total number of mallocs used during MatSetValues calls =0
>             not using I-node (on process 0) routines
>     linear system matrix followed by preconditioner matrix:
>     Mat Object: 1 MPI processes
>       type: mffd
>       rows=961, cols=961
>         Matrix-free approximation:
>           err=1.49012e-08 (relative error in function evaluation)
>           Using wp compute h routine
>               Does not compute normU
>     Mat Object: 1 MPI processes
>       type: mpiaij
>       rows=961, cols=961
>       total: nonzeros=4129, allocated nonzeros=9610
>       total number of mallocs used during MatSetValues calls =0
>         not using I-node (on process 0) routines
>
> Although parallel matrix is used, it runs on a single processor. Because
> of the use of parallel matrices, the overall precondition scheme should be
> bjacobi, and then build a KSP for each block (there is only one block in my
> program). Therefore, there will be a sub KSP object in the PC information.
> But in the above information, a subsubksp is also embedded in the sub KSP
> object. I don't understand the reason for this KSP. Please help me answer.
>

It looks like you are setting the sub PC type to bjacobi and it has a sub
(sub) PC.


>
> 2. Bjacobi is a precondition method in theory. Why is there a subsystem of
> linear equations solver object ksp? In my understanding, precondition is a
> matrix decomposition method, and there is no linear equations,. Therefore,
> there should be no subKSP here.
>
> 3.when I  assemble precondition matrix, I used MatSetVales.
>
> idleft  = row - 1 ;
>
>       idright = row + 1 ;
>
>       idup    = row + 20;
>
>       iddown  = row - 20;
>
>
>
>       /*phi1 field*/
>
>       v[1] = - 1.0 / ( dx * dx / (2* 1.267) + dx * dx / (2* 1.267 ) ) ;
> col[1] = idleft;
>
>       v[2] = - 1.0 / ( dx * dx / (2* 1.267) + dx * dx / (2* 1.267 ) ) ;
> col[2] = idright;
>
>       v[3] = - 1.0/  ( dy * dy / (2* 1.267) + dy * dy / (2* 1.267 ) ) ;
> col[3] = iddown;
>
>       v[4] = - 1.0 / ( dy * dy / (2* 1.267) + dy * dy / (2* 1.267 ) ) ;
> col[4] = idup;
>
>       v[0] = v[1] + v[2] + v[3] + v[4] + v[0] ;   col[0] = row;
>
>       ierr  = MatSetValues(B,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
>
>
There must be junk in "col" somehow. I don't see anything wrong with this
code.  Print out the values of "col" before the MatSetValues call.


> The wrong information is:
>
> [0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> [0]PETSC ERROR: Argument out of range
> [0]PETSC ERROR: Column too large: col 1077215232 max 960
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html
> for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.10.1, Sep, 26, 2018
> [0]PETSC ERROR: ./ex217 on a arch-linux2-c-debug named yjwu-XPS-8910 by
> yjwu Wed Nov  7 09:52:22 2018
> [0]PETSC ERROR: Configure options --with-cc=gcc --with-cxx=g++
> --with-fc=gfortran --download-mpich --download-fblaslapack
> [0]PETSC ERROR: #1 MatSetValues_SeqAIJ() line 442 in
> /home/yjwu/petsc-3.10.1/src/mat/impls/aij/seq/aij.c
> [0]PETSC ERROR: #2 MatSetValues() line 1349 in
> /home/yjwu/petsc-3.10.1/src/mat/interface/matrix.c
> [0]PETSC ERROR: #3 FormJacobian() line 272 in
> /home/yjwu/petsc-3.10.1/src/snes/examples/tutorials/ex217.c
> [0]PETSC ERROR: #4 SNESComputeJacobian() line 2555 in
> /home/yjwu/petsc-3.10.1/src/snes/interface/snes.c
> [0]PETSC ERROR: #5 SNESSolve_NEWTONLS() line 222 in
> /home/yjwu/petsc-3.10.1/src/snes/impls/ls/ls.c
> [0]PETSC ERROR: #6 SNESSolve() line 4396 in
> /home/yjwu/petsc-3.10.1/src/snes/interface/snes.c
> [0]PETSC ERROR: #7 main() line 108 in
> /home/yjwu/petsc-3.10.1/src/snes/examples/tutorials/ex217.c
>
> This seems to be a very obvious question, but I don't know how to solve
> it.
>
> Thanks for your continuous help,
> Yingjie
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20181107/4ab165fc/attachment.html>


More information about the petsc-users mailing list