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) { MatAXPY(st->mat,-st->sigma,st->B,st->str); }
142:       else { MatShift(st->mat,-st->sigma); }
143:       KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
144:     } else {
145:       st->mat = PETSC_NULL;
146:       KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
147:     }
148:   }
149:   KSPSetUp(st->ksp);
150:   return(0);
151: }

155: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
156: {
158:   MatStructure   flg;

161:   /* Nothing to be done if STSetUp has not been called yet */
162:   if (!st->setupcalled) return(0);

164:   /* Check if the new KSP matrix has the same zero structure */
165:   if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
166:     flg = DIFFERENT_NONZERO_PATTERN;
167:   } else {
168:     flg = SAME_NONZERO_PATTERN;
169:   }

171:   switch (st->shift_matrix) {
172:   case ST_MATMODE_INPLACE:
173:     /* Undo previous operations */
174:     if (st->sigma != 0.0) {
175:       if (st->B) {
176:         MatAXPY(st->A,st->sigma,st->B,st->str);
177:       } else {
178:         MatShift(st->A,st->sigma);
179:       }
180:     }
181:     /* Apply new shift */
182:     if (newshift != 0.0) {
183:       if (st->B) {
184:         MatAXPY(st->A,-newshift,st->B,st->str);
185:       } else {
186:         MatShift(st->A,-newshift);
187:       }
188:     }
189:     KSPSetOperators(st->ksp,st->A,st->A,flg);
190:     break;
191:   case ST_MATMODE_SHELL:
192:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
193:     break;
194:   default:
195:     if (st->mat) {
196:       MatCopy(st->A,st->mat,DIFFERENT_NONZERO_PATTERN);
197:     } else {
198:       MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
199:     }
200:     if (newshift != 0.0) {
201:       if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
202:       else { MatShift(st->mat,-newshift); }
203:     }
204:     KSPSetOperators(st->ksp,st->mat,st->mat,flg);
205:   }
206:   st->sigma = newshift;
207:   KSPSetUp(st->ksp);
208:   return(0);
209: }

213: PetscErrorCode STSetFromOptions_Sinvert(ST st)
214: {
216:   PC             pc;
217:   const PCType   pctype;
218:   const KSPType  ksptype;

221:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
222:   KSPGetPC(st->ksp,&pc);
223:   KSPGetType(st->ksp,&ksptype);
224:   PCGetType(pc,&pctype);
225:   if (!pctype && !ksptype) {
226:     if (st->shift_matrix == ST_MATMODE_SHELL) {
227:       /* in shell mode use GMRES with Jacobi as the default */
228:       KSPSetType(st->ksp,KSPGMRES);
229:       PCSetType(pc,PCJACOBI);
230:     } else {
231:       /* use direct solver as default */
232:       KSPSetType(st->ksp,KSPPREONLY);
233:       PCSetType(pc,PCREDUNDANT);
234:     }
235:   }
236:   return(0);
237: }

239: EXTERN_C_BEGIN
242: PetscErrorCode STCreate_Sinvert(ST st)
243: {
245:   st->ops->apply           = STApply_Sinvert;
246:   st->ops->getbilinearform = STGetBilinearForm_Default;
247:   st->ops->applytrans      = STApplyTranspose_Sinvert;
248:   st->ops->postsolve       = STPostSolve_Sinvert;
249:   st->ops->backtr          = STBackTransform_Sinvert;
250:   st->ops->setup           = STSetUp_Sinvert;
251:   st->ops->setshift        = STSetShift_Sinvert;
252:   st->ops->setfromoptions  = STSetFromOptions_Sinvert;
253:   st->ops->checknullspace  = STCheckNullSpace_Default;
254:   return(0);
255: }
256: EXTERN_C_END