[petsc-users] problem with CISS in SLEPC

Jose E. Roman jroman at dsic.upv.es
Thu Jul 11 12:56:16 CDT 2019


Internally the solver builds a subspace of dimension L*M, which defaults to 128.
See the details here http://slepc.upv.es/documentation/reports/str11.pdf
I guess the method is not prepared for smallish matrices. Try reducing the values of L and M.

Jose


> El 11 jul 2019, a las 19:43, Matthew Knepley via petsc-users <petsc-users at mcs.anl.gov> escribió:
> 
> On Thu, Jul 11, 2019 at 12:37 PM Povolotskyi, Mykhailo via petsc-users <petsc-users at mcs.anl.gov> wrote:
> Hello,
> 
> I want to use CISS for my eigenvalue problem. To test it I tried a 
> simple (2x2) matrix case, but the code failed with the following message:
> 
> [0]PETSC ERROR: --------------------- Error Message 
> --------------------------------------------------------------
> [0]PETSC ERROR: Error in external library
> [0]PETSC ERROR: Error in LAPACK subroutine gesvd: info=127
> 
> This does not make sense. The error is supposed to be the number of non-converged diagonals, but
> you only have a 2x2 problem, so something else is wrong.
> 
>    Matt
>  
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html 
> for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.8.4, Mar, 24, 2018
> [0]PETSC ERROR: var1 on a linux-complex named brown-a027.rcac.purdue.edu 
> by mpovolot Thu Jul 11 13:28:21 2019
> [0]PETSC ERROR: Configure options --with-scalar-type=complex --with-x=0 
> --with-hdf5 --download-hdf5=1 --with-single-library=1 --with-pic=1 
> --with-shared-libraries=0 --with-log=0 --with-clanguage=C++ 
> --CXXFLAGS="-fopenmp -fPIC" --CFLAGS="-fopenmp -fPIC" --with-fortran=0 
> --FFLAGS="-fopenmp -fPIC" --with-debugging=0 --with-cc=mpicc 
> --with-fc=mpif90 --with-cxx=mpicxx COPTFLAGS= CXXOPTFLAGS= FOPTFLAGS= 
> --download-metis=1 --download-parmetis=1 
> --with-valgrind-dir=/apps/brown/valgrind/3.13.0_gcc-4.8.5 
> --download-mumps=1 --with-fortran-kernels=0 --download-superlu_dist=1 
> --with-blaslapack-lib="-L/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/lib/intel64 
> -lmkl_intel_lp64 -lmkl_gnu_thread -lmkl_core " 
> --with-blacs-lib=/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/lib/intel64/libmkl_blacs_intelmpi_lp64.so 
> --with-blacs-include=/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/include 
> --with-scalapack-lib="-Wl,-rpath,/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/lib/intel64 
> -L/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/lib/intel64 
> -lmkl_gf_lp64 -lmkl_gnu_thread -lmkl_core  -lpthread 
> -L/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/lib/intel64 
> -lmkl_scalapack_lp64 -lmkl_blacs_intelmpi_lp64" 
> --with-scalapack-include=/apps/cent7/intel/compilers_and_libraries_2017.1.132/linux/mkl/include
> [0]PETSC ERROR: #1 SVD_S() line 611 in 
> /depot/kildisha/apps/brown/nemo5/libs/slepc/build-cplx/src/eps/impls/ciss/ciss.c
> [0]PETSC ERROR: #2 EPSSolve_CISS() line 1028 in 
> /depot/kildisha/apps/brown/nemo5/libs/slepc/build-cplx/src/eps/impls/ciss/ciss.c
> [0]PETSC ERROR: #3 EPSSolve() line 147 in 
> /depot/kildisha/apps/brown/nemo5/libs/slepc/build-cplx/src/eps/interface/epssolve.c
> 
> 
> Here is the code:
> 
> void test1()
> {
>    Mat matrix_petsc;
> MatCreateDense(MPI_COMM_SELF,2,2,PETSC_DECIDE,PETSC_DECIDE,NULL,&matrix_petsc);
> MatSetValue(matrix_petsc,0,0,complex<double>(1,0),INSERT_VALUES);
> MatSetValue(matrix_petsc,1,1,complex<double>(2,0),INSERT_VALUES);
> MatSetValue(matrix_petsc,0,1,complex<double>(0,0),INSERT_VALUES);
> MatSetValue(matrix_petsc,1,0,complex<double>(0,0),INSERT_VALUES);
>    MatAssemblyBegin(matrix_petsc, MAT_FINAL_ASSEMBLY);
>    MatAssemblyEnd(matrix_petsc, MAT_FINAL_ASSEMBLY);
> 
>   {
>        EPS            eps;
>        EPSType        type;
>        RG             rg;
> 
>        EPSCreate(MPI_COMM_SELF,&eps);
>        EPSSetOperators( eps,matrix_petsc,NULL);
>        EPSSetType(eps,EPSCISS);
>        EPSSetProblemType(eps, EPS_NHEP);
> 
>        double real_min = -10;
>        double real_max =  10;
> 
>        double imag_min = -10;
>        double imag_max = 10;
> 
>        EPSGetRG(eps,&rg);
>        RGSetType(rg,RGINTERVAL);
> RGIntervalSetEndpoints(rg,real_min,real_max,imag_min,imag_max);
> 
>        EPSSolve(eps);
> 
>        EPSDestroy(&eps);
>    }
> 
> 
>    MatDestroy(&matrix_petsc);
> 
> }
> 
> 
> Do you see anything wrong in the way I'm calling SLEPC functions?
> 
> Thank you,
> 
> Michael.
> 
> 
> 
> -- 
> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
> -- Norbert Wiener
> 
> https://www.cse.buffalo.edu/~knepley/



More information about the petsc-users mailing list