Actual source code: test10.c
slepc-3.18.0 2022-10-01
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: return 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;
38: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
41: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
42: if (!flag) m=n;
43: N = n*m;
44: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m);
46: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
48: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50: /* K is the 2-D Laplacian */
51: MatCreate(PETSC_COMM_WORLD,&K);
52: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
53: MatSetFromOptions(K);
54: MatSetUp(K);
55: MatGetOwnershipRange(K,&Istart,&Iend);
56: for (II=Istart;II<Iend;II++) {
57: i = II/n; j = II-i*n;
58: if (i>0) MatSetValue(K,II,II-n,-1.0,INSERT_VALUES);
59: if (i<m-1) MatSetValue(K,II,II+n,-1.0,INSERT_VALUES);
60: if (j>0) MatSetValue(K,II,II-1,-1.0,INSERT_VALUES);
61: if (j<n-1) MatSetValue(K,II,II+1,-1.0,INSERT_VALUES);
62: MatSetValue(K,II,II,4.0,INSERT_VALUES);
63: }
64: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
67: /* C is the 1-D Laplacian on horizontal lines */
68: MatCreate(PETSC_COMM_WORLD,&C);
69: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
70: MatSetFromOptions(C);
71: MatSetUp(C);
72: MatGetOwnershipRange(C,&Istart,&Iend);
73: for (II=Istart;II<Iend;II++) {
74: i = II/n; j = II-i*n;
75: if (j>0) MatSetValue(C,II,II-1,-1.0,INSERT_VALUES);
76: if (j<n-1) MatSetValue(C,II,II+1,-1.0,INSERT_VALUES);
77: MatSetValue(C,II,II,2.0,INSERT_VALUES);
78: }
79: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
80: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
82: /* M is a diagonal matrix */
83: MatCreate(PETSC_COMM_WORLD,&M);
84: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
85: MatSetFromOptions(M);
86: MatSetUp(M);
87: MatGetOwnershipRange(M,&Istart,&Iend);
88: for (II=Istart;II<Iend;II++) MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
89: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Create the eigensolver and set various options
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: PEPCreate(PETSC_COMM_WORLD,&pep);
97: A[0] = K; A[1] = C; A[2] = M;
98: PEPSetOperators(pep,3,A);
99: PEPSetProblemType(pep,PEP_HERMITIAN);
100: PEPSetDimensions(pep,4,20,PETSC_DEFAULT);
102: /* setup convergence test relative to the norm of M */
103: MatNorm(M,NORM_1,&norm);
104: PEPSetConvergenceTestFunction(pep,MyConvergedRel,&norm,NULL);
105: PEPSetFromOptions(pep);
107: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108: Solve the eigensystem
109: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111: PEPSolve(pep);
112: PEPGetDimensions(pep,&nev,NULL,NULL);
113: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Display solution and clean up
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
120: PEPDestroy(&pep);
121: MatDestroy(&M);
122: MatDestroy(&C);
123: MatDestroy(&K);
124: SlepcFinalize();
125: return 0;
126: }
128: /*TEST
130: testset:
131: requires: double
132: suffix: 1
134: TEST*/