Actual source code: test16.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

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

 13: #include <slepceps.h>

 15: /*
 16:   MyConvergedAbsolute - Bizarre convergence test that requires more accuracy
 17:   to positive eigenvalues compared to negative ones.
 18: */
 19: PetscErrorCode MyConvergedAbsolute(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 20: {
 21:   *errest = (PetscRealPart(eigr)<0.0)?res:100*res;
 22:   return 0;
 23: }

 25: int main(int argc,char **argv)
 26: {
 27:   Mat            A;           /* problem matrix */
 28:   EPS            eps;         /* eigenproblem solver context */
 29:   PetscInt       n=30,i,Istart,Iend;

 32:   SlepcInitialize(&argc,&argv,(char*)0,help);
 33:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 34:   PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal Eigenproblem, n=%" PetscInt_FMT "\n\n",n);

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 37:      Compute the operator matrix that defines the eigensystem, Ax=kx
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   MatCreate(PETSC_COMM_WORLD,&A);
 41:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 42:   MatSetFromOptions(A);
 43:   MatSetUp(A);

 45:   MatGetOwnershipRange(A,&Istart,&Iend);
 46:   for (i=Istart;i<Iend;i++) {
 47:     if (i>0) MatSetValue(A,i,i-1,-1.0,INSERT_VALUES);
 48:     if (i<n-1) MatSetValue(A,i,i+1,-1.0,INSERT_VALUES);
 49:   }
 50:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 51:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 52:   MatShift(A,-1e-3);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:                         Create the eigensolver
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 57:   EPSCreate(PETSC_COMM_WORLD,&eps);
 58:   EPSSetOperators(eps,A,NULL);
 59:   EPSSetProblemType(eps,EPS_HEP);
 60:   /* set user-defined convergence test */
 61:   EPSSetConvergenceTestFunction(eps,MyConvergedAbsolute,NULL,NULL);
 62:   EPSSetFromOptions(eps);

 64:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 65:                           Solve the problem
 66:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 67:   EPSSolve(eps);
 68:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

 70:   EPSDestroy(&eps);
 71:   MatDestroy(&A);
 72:   SlepcFinalize();
 73:   return 0;
 74: }

 76: /*TEST

 78:    test:
 79:       suffix: 1
 80:       args: -n 200 -eps_nev 6 -eps_ncv 24 -eps_smallest_magnitude
 81:       requires: !single

 83: TEST*/