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-2013, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 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 <slepc-private/stimpl.h>          /*I "slepcst.h" I*/

 28: PetscErrorCode STApply_Sinvert(ST st,Vec x,Vec y)
 29: {

 33:   if (st->nmat>1) {
 34:     /* generalized eigenproblem: y = (A - sB)^-1 B x */
 35:     MatMult(st->T[0],x,st->w);
 36:     STMatSolve(st,1,st->w,y);
 37:   } else {
 38:     /* standard eigenproblem: y = (A - sI)^-1 x */
 39:     STMatSolve(st,1,x,y);
 40:   }
 41:   return(0);
 42: }

 46: PetscErrorCode STApplyTranspose_Sinvert(ST st,Vec x,Vec y)
 47: {

 51:   if (st->nmat>1) {
 52:     /* generalized eigenproblem: y = B^T (A - sB)^-T x */
 53:     STMatSolveTranspose(st,1,x,st->w);
 54:     MatMultTranspose(st->T[0],st->w,y);
 55:   } else {
 56:     /* standard eigenproblem: y = (A - sI)^-T x */
 57:     STMatSolveTranspose(st,1,x,y);
 58:   }
 59:   return(0);
 60: }

 64: PetscErrorCode STBackTransform_Sinvert(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 65: {
 66:   PetscInt    j;
 67: #if !defined(PETSC_USE_COMPLEX)
 68:   PetscScalar t;
 69: #endif

 72: #if !defined(PETSC_USE_COMPLEX)
 73:   for (j=0;j<n;j++) {
 74:     if (eigi[j] == 0) eigr[j] = 1.0 / eigr[j] + st->sigma;
 75:     else {
 76:       t = eigr[j] * eigr[j] + eigi[j] * eigi[j];
 77:       eigr[j] = eigr[j] / t + st->sigma;
 78:       eigi[j] = - eigi[j] / t;
 79:     }
 80:   }
 81: #else
 82:   for (j=0;j<n;j++) {
 83:     eigr[j] = 1.0 / eigr[j] + st->sigma;
 84:   }
 85: #endif
 86:   return(0);
 87: }

 91: PetscErrorCode STPostSolve_Sinvert(ST st)
 92: {
 94:   PetscScalar    s;

 97:   if (st->shift_matrix == ST_MATMODE_INPLACE) {
 98:     if (st->nmat>1) {
 99:       if (st->nmat==3) {
100:         MatAXPY(st->A[0],-st->sigma*st->sigma,st->A[2],st->str);
101:         MatAXPY(st->A[1],-2.0*st->sigma,st->A[2],st->str);
102:         s = -st->sigma;
103:       } else s = st->sigma;
104:       MatAXPY(st->A[0],s,st->A[1],st->str);
105:     } else {
106:       MatShift(st->A[0],st->sigma);
107:     }
108:     st->Astate[0] = ((PetscObject)st->A[0])->state;
109:     st->setupcalled = 0;
110:   }
111:   return(0);
112: }

116: PetscErrorCode STSetUp_Sinvert(ST st)
117: {
119:   PetscScalar    gamma;

122:   /* if the user did not set the shift, use the target value */
123:   if (!st->sigma_set) st->sigma = st->defsigma;
124:   if (st->nmat<3) {
125:     /* T[0] = B */
126:     if (st->nmat>1) { PetscObjectReference((PetscObject)st->A[1]); }
127:     st->T[0] = st->A[1];
128:     gamma = -st->sigma;
129:   } else {
130:     /* T[0] = C */
131:     PetscObjectReference((PetscObject)st->A[2]);
132:     st->T[0] = st->A[2];
133:     /* T[2] = A+sigma*B+sigma*sigma*C */
134:     STMatGAXPY_Private(st,st->sigma,0.0,2,2,PETSC_TRUE);
135:     gamma = 2.0*st->sigma;
136:   }
137:   /* T[1] = A-sigma*B or B+2*sigma*C  */
138:   STMatGAXPY_Private(st,gamma,0.0,1,1,PETSC_TRUE);
139:   if (st->nmat<3) {
140:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
141:     KSPSetOperators(st->ksp,st->T[1],st->T[1],DIFFERENT_NONZERO_PATTERN);
142:     KSPSetUp(st->ksp);
143:     st->kspidx = 1;
144:   }
145:   return(0);
146: }

150: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
151: {
153:   MatStructure   flg;
154:   PetscScalar    alpha,beta;

157:   /* Nothing to be done if STSetUp has not been called yet */
158:   if (!st->setupcalled) return(0);
159:   if (st->nmat<3) {
160:     alpha = -newshift; beta = -st->sigma;
161:   } else {
162:     STMatGAXPY_Private(st,newshift,st->sigma,2,2,PETSC_FALSE);
163:     alpha = 2.0*newshift; beta = 2.0*st->sigma;
164:   }
165:   STMatGAXPY_Private(st,alpha,beta,1,1,PETSC_FALSE);
166:   if (st->kspidx==1 || (st->nmat==3 && st->kspidx==2)) {  /* Update KSP operator */
167:     /* Check if the new KSP matrix has the same zero structure */
168:     if (st->nmat>1 && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) flg = DIFFERENT_NONZERO_PATTERN;
169:     else flg = SAME_NONZERO_PATTERN;
170:     KSPSetOperators(st->ksp,st->T[st->kspidx],st->T[st->kspidx],flg);
171:     KSPSetUp(st->ksp);
172:   }
173:   return(0);
174: }

178: PetscErrorCode STSetFromOptions_Sinvert(ST st)
179: {
181:   PC             pc;
182:   PCType         pctype;
183:   KSPType        ksptype;

186:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
187:   KSPGetPC(st->ksp,&pc);
188:   KSPGetType(st->ksp,&ksptype);
189:   PCGetType(pc,&pctype);
190:   if (!pctype && !ksptype) {
191:     if (st->shift_matrix == ST_MATMODE_SHELL) {
192:       /* in shell mode use GMRES with Jacobi as the default */
193:       KSPSetType(st->ksp,KSPGMRES);
194:       PCSetType(pc,PCJACOBI);
195:     } else {
196:       /* use direct solver as default */
197:       KSPSetType(st->ksp,KSPPREONLY);
198:       PCSetType(pc,PCREDUNDANT);
199:     }
200:   }
201:   return(0);
202: }

206: PETSC_EXTERN PetscErrorCode STCreate_Sinvert(ST st)
207: {
209:   st->ops->apply           = STApply_Sinvert;
210:   st->ops->getbilinearform = STGetBilinearForm_Default;
211:   st->ops->applytrans      = STApplyTranspose_Sinvert;
212:   st->ops->postsolve       = STPostSolve_Sinvert;
213:   st->ops->backtransform   = STBackTransform_Sinvert;
214:   st->ops->setup           = STSetUp_Sinvert;
215:   st->ops->setshift        = STSetShift_Sinvert;
216:   st->ops->setfromoptions  = STSetFromOptions_Sinvert;
217:   st->ops->checknullspace  = STCheckNullSpace_Default;
218:   return(0);
219: }