Actual source code: test14.c
slepc-3.16.3 2022-04-11
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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;
30: *errest = res/norm;
31: return(0);
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;
40: PetscBool terse;
41: PetscReal norm;
42: FN f[2];
43: PetscScalar coeffs;
45: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
46: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
47: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%D\n\n",n);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Create nonlinear eigensolver, define problem in split form
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: NEPCreate(PETSC_COMM_WORLD,&nep);
55: /* Create matrices */
56: MatCreate(PETSC_COMM_WORLD,&A[0]);
57: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
58: MatSetFromOptions(A[0]);
59: MatSetUp(A[0]);
60: MatGetOwnershipRange(A[0],&Istart,&Iend);
61: for (i=Istart;i<Iend;i++) {
62: if (i>0) { MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES); }
63: if (i<n-1) { MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES); }
64: MatSetValue(A[0],i,i,-2.0,INSERT_VALUES);
65: }
66: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
69: MatCreate(PETSC_COMM_WORLD,&A[1]);
70: MatSetSizes(A[1],PETSC_DECIDE,PETSC_DECIDE,n,n);
71: MatSetFromOptions(A[1]);
72: MatSetUp(A[1]);
73: MatAssemblyBegin(A[1],MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(A[1],MAT_FINAL_ASSEMBLY);
75: MatShift(A[1],1.0);
77: /* Define functions */
78: FNCreate(PETSC_COMM_WORLD,&f[0]);
79: FNSetType(f[0],FNRATIONAL);
80: coeffs = 1.0;
81: FNRationalSetNumerator(f[0],1,&coeffs);
82: FNCreate(PETSC_COMM_WORLD,&f[1]);
83: FNSetType(f[1],FNSQRT);
84: NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Set some options and solve
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: NEPSetTarget(nep,1.1);
92: /* setup convergence test relative to the norm of D */
93: MatNorm(A[0],NORM_1,&norm);
94: NEPSetConvergenceTestFunction(nep,MyConvergedRel,&norm,NULL);
96: NEPSetFromOptions(nep);
97: NEPSolve(nep);
99: /* show detailed info unless -terse option is given by user */
100: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
101: if (terse) {
102: NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
103: } else {
104: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
105: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
106: NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
107: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
108: }
109: NEPDestroy(&nep);
110: MatDestroy(&A[0]);
111: MatDestroy(&A[1]);
112: FNDestroy(&f[0]);
113: FNDestroy(&f[1]);
114: SlepcFinalize();
115: return ierr;
116: }
118: /*TEST
120: test:
121: suffix: 1
122: args: -nep_type slp -nep_nev 2 -terse
124: TEST*/