Actual source code: sinvert.c
1: /*
2: Implements the shift-and-invert technique for eigenvalue problems.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include private/stimpl.h
28: PetscErrorCode STApply_Sinvert(ST st,Vec x,Vec y)
29: {
33: if (st->B) {
34: /* generalized eigenproblem: y = (A - sB)^-1 B x */
35: MatMult(st->B,x,st->w);
36: STAssociatedKSPSolve(st,st->w,y);
37: }
38: else {
39: /* standard eigenproblem: y = (A - sI)^-1 x */
40: STAssociatedKSPSolve(st,x,y);
41: }
42: return(0);
43: }
47: PetscErrorCode STApplyTranspose_Sinvert(ST st,Vec x,Vec y)
48: {
52: if (st->B) {
53: /* generalized eigenproblem: y = B^T (A - sB)^-T x */
54: STAssociatedKSPSolveTranspose(st,x,st->w);
55: MatMultTranspose(st->B,st->w,y);
56: }
57: else {
58: /* standard eigenproblem: y = (A - sI)^-T x */
59: STAssociatedKSPSolveTranspose(st,x,y);
60: }
61: return(0);
62: }
66: PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
67: {
68: PetscInt j;
69: #ifndef PETSC_USE_COMPLEX
70: PetscScalar t;
74: for (j=0;j<n;j++) {
75: if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
76: else {
77: t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
78: eigr[j] = eigr[j] / t + st->sigma;
79: eigi[j] = - eigi[j] / t;
80: }
81: }
82: #else
85: for (j=0;j<n;j++) {
86: eigr[j] = 1.0 / eigr[j] + st->sigma;
87: }
88: #endif
89: return(0);
90: }
94: PetscErrorCode STPostSolve_Sinvert(ST st)
95: {
99: if (st->shift_matrix == ST_MATMODE_INPLACE) {
100: if( st->B ) {
101: MatAXPY(st->A,st->sigma,st->B,st->str);
102: } else {
103: MatShift(st->A,st->sigma);
104: }
105: st->setupcalled = 0;
106: }
107: return(0);
108: }
112: PetscErrorCode STSetUp_Sinvert(ST st)
113: {
117: if (st->mat) { MatDestroy(st->mat); }
119: /* if the user did not set the shift, use the target value */
120: if (!st->sigma_set) st->sigma = st->defsigma;
122: switch (st->shift_matrix) {
123: case ST_MATMODE_INPLACE:
124: st->mat = PETSC_NULL;
125: if (st->sigma != 0.0) {
126: if (st->B) {
127: MatAXPY(st->A,-st->sigma,st->B,st->str);
128: } else {
129: MatShift(st->A,-st->sigma);
130: }
131: }
132: KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
133: break;
134: case ST_MATMODE_SHELL:
135: STMatShellCreate(st,&st->mat);
136: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
137: break;
138: default:
139: if (st->sigma != 0.0) {
140: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
141: if (st->B) {
142: MatAXPY(st->mat,-st->sigma,st->B,st->str);
143: } else {
144: MatShift(st->mat,-st->sigma);
145: }
146: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
147: } else {
148: st->mat = PETSC_NULL;
149: KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
150: }
151: }
153: KSPSetUp(st->ksp);
154: return(0);
155: }
159: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
160: {
162: MatStructure flg;
166: /* Nothing to be done if STSetUp has not been called yet */
167: if (!st->setupcalled) return(0);
168:
169: /* Check if the new KSP matrix has the same zero structure */
170: if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
171: flg = DIFFERENT_NONZERO_PATTERN;
172: } else {
173: flg = SAME_NONZERO_PATTERN;
174: }
176: switch (st->shift_matrix) {
177: case ST_MATMODE_INPLACE:
178: /* Undo previous operations */
179: if (st->sigma != 0.0) {
180: if (st->B) {
181: MatAXPY(st->A,st->sigma,st->B,st->str);
182: } else {
183: MatShift(st->A,st->sigma);
184: }
185: }
186: /* Apply new shift */
187: if (newshift != 0.0) {
188: if (st->B) {
189: MatAXPY(st->A,-newshift,st->B,st->str);
190: } else {
191: MatShift(st->A,-newshift);
192: }
193: }
194: KSPSetOperators(st->ksp,st->A,st->A,flg);
195: break;
196: case ST_MATMODE_SHELL:
197: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
198: break;
199: default:
200: if (st->mat) {
201: MatCopy(st->A,st->mat,SUBSET_NONZERO_PATTERN);
202: } else {
203: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
204: }
205: if (newshift != 0.0) {
206: if (st->B) {
207: MatAXPY(st->mat,-newshift,st->B,st->str);
208: } else {
209: MatShift(st->mat,-newshift);
210: }
211: }
212: KSPSetOperators(st->ksp,st->mat,st->mat,flg);
213: }
214: st->sigma = newshift;
215: KSPSetUp(st->ksp);
216: return(0);
217: }
221: PetscErrorCode STSetFromOptions_Sinvert(ST st)
222: {
224: PC pc;
225: const PCType pctype;
226: const KSPType ksptype;
230: KSPGetPC(st->ksp,&pc);
231: KSPGetType(st->ksp,&ksptype);
232: PCGetType(pc,&pctype);
233: if (!pctype && !ksptype) {
234: if (st->shift_matrix == ST_MATMODE_SHELL) {
235: /* in shell mode use GMRES with Jacobi as the default */
236: KSPSetType(st->ksp,KSPGMRES);
237: PCSetType(pc,PCJACOBI);
238: } else {
239: /* use direct solver as default */
240: KSPSetType(st->ksp,KSPPREONLY);
241: PCSetType(pc,PCREDUNDANT);
242: }
243: }
245: return(0);
246: }
251: PetscErrorCode STCreate_Sinvert(ST st)
252: {
254: st->data = 0;
256: st->ops->apply = STApply_Sinvert;
257: st->ops->getbilinearform = STGetBilinearForm_Default;
258: st->ops->applytrans = STApplyTranspose_Sinvert;
259: st->ops->postsolve = STPostSolve_Sinvert;
260: st->ops->backtr = STBackTransform_Sinvert;
261: st->ops->setup = STSetUp_Sinvert;
262: st->ops->setshift = STSetShift_Sinvert;
263: st->ops->view = STView_Default;
264: st->ops->setfromoptions = STSetFromOptions_Sinvert;
265:
266: st->checknullspace = STCheckNullSpace_Default;
268: return(0);
269: }