Actual source code: test10.c
slepc-3.13.2 2020-05-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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 PEP (based on ex16.c).\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
14: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
16: #include <slepcpep.h>
18: /*
19: MyConvergedRel - Convergence test relative to the norm of M (given in ctx).
20: */
21: PetscErrorCode MyConvergedRel(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
22: {
23: PetscReal norm = *(PetscReal*)ctx;
26: *errest = res/norm;
27: return(0);
28: }
30: int main(int argc,char **argv)
31: {
32: Mat M,C,K,A[3]; /* problem matrices */
33: PEP pep; /* polynomial eigenproblem solver context */
34: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
35: PetscBool flag;
36: PetscReal norm;
39: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
41: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
42: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
43: if (!flag) m=n;
44: N = n*m;
45: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: /* K is the 2-D Laplacian */
52: MatCreate(PETSC_COMM_WORLD,&K);
53: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
54: MatSetFromOptions(K);
55: MatSetUp(K);
56: MatGetOwnershipRange(K,&Istart,&Iend);
57: for (II=Istart;II<Iend;II++) {
58: i = II/n; j = II-i*n;
59: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
60: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
61: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
62: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
63: MatSetValue(K,II,II,4.0,INSERT_VALUES);
64: }
65: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
66: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
68: /* C is the 1-D Laplacian on horizontal lines */
69: MatCreate(PETSC_COMM_WORLD,&C);
70: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
71: MatSetFromOptions(C);
72: MatSetUp(C);
73: MatGetOwnershipRange(C,&Istart,&Iend);
74: for (II=Istart;II<Iend;II++) {
75: i = II/n; j = II-i*n;
76: if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
77: if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
78: MatSetValue(C,II,II,2.0,INSERT_VALUES);
79: }
80: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
81: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
83: /* M is a diagonal matrix */
84: MatCreate(PETSC_COMM_WORLD,&M);
85: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
86: MatSetFromOptions(M);
87: MatSetUp(M);
88: MatGetOwnershipRange(M,&Istart,&Iend);
89: for (II=Istart;II<Iend;II++) {
90: MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
91: }
92: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Create the eigensolver and set various options
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: PEPCreate(PETSC_COMM_WORLD,&pep);
100: A[0] = K; A[1] = C; A[2] = M;
101: PEPSetOperators(pep,3,A);
102: PEPSetProblemType(pep,PEP_HERMITIAN);
103: PEPSetDimensions(pep,4,20,PETSC_DEFAULT);
105: /* setup convergence test relative to the norm of M */
106: MatNorm(M,NORM_1,&norm);
107: PEPSetConvergenceTestFunction(pep,MyConvergedRel,&norm,NULL);
108: PEPSetFromOptions(pep);
110: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111: Solve the eigensystem
112: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
114: PEPSolve(pep);
115: PEPGetDimensions(pep,&nev,NULL,NULL);
116: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
118: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119: Display solution and clean up
120: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
123: PEPDestroy(&pep);
124: MatDestroy(&M);
125: MatDestroy(&C);
126: MatDestroy(&K);
127: SlepcFinalize();
128: return ierr;
129: }
131: /*TEST
133: testset:
134: requires: double
135: suffix: 1
137: TEST*/