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-2010, Universidad Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.
 10:       
 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 private/stimpl.h

 29: PetscErrorCode STApply_Shift(ST st,Vec x,Vec y)
 30: {

 34:   if (st->B) {
 35:     /* generalized eigenproblem: y = (B^-1 A + sI) x */
 36:     MatMult(st->A,x,st->w);
 37:     STAssociatedKSPSolve(st,st->w,y);
 38:   }
 39:   else {
 40:     /* standard eigenproblem: y = (A + sI) x */
 41:     MatMult(st->A,x,y);
 42:   }
 43:   if (st->sigma != 0.0) {
 44:     VecAXPY(y,st->sigma,x);
 45:   }
 46:   return(0);
 47: }

 51: PetscErrorCode STApplyTranspose_Shift(ST st,Vec x,Vec y)
 52: {

 56:   if (st->B) {
 57:     /* generalized eigenproblem: y = (A^T B^-T + sI) x */
 58:     STAssociatedKSPSolveTranspose(st,x,st->w);
 59:     MatMultTranspose(st->A,st->w,y);
 60:   }
 61:   else {
 62:     /* standard eigenproblem: y = (A^T + sI) x */
 63:     MatMultTranspose(st->A,x,y);
 64:   }
 65:   if (st->sigma != 0.0) {
 66:     VecAXPY(y,st->sigma,x);
 67:   }
 68:   return(0);
 69: }

 73: PetscErrorCode STBackTransform_Shift(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 74: {
 75:   PetscInt j;
 78:   for (j=0;j<n;j++) {
 79:     eigr[j] -= st->sigma;
 80:   }
 81:   return(0);
 82: }

 86: PetscErrorCode STSetUp_Shift(ST st)
 87: {

 91:   if (st->B) {
 92:     KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);
 93:     KSPSetUp(st->ksp);
 94:   }
 95:   return(0);
 96: }

100: PetscErrorCode STView_Shift(ST st,PetscViewer viewer)
101: {

105:   if (st->B) {
106:     STView_Default(st,viewer);
107:   }
108:   return(0);
109: }

113: PetscErrorCode STSetFromOptions_Shift(ST st)
114: {
116:   PC             pc;
117:   const PCType   pctype;
118:   const KSPType  ksptype;


122:   KSPGetPC(st->ksp,&pc);
123:   KSPGetType(st->ksp,&ksptype);
124:   PCGetType(pc,&pctype);
125:   if (!pctype && !ksptype) {
126:     if (st->shift_matrix == ST_MATMODE_SHELL) {
127:       /* in shell mode use GMRES with Jacobi as the default */
128:       KSPSetType(st->ksp,KSPGMRES);
129:       PCSetType(pc,PCJACOBI);
130:     } else {
131:       /* use direct solver as default */
132:       KSPSetType(st->ksp,KSPPREONLY);
133:       PCSetType(pc,PCREDUNDANT);
134:     }
135:   }

137:   return(0);
138: }

143: PetscErrorCode STCreate_Shift(ST st)
144: {
146:   st->ops->apply           = STApply_Shift;
147:   st->ops->getbilinearform = STGetBilinearForm_Default;
148:   st->ops->applytrans      = STApplyTranspose_Shift;
149:   st->ops->backtr          = STBackTransform_Shift;
150:   st->ops->setfromoptions  = STSetFromOptions_Shift;
151:   st->ops->setup           = STSetUp_Shift;
152:   st->ops->view            = STView_Shift;
153:   st->checknullspace       = 0;
154:   return(0);
155: }