[petsc-users] Use block Jacobi preconditioner with SNES

Smith, Barry F. bsmith at mcs.anl.gov
Wed Aug 29 10:10:12 CDT 2018


   You need to add the line 

  ierr = MatSetVariableBlockSizes(J,3,lens);CHKERRQ(ierr);

with the number of blocks replacing 3 and the sizes of each block in lens

   Barry


> On Aug 29, 2018, at 4:27 AM, Ali Reza Khaz'ali <arkhazali at cc.iut.ac.ir> wrote:
> 
> Barry,
> 
> Thanks a lot for your efforts. Using PCVPBJACOBI causes SNESSolve to throw the following error (the code and the system being simulated is the same, just PCBJACOBI changed to PCVPBJACOBI):
> 
> 
> [0]PETSC ERROR: --------------------- Error Message ---------------------------------------
> -----------------------
> [0]PETSC ERROR: Nonconforming object sizes
> [0]PETSC ERROR: Total blocksizes 0 doesn't match number matrix rows 108000
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Development GIT revision: v3.9.3-1238-g434eb0c8b5  GIT Date: 2018-08-28 19:19:26 -0500
> [0]PETSC ERROR: E:\Documents\Visual Studio 2015\Projects\compsim\x64\Release\compsim.exe on a  named ALIREZA-PC by AliReza Wed Aug 29 13:46:38 2018
> [0]PETSC ERROR: Configure options --prefix=/home/alireza/PetscGit --with-mkl_pardiso-dir=/cygdrive/E/Program_Files_x86/IntelSWTools/compilers_and_libraries/windows/mkl --
> with-hypre-include=/cygdrive/E/hypre-2.11.2/Builds/Bins/include --with-hypre-lib=/cygdrive/E/hypre-2.11.2/Builds/Bins/lib/HYPRE.lib --with-ml-include=/cygdrive/E/Trilinos
> -master/Bins/include --with-ml-lib=/cygdrive/E/Trilinos-master/Bins/lib/ml.lib ظ€ôwith-openmp --with-cc="win32fe icl" --with-fc="win32fe ifort" --with-mpi-include=/cygdri
> ve/E/Program_Files_x86/IntelSWTools/compilers_and_libraries/windows/mpi/intel64/include --with-mpi-lib=/cygdrive/E/Program_Files_x86/IntelSWTools/compilers_and_libraries/
> windows/mpi/intel64/lib/impi.lib --with-mpi-mpiexec=/cygdrive/E/Program_Files_x86/IntelSWTools/compilers_and_libraries/windows/mpi/intel64/bin/mpiexec.exe --with-debuggin
> g=0 --with-blas-lib=/cygdrive/E/Program_Files_x86/IntelSWTools/compilers_and_libraries/windows/mkl/lib/intel64_win/mkl_rt.lib --with-lapack-lib=/cygdrive/E/Program_Files_
> x86/IntelSWTools/compilers_and_libraries/windows/mkl/lib/intel64_win/mkl_rt.lib -CFLAGS="-O2 -MT -wd4996 -Qopenmp" -CXXFLAGS="-O2 -MT -wd4996 -Qopenmp" -FFLAGS="-MT -O2 -
> Qopenmp"
> [0]PETSC ERROR: #1 MatInvertVariableBlockDiagonal_SeqAIJ() line 1569 in E:\cygwin64\home\alireza\PETSc\src\mat\impls\aij\seq\aij.c
> [0]PETSC ERROR: #2 MatInvertVariableBlockDiagonal() line 10482 in E:\cygwin64\home\alireza\PETSc\src\mat\interface\matrix.c
> [0]PETSC ERROR: #3 PCSetUp_VPBJacobi() line 120 in E:\cygwin64\home\alireza\PETSc\src\ksp\pc\impls\vpbjacobi\vpbjacobi.c
> [0]PETSC ERROR: #4 PCSetUp() line 932 in E:\cygwin64\home\alireza\PETSc\src\ksp\pc\interface\precon.c
> [0]PETSC ERROR: #5 KSPSetUp() line 381 in E:\cygwin64\home\alireza\PETSc\src\ksp\ksp\interface\itfunc.c
> [0]PETSC ERROR: #6 KSPSolve() line 612 in E:\cygwin64\home\alireza\PETSc\src\ksp\ksp\interface\itfunc.c
> [0]PETSC ERROR: #7 SNESSolve_NEWTONLS() line 224 in E:\cygwin64\home\alireza\PETSc\src\snes\impls\ls\ls.c
> [0]PETSC ERROR: #8 SNESSolve() line 4355 in E:\cygwin64\home\alireza\PETSc\src\snes\interface\snes.c
> 
> Many thanks,
> Ali
> 
> On 8/29/2018 2:24 AM, Smith, Barry F. wrote:
>>   Ali,
>> 
>>     In the branch barry/feature-PCVPBJACOBI (see src/snes/examples/tutorials/ex5.c) I have implemented PCVPBJACOBI that is a point block Jacobi with variable size blocks. It has very little testing.
>> 
>>     Could you please try your code with it; you should see very very similar convergence history with this branch and the previous approach you used (the only numerical difference between the two approaches is that this one uses a dense LU factorization for each small block while the previous used a sparse factorization). This new one should also be slightly faster in the KSPSolve().
>> 
>>    Thanks
>> 
>>     Barry
>> 
>> 
>> 
>> 
>>> On Aug 27, 2018, at 4:04 PM, Jed Brown <jed at jedbrown.org> wrote:
>>> 
>>> "Smith, Barry F." <bsmith at mcs.anl.gov> writes:
>>> 
>>>>    I have added new functionality within SNES to allow one to do as you desire in the branch barry/feature-snessetkspsetupcallback. Please see src/snes/examples/tests/ex3.c
>>> Note that this example is block Jacobi with O(1) sparse blocks per
>>> process, not variable-sized point-block Jacobi which I think is what Ali
>>> had in mind.
>>> 
>>>>    Please let us know if it works for you and if you have any questions or problems.
>>>> 
>>>>    Barry
>>>> 
>>>>   Note that it has to be handled by a callback called from within the SNES solver because that is the first time the matrix exists in a form that the block sizes may be set.
>>>> 
>>>>    Sorry for the runaround with so many emails.
>>>> 
>>>> 
>>>>> On Aug 24, 2018, at 3:48 PM, Ali Reza Khaz'ali <arkhazali at cc.iut.ac.ir> wrote:
>>>>> 
>>>>> Hi,
>>>>> 
>>>>> I am trying to use block Jacobi preconditioner in SNES (SNESNEWTONLS). However, PCBJacobiGetSubKSP function returns an error stating "Object is in wrong state, Must call KSPSetUp() or PCSetUp() first". When I add KSPSetUp, I got and error from them as: "Must call DMShellSetGlobalVector() or DMShellSetCreateGlobalVector()", and if PCSetUp is added, "Object is in wrong state, Matrix must be set first" error is printed.
>>>>> 
>>>>> Below is a part of my code. It is run serially. Any help is much appreciated.
>>>>> 
>>>>>    ierr = SNESGetKSP(snes, &Petsc_ksp);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = KSPGetPC(Petsc_ksp, &Petsc_pc);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = KSPSetTolerances(Petsc_ksp, 1.e-3, 1e-3, PETSC_DEFAULT, 20000);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = SNESSetTolerances(snes, 1e-1, 1e-1, 1e-1, 2000, 2000);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = SNESSetType(snes, SNESNEWTONLS);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = KSPSetType(Petsc_ksp, KSPGMRES);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = PCSetType(Petsc_pc, PCBJACOBI);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = PCSetType(Petsc_pc, PCBJACOBI);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = PCBJacobiSetTotalBlocks(Petsc_pc, 2*Nx*Ny*Nz, SadeqSize);
>>>>>    CHKERRQ(ierr);
>>>>> 
>>>>>    SNESSetUp(snes);
>>>>>    CHKERRQ(ierr);
>>>>>    ierr = PCBJacobiGetSubKSP(Petsc_pc, &nLocal, &Firstly, &subKSP);
>>>>>    CHKERRQ(ierr);
>>>>> 
>>>>>    for (i = 0; i < nLocal; i++) {
>>>>>        ierr = KSPGetPC(subKSP[i], &SubPc);
>>>>>        CHKERRQ(ierr);
>>>>>        ierr = PCSetType(SubPc, PCLU);
>>>>>        CHKERRQ(ierr);
>>>>>        ierr = PCFactorSetMatSolverPackage(SubPc, "mkl_pardiso");
>>>>>        CHKERRQ(ierr);
>>>>>        ierr = KSPSetType(subKSP[i], KSPPREONLY);
>>>>>        CHKERRQ(ierr);
>>>>>        ierr = KSPSetTolerances(subKSP[i], 1.e-6, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
>>>>>        CHKERRQ(ierr);
>>>>>    }
>>>>>    ierr = SNESSolve(snes, NULL, Petsc_X);
>>>>>    CHKERRQ(ierr);
>>>>> 
>>>>> 
>>>>> 
>>>>> -- 
>>>>> Ali Reza Khaz’ali
>>>>> Assistant Professor of Petroleum Engineering,
>>>>> Department of Chemical Engineering
>>>>> Isfahan University of Technology
>>>>> Isfahan, Iran
>>>>> 
> 
> -- 
> Ali Reza Khaz’ali
> Assistant Professor of Petroleum Engineering,
> Department of Chemical Engineering
> Isfahan University of Technology
> Isfahan, Iran
> 



More information about the petsc-users mailing list