Actual source code: test14.c

slepc-3.17.0 2022-03-31
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:   *errest = res/norm;
 30:   PetscFunctionReturn(0);
 31: }

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

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

 47:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 48:      Create nonlinear eigensolver, define problem in split form
 49:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 51:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

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

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

 82:   NEPSetTarget(nep,1.1);

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

 88:   NEPSetFromOptions(nep);
 89:   NEPSolve(nep);

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

109: /*TEST

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

115: TEST*/