Actual source code: test16.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: static char help[] = "Tests a user-defined convergence test.\n\n";

 24: #include <slepceps.h>


 29: /*
 30:   MyConvergedAbsolute - Bizarre convergence test that requires more accuracy
 31:   to positive eigenvalues compared to negative ones.
 32: */
 33: PetscErrorCode MyConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 34: {
 36:   *errest = (PetscRealPart(eigr)<0.0)?res:100*res;
 37:   return(0);
 38: }

 42: int main(int argc,char **argv)
 43: {
 44:   Mat            A;           /* problem matrix */
 45:   EPS            eps;         /* eigenproblem solver context */
 46:   PetscInt       n=30,i,Istart,Iend;

 49:   SlepcInitialize(&argc,&argv,(char*)0,help);
 50:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 51:   PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal Eigenproblem, n=%D\n\n",n);

 53:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 54:      Compute the operator matrix that defines the eigensystem, Ax=kx
 55:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 57:   MatCreate(PETSC_COMM_WORLD,&A);
 58:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 59:   MatSetFromOptions(A);
 60:   MatSetUp(A);

 62:   MatGetOwnershipRange(A,&Istart,&Iend);
 63:   for (i=Istart;i<Iend;i++) {
 64:     if (i>0) { MatSetValue(A,i,i-1,-1.0,INSERT_VALUES); }
 65:     if (i<n-1) { MatSetValue(A,i,i+1,-1.0,INSERT_VALUES); }
 66:   }
 67:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 68:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 69:   MatShift(A,-1e-3);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:                         Create the eigensolver
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 74:   EPSCreate(PETSC_COMM_WORLD,&eps);
 75:   EPSSetOperators(eps,A,NULL);
 76:   EPSSetProblemType(eps,EPS_HEP);
 77:   /* set user-defined convergence test */
 78:   EPSSetConvergenceTestFunction(eps,MyConvergedAbsolute,NULL,NULL);
 79:   EPSSetFromOptions(eps);

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:                           Solve the problem
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 84:   EPSSolve(eps);
 85:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

 87:   EPSDestroy(&eps);
 88:   MatDestroy(&A);
 89:   SlepcFinalize();
 90:   return ierr;
 91: }