Actual source code: ex16.c
slepc-3.7.4 2017-05-17
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Simple quadratic eigenvalue problem.\n\n"
23: "The command line options are:\n"
24: " -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
25: " -m <m>, where <m> = number of grid subdivisions in y dimension.\n\n";
27: #include <slepcpep.h>
31: int main(int argc,char **argv)
32: {
33: Mat M,C,K,A[3]; /* problem matrices */
34: PEP pep; /* polynomial eigenproblem solver context */
35: PEPType type;
36: PetscInt N,n=10,m,Istart,Iend,II,nev,i,j,nconv;
37: PetscBool flag,terse;
38: PetscReal error,re,im;
39: PetscScalar kr,ki;
40: Vec xr,xi;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
45: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
46: PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
47: if (!flag) m=n;
48: N = n*m;
49: PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)\n\n",N,n,m);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /* K is the 2-D Laplacian */
56: MatCreate(PETSC_COMM_WORLD,&K);
57: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
58: MatSetFromOptions(K);
59: MatSetUp(K);
60: MatGetOwnershipRange(K,&Istart,&Iend);
61: for (II=Istart;II<Iend;II++) {
62: i = II/n; j = II-i*n;
63: if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
64: if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
65: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
66: if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
67: MatSetValue(K,II,II,4.0,INSERT_VALUES);
68: }
69: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
70: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
72: /* C is the 1-D Laplacian on horizontal lines */
73: MatCreate(PETSC_COMM_WORLD,&C);
74: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
75: MatSetFromOptions(C);
76: MatSetUp(C);
77: MatGetOwnershipRange(C,&Istart,&Iend);
78: for (II=Istart;II<Iend;II++) {
79: i = II/n; j = II-i*n;
80: if (j>0) { MatSetValue(C,II,II-1,-1.0,INSERT_VALUES); }
81: if (j<n-1) { MatSetValue(C,II,II+1,-1.0,INSERT_VALUES); }
82: MatSetValue(C,II,II,2.0,INSERT_VALUES);
83: }
84: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
87: /* M is a diagonal matrix */
88: MatCreate(PETSC_COMM_WORLD,&M);
89: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
90: MatSetFromOptions(M);
91: MatSetUp(M);
92: MatGetOwnershipRange(M,&Istart,&Iend);
93: for (II=Istart;II<Iend;II++) {
94: MatSetValue(M,II,II,(PetscReal)(II+1),INSERT_VALUES);
95: }
96: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create the eigensolver and set various options
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: /*
104: Create eigensolver context
105: */
106: PEPCreate(PETSC_COMM_WORLD,&pep);
108: /*
109: Set matrices and problem type
110: */
111: A[0] = K; A[1] = C; A[2] = M;
112: PEPSetOperators(pep,3,A);
113: PEPSetProblemType(pep,PEP_HERMITIAN);
115: /*
116: Set solver parameters at runtime
117: */
118: PEPSetFromOptions(pep);
120: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121: Solve the eigensystem
122: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: PEPSolve(pep);
126: /*
127: Optional: Get some information from the solver and display it
128: */
129: PEPGetType(pep,&type);
130: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
131: PEPGetDimensions(pep,&nev,NULL,NULL);
132: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Display solution and clean up
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: /* show detailed info unless -terse option is given by user */
139: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
140: if (terse) {
141: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
142: } else {
143: PEPGetConverged(pep,&nconv);
144: if (nconv>0) {
145: MatCreateVecs(M,&xr,&xi);
146: /* display eigenvalues and relative errors */
147: PetscPrintf(PETSC_COMM_WORLD,
148: "\n k ||P(k)x||/||kx||\n"
149: " ----------------- ------------------\n");
150: for (i=0;i<nconv;i++) {
151: /* get converged eigenpairs */
152: PEPGetEigenpair(pep,i,&kr,&ki,xr,xi);
153: /* compute the relative error associated to each eigenpair */
154: PEPComputeError(pep,i,PEP_ERROR_RELATIVE,&error);
155: #if defined(PETSC_USE_COMPLEX)
156: re = PetscRealPart(kr);
157: im = PetscImaginaryPart(kr);
158: #else
159: re = kr;
160: im = ki;
161: #endif
162: if (im!=0.0) {
163: PetscPrintf(PETSC_COMM_WORLD," %9f%+9fi %12g\n",(double)re,(double)im,(double)error);
164: } else {
165: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",(double)re,(double)error);
166: }
167: }
168: PetscPrintf(PETSC_COMM_WORLD,"\n");
169: VecDestroy(&xr);
170: VecDestroy(&xi);
171: }
172: }
173: PEPDestroy(&pep);
174: MatDestroy(&M);
175: MatDestroy(&C);
176: MatDestroy(&K);
177: SlepcFinalize();
178: return ierr;
179: }