Actual source code: test14.c

slepc-3.21.2 2024-09-25
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 in NEP.\n\n"
 12:   "The command line options are:\n"
 13:   "  -n <n>, where <n> = matrix dimension.\n";

 15: /*
 16:    Solve T(lambda)x=0 with T(lambda) = -D+sqrt(lambda)*I
 17:       where D is the Laplacian operator in 1 dimension
 18: */

 20: #include <slepcnep.h>

 22: /*
 23:   MyConvergedRel - Convergence test relative to the norm of D (given in ctx).
 24: */
 25: PetscErrorCode MyConvergedRel(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
 26: {
 27:   PetscReal norm = *(PetscReal*)ctx;

 29:   PetscFunctionBegin;
 30:   *errest = res/norm;
 31:   PetscFunctionReturn(PETSC_SUCCESS);
 32: }

 34: int main(int argc,char **argv)
 35: {
 36:   NEP            nep;             /* nonlinear eigensolver context */
 37:   Mat            A[2];
 38:   PetscInt       n=100,Istart,Iend,i;
 39:   PetscBool      terse;
 40:   PetscReal      norm;
 41:   FN             f[2];
 42:   PetscScalar    coeffs;

 44:   PetscFunctionBeginUser;
 45:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 46:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 47:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%" PetscInt_FMT "\n\n",n));

 49:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 50:      Create nonlinear eigensolver, define problem in split form
 51:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 53:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));

 55:   /* Create matrices */
 56:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A[0]));
 57:   PetscCall(MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n));
 58:   PetscCall(MatSetFromOptions(A[0]));
 59:   PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
 60:   for (i=Istart;i<Iend;i++) {
 61:     if (i>0) PetscCall(MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES));
 62:     if (i<n-1) PetscCall(MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES));
 63:     PetscCall(MatSetValue(A[0],i,i,-2.0,INSERT_VALUES));
 64:   }
 65:   PetscCall(MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY));
 66:   PetscCall(MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY));

 68:   PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&A[1]));

 70:   /* Define functions */
 71:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[0]));
 72:   PetscCall(FNSetType(f[0],FNRATIONAL));
 73:   coeffs = 1.0;
 74:   PetscCall(FNRationalSetNumerator(f[0],1,&coeffs));
 75:   PetscCall(FNCreate(PETSC_COMM_WORLD,&f[1]));
 76:   PetscCall(FNSetType(f[1],FNSQRT));
 77:   PetscCall(NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN));

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:                    Set some options and solve
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 83:   PetscCall(NEPSetTarget(nep,1.1));

 85:   /* setup convergence test relative to the norm of D */
 86:   PetscCall(MatNorm(A[0],NORM_1,&norm));
 87:   PetscCall(NEPSetConvergenceTestFunction(nep,MyConvergedRel,&norm,NULL));

 89:   PetscCall(NEPSetFromOptions(nep));
 90:   PetscCall(NEPSolve(nep));

 92:   /* show detailed info unless -terse option is given by user */
 93:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
 94:   if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL));
 95:   else {
 96:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
 97:     PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
 98:     PetscCall(NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
 99:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
100:   }
101:   PetscCall(NEPDestroy(&nep));
102:   PetscCall(MatDestroy(&A[0]));
103:   PetscCall(MatDestroy(&A[1]));
104:   PetscCall(FNDestroy(&f[0]));
105:   PetscCall(FNDestroy(&f[1]));
106:   PetscCall(SlepcFinalize());
107:   return 0;
108: }

110: /*TEST

112:    test:
113:       suffix: 1
114:       args: -nep_type slp -nep_nev 2 -terse

116: TEST*/