Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
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";
15 :
16 : #include <slepcpep.h>
17 :
18 : /*
19 : MyConvergedRel - Convergence test relative to the norm of M (given in ctx).
20 : */
21 6 : PetscErrorCode MyConvergedRel(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
22 : {
23 6 : PetscReal norm = *(PetscReal*)ctx;
24 :
25 6 : PetscFunctionBegin;
26 6 : *errest = res/norm;
27 6 : PetscFunctionReturn(PETSC_SUCCESS);
28 : }
29 :
30 1 : int main(int argc,char **argv)
31 : {
32 1 : Mat M,C,K,A[3]; /* problem matrices */
33 1 : PEP pep; /* polynomial eigenproblem solver context */
34 1 : PetscInt N,n=10,m,Istart,Iend,II,nev,i,j;
35 1 : PetscBool flag;
36 1 : PetscReal norm;
37 :
38 1 : PetscFunctionBeginUser;
39 1 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
40 :
41 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
42 1 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
43 1 : if (!flag) m=n;
44 1 : N = n*m;
45 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
46 :
47 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48 : Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
49 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
50 :
51 : /* K is the 2-D Laplacian */
52 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
53 1 : PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N));
54 1 : PetscCall(MatSetFromOptions(K));
55 1 : PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
56 101 : for (II=Istart;II<Iend;II++) {
57 100 : i = II/n; j = II-i*n;
58 100 : if (i>0) PetscCall(MatSetValue(K,II,II-n,-1.0,INSERT_VALUES));
59 100 : if (i<m-1) PetscCall(MatSetValue(K,II,II+n,-1.0,INSERT_VALUES));
60 100 : if (j>0) PetscCall(MatSetValue(K,II,II-1,-1.0,INSERT_VALUES));
61 100 : if (j<n-1) PetscCall(MatSetValue(K,II,II+1,-1.0,INSERT_VALUES));
62 100 : PetscCall(MatSetValue(K,II,II,4.0,INSERT_VALUES));
63 : }
64 1 : PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
65 1 : PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
66 :
67 : /* C is the 1-D Laplacian on horizontal lines */
68 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
69 1 : PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N));
70 1 : PetscCall(MatSetFromOptions(C));
71 1 : PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
72 101 : for (II=Istart;II<Iend;II++) {
73 100 : i = II/n; j = II-i*n;
74 100 : if (j>0) PetscCall(MatSetValue(C,II,II-1,-1.0,INSERT_VALUES));
75 100 : if (j<n-1) PetscCall(MatSetValue(C,II,II+1,-1.0,INSERT_VALUES));
76 100 : PetscCall(MatSetValue(C,II,II,2.0,INSERT_VALUES));
77 : }
78 1 : PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
79 1 : PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
80 :
81 : /* M is a diagonal matrix */
82 1 : PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
83 1 : PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N));
84 1 : PetscCall(MatSetFromOptions(M));
85 1 : PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
86 101 : for (II=Istart;II<Iend;II++) PetscCall(MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES));
87 1 : PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
88 1 : PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
89 :
90 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91 : Create the eigensolver and set various options
92 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
93 :
94 1 : PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
95 1 : A[0] = K; A[1] = C; A[2] = M;
96 1 : PetscCall(PEPSetOperators(pep,3,A));
97 1 : PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
98 1 : PetscCall(PEPSetDimensions(pep,4,20,PETSC_DETERMINE));
99 :
100 : /* setup convergence test relative to the norm of M */
101 1 : PetscCall(MatNorm(M,NORM_1,&norm));
102 1 : PetscCall(PEPSetConvergenceTestFunction(pep,MyConvergedRel,&norm,NULL));
103 1 : PetscCall(PEPSetFromOptions(pep));
104 :
105 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106 : Solve the eigensystem
107 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108 :
109 1 : PetscCall(PEPSolve(pep));
110 1 : PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
111 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
112 :
113 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114 : Display solution and clean up
115 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116 :
117 1 : PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
118 1 : PetscCall(PEPDestroy(&pep));
119 1 : PetscCall(MatDestroy(&M));
120 1 : PetscCall(MatDestroy(&C));
121 1 : PetscCall(MatDestroy(&K));
122 1 : PetscCall(SlepcFinalize());
123 : return 0;
124 : }
125 :
126 : /*TEST
127 :
128 : testset:
129 : requires: double
130 : suffix: 1
131 :
132 : TEST*/
|