Actual source code: shift.c
1: /*
2: Shift spectral transformation, applies (A + sigma I) as operator, or
3: inv(B)(A + sigma B) for generalized problems
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <slepc-private/stimpl.h> /*I "slepcst.h" I*/
29: PetscErrorCode STApply_Shift(ST st,Vec x,Vec y)
30: {
34: if (st->nmat>1) {
35: /* generalized eigenproblem: y = B^-1 (A + sB) x */
36: MatMult(st->T[0],x,st->w);
37: STMatSolve(st,1,st->w,y);
38: } else {
39: /* standard eigenproblem: y = (A + sI) x */
40: MatMult(st->T[0],x,y);
41: }
42: return(0);
43: }
47: PetscErrorCode STApplyTranspose_Shift(ST st,Vec x,Vec y)
48: {
52: if (st->nmat>1) {
53: /* generalized eigenproblem: y = (A + sB)^T B^-T x */
54: STMatSolveTranspose(st,1,x,st->w);
55: MatMultTranspose(st->T[0],st->w,y);
56: } else {
57: /* standard eigenproblem: y = (A^T + sI) x */
58: MatMultTranspose(st->T[0],x,y);
59: }
60: return(0);
61: }
65: PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
66: {
67: PetscInt j;
70: for (j=0;j<n;j++) {
71: eigr[j] -= st->sigma;
72: }
73: return(0);
74: }
78: PetscErrorCode STPostSolve_Shift(ST st)
79: {
81: PetscScalar s;
84: if (st->shift_matrix == ST_MATMODE_INPLACE) {
85: if (st->nmat>1) {
86: if (st->nmat==3) {
87: MatAXPY(st->A[0],-st->sigma*st->sigma,st->A[2],st->str);
88: MatAXPY(st->A[1],2.0*st->sigma,st->A[2],st->str);
89: s = st->sigma;
90: } else s = -st->sigma;
91: MatAXPY(st->A[0],s,st->A[1],st->str);
92: } else {
93: MatShift(st->A[0],st->sigma);
94: }
95: st->Astate[0] = ((PetscObject)st->A[0])->state;
96: st->setupcalled = 0;
97: }
98: return(0);
99: }
103: PetscErrorCode STSetUp_Shift(ST st)
104: {
108: if (st->nmat<3) {
109: /* T[1] = B */
110: if (st->nmat>1) { PetscObjectReference((PetscObject)st->A[1]); }
111: st->T[1] = st->A[1];
112: /* T[0] = A+sigma*B */
113: STMatGAXPY_Private(st,st->sigma,0.0,1,0,PETSC_TRUE);
114: } else {
115: /* T[2] = C */
116: PetscObjectReference((PetscObject)st->A[2]);
117: st->T[2] = st->A[2];
118: /* T[0] = A-sigma*B+sigma*sigma*C */
119: STMatGAXPY_Private(st,-st->sigma,0.0,2,0,PETSC_TRUE);
120: /* T[1] = B-2*sigma*C */
121: STMatGAXPY_Private(st,-2.0*st->sigma,0.0,1,1,PETSC_TRUE);
122: }
123: if (st->nmat==2) {
124: if (!st->ksp) { STGetKSP(st,&st->ksp); }
125: KSPSetOperators(st->ksp,st->T[1],st->T[1],DIFFERENT_NONZERO_PATTERN);
126: KSPSetUp(st->ksp);
127: st->kspidx = 1;
128: }
129: return(0);
130: }
134: PetscErrorCode STSetShift_Shift(ST st,PetscScalar newshift)
135: {
137: MatStructure flg;
140: /* Nothing to be done if STSetUp has not been called yet */
141: if (!st->setupcalled) return(0);
143: if (st->nmat<3) {
144: STMatGAXPY_Private(st,newshift,st->sigma,1,0,PETSC_FALSE);
145: } else {
146: STMatGAXPY_Private(st,-newshift,-st->sigma,2,2,PETSC_FALSE);
147: STMatGAXPY_Private(st,-2.0*newshift,-2.0*st->sigma,1,1,PETSC_FALSE);
148: }
150: if (st->kspidx==0 || (st->nmat==3 && st->kspidx==1)) { /* Update KSP operator */
151: /* Check if the new KSP matrix has the same zero structure */
152: if (st->nmat>1 && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) flg = DIFFERENT_NONZERO_PATTERN;
153: else flg = SAME_NONZERO_PATTERN;
154: KSPSetOperators(st->ksp,st->T[st->kspidx],st->T[st->kspidx],flg);
155: KSPSetUp(st->ksp);
156: }
157: return(0);
158: }
162: PetscErrorCode STSetFromOptions_Shift(ST st)
163: {
165: PC pc;
166: PCType pctype;
167: KSPType ksptype;
170: if (!st->ksp) { STGetKSP(st,&st->ksp); }
171: KSPGetPC(st->ksp,&pc);
172: KSPGetType(st->ksp,&ksptype);
173: PCGetType(pc,&pctype);
174: if (!pctype && !ksptype) {
175: if (st->shift_matrix == ST_MATMODE_SHELL) {
176: /* in shell mode use GMRES with Jacobi as the default */
177: KSPSetType(st->ksp,KSPGMRES);
178: PCSetType(pc,PCJACOBI);
179: } else {
180: /* use direct solver as default */
181: KSPSetType(st->ksp,KSPPREONLY);
182: PCSetType(pc,PCREDUNDANT);
183: }
184: }
185: return(0);
186: }
190: PETSC_EXTERN PetscErrorCode STCreate_Shift(ST st)
191: {
193: st->ops->apply = STApply_Shift;
194: st->ops->getbilinearform = STGetBilinearForm_Default;
195: st->ops->applytrans = STApplyTranspose_Shift;
196: st->ops->postsolve = STPostSolve_Shift;
197: st->ops->backtransform = STBackTransform_Shift;
198: st->ops->setfromoptions = STSetFromOptions_Shift;
199: st->ops->setup = STSetUp_Shift;
200: st->ops->setshift = STSetShift_Shift;
201: return(0);
202: }