[petsc-users] Help with generalized eigenvalue problem
Jose E. Roman
jroman at dsic.upv.es
Tue Jul 7 11:01:49 CDT 2020
I don't know what options you are using, shift-and-invert or not, which, target. Anyway, if A and B have a common null space, then A-sigma*B will always be singular, so the KSP solver will have to deal with a singular system.
If you pass the nullspace vectors to the EPS solver via https://slepc.upv.es/documentation/current/docs/manualpages/EPS/EPSSetDeflationSpace.html then these vectors will also be attached to the KSP solver.
Also, since you have MUMPS, try adding
-st_pc_factor_mat_solver_type mumps
(see section 3.4.1 of the SLEPc users manual). MUMPS should be quite robust with respect to singular linear systems.
Jose
> El 7 jul 2020, a las 17:40, Aulisa, Eugenio <Eugenio.Aulisa at ttu.edu> escribió:
>
> Can somebody help me figure it out the solver options for the attached generalized eigenvalue problem?
>
> A x = lambda B x
>
> where A and B have the same null space, with 5 zero eigenvalues.
>
> The only non indefinite eigenvalue should be 41.1892
>
> The problem is Hermite semi-positive so automatically slepc should purify it, but to be sure I added
>
> EPSSetPurify(eps, PETSC_TRUE);
>
> I guess with my non-solver-options slepc is still trying to do the inverse of the full
> B matrix with an LU decomposition and gets zero pivot.
>
> This is what I get as a PETSC error Message
>
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Zero pivot in LU factorization: https://www.mcs.anl.gov/petsc/documentation/faq.html#zeropivot
> [0]PETSC ERROR: Bad LU factorization
> [0]PETSC ERROR: See https://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Development GIT revision: v3.11.3-2263-gce77f2ed1a GIT Date: 2019-09-26 13:31:14 -0500
> [0]PETSC ERROR: Nitsche_ex4a on a arch-linux2-c-opt named linux-8biu by eaulisa Tue Jul 7 10:17:55 2020
> [0]PETSC ERROR: Configure options --with-debugging=0 --with-x=1 COPTFLAGS="-O3 -march=native -mtune=native" CXXOPTFLAGS="-O3 -march=native -mtune=native" FOPTFLAGS="-O3 -march=native -mtune=native" --download-openmpi=1 --download-fblaslapack=1 --download-hdf5=1 --download-metis=1 --download-parmetis=1 --with-shared-libraries=1 --download-blacs=1 --download-scalapack=1 --download-mumps=1 --download-suitesparse
> [0]PETSC ERROR: #1 MatLUFactor_SeqDense() line 633 in /home/eaulisa/software/petsc/src/mat/impls/dense/seq/dense.c
> [0]PETSC ERROR: #2 MatLUFactorNumeric_SeqDense() line 432 in /home/eaulisa/software/petsc/src/mat/impls/dense/seq/dense.c
> [0]PETSC ERROR: #3 MatLUFactorNumeric() line 3056 in /home/eaulisa/software/petsc/src/mat/interface/matrix.c
> [0]PETSC ERROR: #4 PCSetUp_LU() line 126 in /home/eaulisa/software/petsc/src/ksp/pc/impls/factor/lu/lu.c
> [0]PETSC ERROR: #5 PCSetUp() line 894 in /home/eaulisa/software/petsc/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #6 KSPSetUp() line 377 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #7 STSetUp_Shift() line 119 in /home/eaulisa/software/slepc/src/sys/classes/st/impls/shift/shift.c
> [0]PETSC ERROR: #8 STSetUp() line 271 in /home/eaulisa/software/slepc/src/sys/classes/st/interface/stsolve.c
> [0]PETSC ERROR: #9 EPSSetUp() line 273 in /home/eaulisa/software/slepc/src/eps/interface/epssetup.c
> [0]PETSC ERROR: #10 EPSSolve() line 136 in /home/eaulisa/software/slepc/src/eps/interface/epssolve.c
>
>
>
> Thanks
> Eugenio
>
>
>
>
> Eugenio Aulisa
>
> Department of Mathematics and Statistics,
> Texas Tech University
> Lubbock TX, 79409-1042
> room: 226
> http://www.math.ttu.edu/~eaulisa/
> phone: (806) 834-6684
> fax: (806) 742-1112
>
>
>
> <genEigen.cpp>
More information about the petsc-users
mailing list