Actual source code: test10.c
slepc-3.17.1 2022-04-11
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 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;
25: *errest = res/norm;
26: PetscFunctionReturn(0);
27: }
29: int main(int argc,char **argv)
30: {
31: Mat M,C,K,A[3]; /* problem matrices */
32: PEP pep; /* polynomial eigenproblem solver context */
33: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
34: PetscBool flag;
35: PetscReal norm;
37: SlepcInitialize(&argc,&argv,(char*)0,help);
39: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
41: if (!flag) m=n;
42: N = n*m;
43: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
45: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
47: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49: /* K is the 2-D Laplacian */
50: MatCreate(PETSC_COMM_WORLD,&K);
51: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
52: MatSetFromOptions(K);
53: MatSetUp(K);
54: MatGetOwnershipRange(K,&Istart,&Iend);
55: for (II=Istart;II<Iend;II++) {
56: i = II/n; j = II-i*n;
57: if (i>0) MatSetValue(K,II,II-n,-1.0,INSERT_VALUES);
58: if (i<m-1) MatSetValue(K,II,II+n,-1.0,INSERT_VALUES);
59: if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
60: if (j<n-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
61: MatSetValue(K,II,II,4.0,INSERT_VALUES);
62: }
63: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
64: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
66: /* C is the 1-D Laplacian on horizontal lines */
67: MatCreate(PETSC_COMM_WORLD,&C);
68: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
69: MatSetFromOptions(C);
70: MatSetUp(C);
71: MatGetOwnershipRange(C,&Istart,&Iend);
72: for (II=Istart;II<Iend;II++) {
73: i = II/n; j = II-i*n;
74: if (j>0) MatSetValue(C,II,II-1,-1.0,INSERT_VALUES);
75: if (j<n-1) MatSetValue(C,II,II+1,-1.0,INSERT_VALUES);
76: MatSetValue(C,II,II,2.0,INSERT_VALUES);
77: }
78: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
79: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
81: /* M is a diagonal matrix */
82: MatCreate(PETSC_COMM_WORLD,&M);
83: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
84: MatSetFromOptions(M);
85: MatSetUp(M);
86: MatGetOwnershipRange(M,&Istart,&Iend);
87: for (II=Istart;II<Iend;II++) MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
88: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
89: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Create the eigensolver and set various options
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: PEPCreate(PETSC_COMM_WORLD,&pep);
96: A[0] = K; A[1] = C; A[2] = M;
97: PEPSetOperators(pep,3,A);
98: PEPSetProblemType(pep,PEP_HERMITIAN);
99: PEPSetDimensions(pep,4,20,PETSC_DEFAULT);
101: /* setup convergence test relative to the norm of M */
102: MatNorm(M,NORM_1,&norm);
103: PEPSetConvergenceTestFunction(pep,MyConvergedRel,&norm,NULL);
104: PEPSetFromOptions(pep);
106: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: Solve the eigensystem
108: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110: PEPSolve(pep);
111: PEPGetDimensions(pep,&nev,NULL,NULL);
112: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Display solution and clean up
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
119: PEPDestroy(&pep);
120: MatDestroy(&M);
121: MatDestroy(&C);
122: MatDestroy(&K);
123: SlepcFinalize();
124: return 0;
125: }
127: /*TEST
129: testset:
130: requires: double
131: suffix: 1
133: TEST*/