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-2011, Universitat 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> /*I "slepcst.h" I*/
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: #if !defined(PETSC_USE_COMPLEX)
70: PetscScalar t;
71: #endif
74: #if !defined(PETSC_USE_COMPLEX)
75: for (j=0;j<n;j++) {
76: if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
77: else {
78: t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
79: eigr[j] = eigr[j] / t + st->sigma;
80: eigi[j] = - eigi[j] / t;
81: }
82: }
83: #else
84: for (j=0;j<n;j++) {
85: eigr[j] = 1.0 / eigr[j] + st->sigma;
86: }
87: #endif
88: return(0);
89: }
93: PetscErrorCode STPostSolve_Sinvert(ST st)
94: {
98: if (st->shift_matrix == ST_MATMODE_INPLACE) {
99: if (st->B) {
100: MatAXPY(st->A,st->sigma,st->B,st->str);
101: } else {
102: MatShift(st->A,st->sigma);
103: }
104: st->setupcalled = 0;
105: }
106: return(0);
107: }
111: PetscErrorCode STSetUp_Sinvert(ST st)
112: {
116: MatDestroy(&st->mat);
118: /* if the user did not set the shift, use the target value */
119: if (!st->sigma_set) st->sigma = st->defsigma;
121: if (!st->ksp) { STGetKSP(st,&st->ksp); }
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;
165: /* Nothing to be done if STSetUp has not been called yet */
166: if (!st->setupcalled) return(0);
167:
168: /* Check if the new KSP matrix has the same zero structure */
169: if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
170: flg = DIFFERENT_NONZERO_PATTERN;
171: } else {
172: flg = SAME_NONZERO_PATTERN;
173: }
175: switch (st->shift_matrix) {
176: case ST_MATMODE_INPLACE:
177: /* Undo previous operations */
178: if (st->sigma != 0.0) {
179: if (st->B) {
180: MatAXPY(st->A,st->sigma,st->B,st->str);
181: } else {
182: MatShift(st->A,st->sigma);
183: }
184: }
185: /* Apply new shift */
186: if (newshift != 0.0) {
187: if (st->B) {
188: MatAXPY(st->A,-newshift,st->B,st->str);
189: } else {
190: MatShift(st->A,-newshift);
191: }
192: }
193: KSPSetOperators(st->ksp,st->A,st->A,flg);
194: break;
195: case ST_MATMODE_SHELL:
196: KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
197: break;
198: default:
199: if (st->mat) {
200: MatCopy(st->A,st->mat,SUBSET_NONZERO_PATTERN);
201: } else {
202: MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
203: }
204: if (newshift != 0.0) {
205: if (st->B) {
206: MatAXPY(st->mat,-newshift,st->B,st->str);
207: } else {
208: MatShift(st->mat,-newshift);
209: }
210: }
211: KSPSetOperators(st->ksp,st->mat,st->mat,flg);
212: }
213: st->sigma = newshift;
214: KSPSetUp(st->ksp);
215: return(0);
216: }
220: PetscErrorCode STSetFromOptions_Sinvert(ST st)
221: {
223: PC pc;
224: const PCType pctype;
225: const KSPType ksptype;
228: if (!st->ksp) { STGetKSP(st,&st->ksp); }
229: KSPGetPC(st->ksp,&pc);
230: KSPGetType(st->ksp,&ksptype);
231: PCGetType(pc,&pctype);
232: if (!pctype && !ksptype) {
233: if (st->shift_matrix == ST_MATMODE_SHELL) {
234: /* in shell mode use GMRES with Jacobi as the default */
235: KSPSetType(st->ksp,KSPGMRES);
236: PCSetType(pc,PCJACOBI);
237: } else {
238: /* use direct solver as default */
239: KSPSetType(st->ksp,KSPPREONLY);
240: PCSetType(pc,PCREDUNDANT);
241: }
242: }
243: return(0);
244: }
249: PetscErrorCode STCreate_Sinvert(ST st)
250: {
252: st->ops->apply = STApply_Sinvert;
253: st->ops->getbilinearform = STGetBilinearForm_Default;
254: st->ops->applytrans = STApplyTranspose_Sinvert;
255: st->ops->postsolve = STPostSolve_Sinvert;
256: st->ops->backtr = STBackTransform_Sinvert;
257: st->ops->setup = STSetUp_Sinvert;
258: st->ops->setshift = STSetShift_Sinvert;
259: st->ops->setfromoptions = STSetFromOptions_Sinvert;
260: st->ops->checknullspace = STCheckNullSpace_Default;
261: return(0);
262: }