[petsc-users] slepc4py - Using shell matrix and explicit preconditioning matrix to solve EPS

Jose E. Roman jroman at dsic.upv.es
Wed Aug 1 04:00:23 CDT 2018


If Krylov methods work for you then no need to worry about Davidson methods. JD has many parameters that influence convergence. Also, it relies on having a good preconditioner, maybe ASM is not enough for your problem. Some time ago with we used the pARMS preconditioner, see https://doi.org/10.1016/j.cpc.2011.12.018 - but the results are problem dependent. Also, you can try GD instead of JD, which is simpler and often gives better performance. See a detailed explanation here: https://doi.org/10.1145/2543696

Jose

> El 1 ago 2018, a las 10:43, Michael Werner <michael.werner at dlr.de> escribió:
> 
> 
> Thanks for the quick reply, your suggestion worked perfect. Concerning your remark about Davidson-type methods: I already tried that, however I had problems getting them to converge to the expected eigenvalues. Even when I was using the full explicit system matrix and the options suggested in the slepc manual (-st_type precond -eps_type jd -st_pc_type asm -st_ksp_type bcgsl -st_ksp_rtol 0.001 -eps_target 0.01 -eps_target_real), the results were always "hit-and-miss", the results didn't reliably converge towards the eigenvalues suggested by the exact shift-and-invert.
> 
> Kind regards,
> Michael
> 
> Jose E. Roman writes:
> 
>> Try calling E.setUp() before replacing the matrix. But note that E.setUp() will try to build the preconditioner, so you will have to set PC=none and later change it to lu.
>> 
>> Also note that in Krylov methods with exact shift-and-invert (LU) the preconditioner matrix should be exactly the same as A-sigma*B, otherwise you may get unexpected results. Davidson-type methods allow using a different preconditioner.
>> 
>> Jose
>> 
>>> El 1 ago 2018, a las 10:11, Michael Werner <michael.werner at dlr.de> escribió:
>>> Hello,
>>> I'm trying to find the smallest eigenvalues of a linear system created by CFD simulations. To reduce memory requirements, I want to use a shell matrix (A_Shell) to provide the matrix-vector product, and a lower-order explicit matrix (P) as preconditioner. As I'm solving a generalized problem, I'm also passing a second explicit matrix (B) to the EPS. However, I don't see a way to provide the explicit preconditioning matrix to the EPS.
>>> I'm using the following options: -st_type sinvert -st_pc_type lu -st_matmode shell -eps_target 0.01 -eps_target_real  My Python code looks like this:
>>> E = SLEPc.EPS().create(comm=PETSc.COMM_WORLD)
>>> E.setOperators(A_Shell, B)
>>> E.setProblemType(SLEPc.EPS.ProblemType.GNHEP)
>>> st = E.getST()
>>> ksp = st.getKSP()
>>> pc = ksp.getPC()
>>> pc.setOperators(A_Shell, P)
>>> E.setFromOptions()
>>> E.solve()
>>> If I run this, I get the following error:
>>>> Error: error code 92
>>>> [2] EPSSolve() line 135 in /home/conda/feedstock_root/build_artifacts/slepc_1531127798016/work/src/eps/interface/epssolve.c
>>>> [2] EPSSetUp() line 263 in /home/conda/feedstock_root/build_artifacts/slepc_1531127798016/work/src/eps/interface/epssetup.c
>>>> [2] STSetUp() line 271 in /home/conda/feedstock_root/build_artifacts/slepc_1531127798016/work/src/sys/classes/st/interface/stsolve.c
>>>> [2] STSetUp_Sinvert() line 132 in /home/conda/feedstock_root/build_artifacts/slepc_1531127798016/work/src/sys/classes/st/impls/sinvert/sinvert.c
>>>> [2] KSPSetUp() line 381 in /home/conda/feedstock_root/build_artifacts/petsc_1530637062977/work/src/ksp/ksp/interface/itfunc.c
>>>> [2] PCSetUp() line 923 in /home/conda/feedstock_root/build_artifacts/petsc_1530637062977/work/src/ksp/pc/interface/precon.c
>>>> [2] PCSetUp_LU() line 113 in /home/conda/feedstock_root/build_artifacts/petsc_1530637062977/work/src/ksp/pc/impls/factor/lu/lu.c
>>>> [2] MatGetFactor() line 4332 in /home/conda/feedstock_root/build_artifacts/petsc_1530637062977/work/src/mat/interface/matrix.c
>>>> [2] See http://www.mcs.anl.gov/petsc/documentation/linearsolvertable.html for possible LU and Cholesky solvers
>>>> [2] MatSolverType mumps does not support matrix type shell
>>> I also tried setting up the preconditioner in advance, but that didn't work either and produced the same errors as above:
>>> pc = PETSc.PC().create(comm=PETSc.COMM_WORLD)
>>> pc.setOperators(A_Shell, P)
>>> pc.setFromOptions()
>>> pc.setUp()
>>> E = SLEPc.EPS().create(comm=PETSc.COMM_WORLD)
>>> E.setOperators(A_Shell, B)
>>> E.setProblemType(SLEPc.EPS.ProblemType.GNHEP)
>>> ksp = E.getKSP()
>>> ksp.setPC(pc)
>>> E.setFromOptions()
>>> E.solve()
>>> If I simply try to solve a KSP with those two matrices, it works without any problems.
> 



More information about the petsc-users mailing list