[petsc-users] SLEPc fails with POWER ITERATION method

Xujun Zhao xzhao99 at gmail.com
Fri Aug 7 11:21:40 CDT 2015


Hi all,

I am solving the max eigenvalue of a Shell matrix using SLEPc. the Shell
operation is set MATOP_MULT with user-defined function u = M*f. It works
with the Krylov-Schur and Arnoldi method, but fails when I use Power
Iteration method and several others. This is strange, because some of those
are supposed to work with any type of problem.

I also wrote a power iteration algorithm by myself, and it works well and
obtains the same results with that from Krylov_Schur. I am curious why this
doesn't work is SLEPc. The following are my code and error messages:

----------------------------------------------------------------------------------
  ierr = MatSetFromOptions(M);      CHKERRQ(ierr);
  ierr =
MatShellSetOperation(M,MATOP_MULT,(void(*)())_MatMult_Stokes);CHKERRQ(ierr);

  ierr = PetscPrintf(PETSC_COMM_WORLD,"\n");
  printf("--->test: n = %d, N = %d, rank = %d\n",n, N, (int)rank);


  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   Create the eigensolver and set various options
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);  CHKERRQ(ierr);
  ierr = EPSSetOperators(eps,M,NULL);       CHKERRQ(ierr);
  ierr = EPSSetProblemType(eps,EPS_HEP);   CHKERRQ(ierr);


  // EPSKRYLOVSCHUR(Default)/EPSARNOLDI/
  // does NOT work: EPSPOWER/EPSLANCZOS/EPSSUBSPACE
  ierr = EPSSetType(eps,EPSPOWER);          CHKERRQ(ierr);


  /*  Select portion of spectrum */

  if(option=="smallest")  // LOBPCG for smallest eigenvalue problem!
  { ierr = EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);CHKERRQ(ierr);  }
  else if(option=="largest")
  { ierr = EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); CHKERRQ(ierr);  }
  else
  { ierr = EPSSetFromOptions(eps);  CHKERRQ(ierr);  }
  // end if-else


  // Set the tolerance and maximum iteration
  ierr = EPSSetTolerances(eps,  tol,  maxits);  CHKERRQ(ierr);
  ierr = PetscPrintf(PETSC_COMM_WORLD,"EPS tol = %f, maxits =
%d\n",tol,maxits);CHKERRQ(ierr);


  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
   Solve the eigensystem and get the solution
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

  ierr = PetscPrintf(PETSC_COMM_WORLD,"EPS solve starts
...\n");CHKERRQ(ierr);
  ierr = EPSSolve(eps);CHKERRQ(ierr);

----------------------------------------------------------------------------------


EPS tol = 0.000001, maxits = 100

EPS solve starts ...

[0]PETSC ERROR: --------------------- Error Message
--------------------------------------------------------------

[0]PETSC ERROR: Wrong value of eps->which

[0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for
trouble shooting.

[0]PETSC ERROR: Petsc Release Version 3.5.4, May, 23, 2015

[0]PETSC ERROR: ./example-opt on a arch-darwin-c-opt named
mcswl121.mcs.anl.gov by xzhao Fri Aug  7 10:31:38 2015

[0]PETSC ERROR: Configure options --with-cc=gcc-4.9 --with-cxx=g++-4.9
--with-fc=gfortran-4.9 --with-cxx-dialect=C++11 --download-mpich
--download-fblaslapack --download-scalapack --download-mumps
--download-superlu_dist --download-hypre --download-ml --download-parmetis
--download-metis --download-triangle --download-chaco --download-elemental
--with-debugging=0

[0]PETSC ERROR: #1 EPSSetUp_Power() line 64 in
/Users/xzhao/software/slepc/slepc-3.5.4/src/eps/impls/power/power.c

[0]PETSC ERROR: #2 EPSSetUp() line 120 in
/Users/xzhao/software/slepc/slepc-3.5.4/src/eps/interface/epssetup.c

[0]PETSC ERROR: #3 EPSSolve() line 88 in
/Users/xzhao/software/slepc/slepc-3.5.4/src/eps/interface/epssolve.c

[0]PETSC ERROR: #4 compute_eigenvalue() line 318 in brownian_system.C
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20150807/8e0de2f2/attachment.html>


More information about the petsc-users mailing list