Actual source code: test14.c
slepc-3.18.2 2023-01-26
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: return 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;
44: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%" PetscInt_FMT "\n\n",n);
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Create nonlinear eigensolver, define problem in split form
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: NEPCreate(PETSC_COMM_WORLD,&nep);
54: /* Create matrices */
55: MatCreate(PETSC_COMM_WORLD,&A[0]);
56: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
57: MatSetFromOptions(A[0]);
58: MatSetUp(A[0]);
59: MatGetOwnershipRange(A[0],&Istart,&Iend);
60: for (i=Istart;i<Iend;i++) {
61: if (i>0) MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES);
62: if (i<n-1) MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES);
63: MatSetValue(A[0],i,i,-2.0,INSERT_VALUES);
64: }
65: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
68: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&A[1]);
70: /* Define functions */
71: FNCreate(PETSC_COMM_WORLD,&f[0]);
72: FNSetType(f[0],FNRATIONAL);
73: coeffs = 1.0;
74: FNRationalSetNumerator(f[0],1,&coeffs);
75: FNCreate(PETSC_COMM_WORLD,&f[1]);
76: FNSetType(f[1],FNSQRT);
77: NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Set some options and solve
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83: NEPSetTarget(nep,1.1);
85: /* setup convergence test relative to the norm of D */
86: MatNorm(A[0],NORM_1,&norm);
87: NEPSetConvergenceTestFunction(nep,MyConvergedRel,&norm,NULL);
89: NEPSetFromOptions(nep);
90: NEPSolve(nep);
92: /* show detailed info unless -terse option is given by user */
93: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
94: if (terse) NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
95: else {
96: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
97: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
98: NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
99: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
100: }
101: NEPDestroy(&nep);
102: MatDestroy(&A[0]);
103: MatDestroy(&A[1]);
104: FNDestroy(&f[0]);
105: FNDestroy(&f[1]);
106: SlepcFinalize();
107: return 0;
108: }
110: /*TEST
112: test:
113: suffix: 1
114: args: -nep_type slp -nep_nev 2 -terse
116: TEST*/