Actual source code: test16.c

slepc-3.22.1 2024-10-28
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:   PetscFunctionBegin;
 22:   *errest = (PetscRealPart(eigr)<0.0)?res:100*res;
 23:   PetscFunctionReturn(PETSC_SUCCESS);
 24: }

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

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

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

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

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

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

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

 70:   PetscCall(EPSDestroy(&eps));
 71:   PetscCall(MatDestroy(&A));
 72:   PetscCall(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*/