| 1376 |
slepc |
1 |
/*
|
|
|
2 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
| 1672 |
slepc |
3 |
SLEPc - Scalable Library for Eigenvalue Problem Computations
|
| 2116 |
eromero |
4 |
Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
|
| 1242 |
slepc |
5 |
|
| 1672 |
slepc |
6 |
This file is part of SLEPc.
|
|
|
7 |
|
|
|
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.
|
|
|
11 |
|
|
|
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.
|
|
|
16 |
|
|
|
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/>.
|
| 1376 |
slepc |
19 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
20 |
*/
|
|
|
21 |
|
| 1242 |
slepc |
22 |
static char help[] = "Generalized Symmetric eigenproblem.\n\n"
|
|
|
23 |
"The problem is Ax = lambda Bx, with:\n"
|
|
|
24 |
" A = Laplacian operator in 2-D\n"
|
|
|
25 |
" B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
|
|
|
26 |
"The command line options are:\n"
|
|
|
27 |
" -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
|
|
|
28 |
" -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
|
|
|
29 |
" -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";
|
|
|
30 |
|
| 2283 |
jroman |
31 |
#include <slepceps.h>
|
| 1242 |
slepc |
32 |
|
|
|
33 |
#undef __FUNCT__
|
|
|
34 |
#define __FUNCT__ "main"
|
|
|
35 |
int main( int argc, char **argv )
|
|
|
36 |
{
|
|
|
37 |
Mat A, B; /* matrices */
|
|
|
38 |
EPS eps; /* eigenproblem solver context */
|
| 1827 |
jroman |
39 |
ST st; /* spectral transformation context */
|
| 1507 |
slepc |
40 |
const EPSType type;
|
| 1242 |
slepc |
41 |
PetscReal error, tol, re, im;
|
|
|
42 |
PetscScalar kr, ki;
|
|
|
43 |
PetscErrorCode ierr;
|
| 1961 |
antodo |
44 |
PetscInt N, n=10, m, Istart, Iend, II, nev, maxit, i, j, its, nconv, nulldim=0;
|
| 2216 |
jroman |
45 |
PetscBool flag;
|
| 1242 |
slepc |
46 |
|
|
|
47 |
SlepcInitialize(&argc,&argv,(char*)0,help);
|
|
|
48 |
|
|
|
49 |
ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
|
|
|
50 |
ierr = PetscOptionsGetInt(PETSC_NULL,"-m",&m,&flag);CHKERRQ(ierr);
|
| 2177 |
jroman |
51 |
if(!flag) m=n;
|
| 1242 |
slepc |
52 |
N = n*m;
|
|
|
53 |
ierr = PetscOptionsGetInt(PETSC_NULL,"-nulldim",&nulldim,PETSC_NULL);CHKERRQ(ierr);
|
|
|
54 |
ierr = PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%d (%dx%d grid), null(B)=%d\n\n",N,n,m,nulldim);CHKERRQ(ierr);
|
|
|
55 |
|
|
|
56 |
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
57 |
Compute the matrices that define the eigensystem, Ax=kBx
|
|
|
58 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
|
|
|
59 |
|
|
|
60 |
ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
|
|
|
61 |
ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
|
|
|
62 |
ierr = MatSetFromOptions(A);CHKERRQ(ierr);
|
|
|
63 |
|
|
|
64 |
ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
|
|
|
65 |
ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
|
|
|
66 |
ierr = MatSetFromOptions(B);CHKERRQ(ierr);
|
|
|
67 |
|
|
|
68 |
ierr = MatGetOwnershipRange(A,&Istart,&Iend);CHKERRQ(ierr);
|
|
|
69 |
for( II=Istart; II<Iend; II++ ) {
|
| 1961 |
antodo |
70 |
i = II/n; j = II-i*n;
|
|
|
71 |
if(i>0) { ierr = MatSetValue(A,II,II-n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
|
|
|
72 |
if(i<m-1) { ierr = MatSetValue(A,II,II+n,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
|
|
|
73 |
if(j>0) { ierr = MatSetValue(A,II,II-1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
|
|
|
74 |
if(j<n-1) { ierr = MatSetValue(A,II,II+1,-1.0,INSERT_VALUES);CHKERRQ(ierr); }
|
|
|
75 |
ierr = MatSetValue(A,II,II,4.0,INSERT_VALUES);CHKERRQ(ierr);
|
|
|
76 |
if (II>=nulldim) { ierr = MatSetValue(B,II,II,4.0,INSERT_VALUES);CHKERRQ(ierr); }
|
| 1242 |
slepc |
77 |
}
|
|
|
78 |
|
|
|
79 |
ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
|
|
|
80 |
ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
|
|
|
81 |
ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
|
|
|
82 |
ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
|
|
|
83 |
|
|
|
84 |
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
85 |
Create the eigensolver and set various options
|
|
|
86 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
|
|
|
87 |
|
|
|
88 |
/*
|
|
|
89 |
Create eigensolver context
|
|
|
90 |
*/
|
|
|
91 |
ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
|
|
|
92 |
|
|
|
93 |
/*
|
|
|
94 |
Set operators. In this case, it is a generalized eigenvalue problem
|
|
|
95 |
*/
|
|
|
96 |
ierr = EPSSetOperators(eps,A,B);CHKERRQ(ierr);
|
|
|
97 |
ierr = EPSSetProblemType(eps,EPS_GHEP);CHKERRQ(ierr);
|
|
|
98 |
|
| 1827 |
jroman |
99 |
/*
|
|
|
100 |
Use shift-and-invert to avoid solving linear systems with a singular B
|
|
|
101 |
in case nulldim>0
|
|
|
102 |
*/
|
|
|
103 |
ierr = EPSGetST(eps,&st);CHKERRQ(ierr);
|
| 2122 |
jroman |
104 |
ierr = EPSSetTarget(eps,0.0);CHKERRQ(ierr);
|
|
|
105 |
ierr = EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);CHKERRQ(ierr);
|
| 2092 |
jroman |
106 |
ierr = STSetType(st,STSINVERT);CHKERRQ(ierr);
|
| 1827 |
jroman |
107 |
|
| 1242 |
slepc |
108 |
/*
|
|
|
109 |
Set solver parameters at runtime
|
|
|
110 |
*/
|
|
|
111 |
ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
|
|
|
112 |
|
|
|
113 |
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
114 |
Solve the eigensystem
|
|
|
115 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
|
|
|
116 |
|
|
|
117 |
ierr = EPSSolve(eps);CHKERRQ(ierr);
|
|
|
118 |
ierr = EPSGetIterationNumber(eps, &its);CHKERRQ(ierr);
|
|
|
119 |
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);CHKERRQ(ierr);
|
|
|
120 |
|
|
|
121 |
/*
|
|
|
122 |
Optional: Get some information from the solver and display it
|
|
|
123 |
*/
|
|
|
124 |
ierr = EPSGetType(eps,&type);CHKERRQ(ierr);
|
|
|
125 |
ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
|
| 1575 |
slepc |
126 |
ierr = EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
|
| 1242 |
slepc |
127 |
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);CHKERRQ(ierr);
|
|
|
128 |
ierr = EPSGetTolerances(eps,&tol,&maxit);CHKERRQ(ierr);
|
|
|
129 |
ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);CHKERRQ(ierr);
|
|
|
130 |
|
|
|
131 |
/* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
132 |
Display solution and clean up
|
|
|
133 |
- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
|
|
|
134 |
|
|
|
135 |
/*
|
|
|
136 |
Get number of converged approximate eigenpairs
|
|
|
137 |
*/
|
|
|
138 |
ierr = EPSGetConverged(eps,&nconv);CHKERRQ(ierr);
|
|
|
139 |
ierr = PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
|
|
|
140 |
CHKERRQ(ierr);
|
|
|
141 |
|
|
|
142 |
if (nconv>0) {
|
|
|
143 |
/*
|
|
|
144 |
Display eigenvalues and relative errors
|
|
|
145 |
*/
|
|
|
146 |
ierr = PetscPrintf(PETSC_COMM_WORLD,
|
|
|
147 |
" k ||Ax-kBx||/||kx||\n"
|
|
|
148 |
" ----------------- ------------------\n" );CHKERRQ(ierr);
|
|
|
149 |
|
|
|
150 |
for( i=0; i<nconv; i++ ) {
|
|
|
151 |
/*
|
|
|
152 |
Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
|
|
|
153 |
ki (imaginary part)
|
|
|
154 |
*/
|
|
|
155 |
ierr = EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
|
|
|
156 |
/*
|
|
|
157 |
Compute the relative error associated to each eigenpair
|
|
|
158 |
*/
|
|
|
159 |
ierr = EPSComputeRelativeError(eps,i,&error);CHKERRQ(ierr);
|
|
|
160 |
|
|
|
161 |
#ifdef PETSC_USE_COMPLEX
|
|
|
162 |
re = PetscRealPart(kr);
|
|
|
163 |
im = PetscImaginaryPart(kr);
|
|
|
164 |
#else
|
|
|
165 |
re = kr;
|
|
|
166 |
im = ki;
|
|
|
167 |
#endif
|
|
|
168 |
if (im!=0.0) {
|
| 1892 |
jroman |
169 |
ierr = PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);CHKERRQ(ierr);
|
| 1242 |
slepc |
170 |
} else {
|
| 1892 |
jroman |
171 |
ierr = PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);CHKERRQ(ierr);
|
| 1242 |
slepc |
172 |
}
|
|
|
173 |
}
|
|
|
174 |
ierr = PetscPrintf(PETSC_COMM_WORLD,"\n" );CHKERRQ(ierr);
|
|
|
175 |
}
|
|
|
176 |
|
|
|
177 |
/*
|
|
|
178 |
Free work space
|
|
|
179 |
*/
|
|
|
180 |
ierr = EPSDestroy(eps);CHKERRQ(ierr);
|
| 2305 |
jroman |
181 |
ierr = MatDestroy(&A);CHKERRQ(ierr);
|
|
|
182 |
ierr = MatDestroy(&B);CHKERRQ(ierr);
|
| 1242 |
slepc |
183 |
ierr = SlepcFinalize();CHKERRQ(ierr);
|
|
|
184 |
return 0;
|
|
|
185 |
}
|
|
|
186 |
|