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-2009, Universidad 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

 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,PetscScalar *eigr,PetscScalar *eigi)
 67: {
 68: #ifndef PETSC_USE_COMPLEX
 69:   PetscScalar t;
 73:   if (*eigi == 0) *eigr = 1.0 / *eigr + st->sigma;
 74:   else {
 75:     t = *eigr * *eigr + *eigi * *eigi;
 76:     *eigr = *eigr / t + st->sigma;
 77:     *eigi = - *eigi / t;
 78:   }
 79: #else
 82:   *eigr = 1.0 / *eigr + st->sigma;
 83: #endif
 84:   return(0);
 85: }

 89: PetscErrorCode STPostSolve_Sinvert(ST st)
 90: {

 94:   if (st->shift_matrix == STMATMODE_INPLACE) {
 95:     if( st->B ) {
 96:       MatAXPY(st->A,st->sigma,st->B,st->str);
 97:     } else {
 98:       MatShift(st->A,st->sigma);
 99:     }
100:     st->setupcalled = 0;
101:   }
102:   return(0);
103: }

107: PetscErrorCode STSetUp_Sinvert(ST st)
108: {

112:   if (st->mat) { MatDestroy(st->mat); }

114:   switch (st->shift_matrix) {
115:   case STMATMODE_INPLACE:
116:     st->mat = PETSC_NULL;
117:     if (st->sigma != 0.0) {
118:       if (st->B) {
119:         MatAXPY(st->A,-st->sigma,st->B,st->str);
120:       } else {
121:         MatShift(st->A,-st->sigma);
122:       }
123:     }
124:     KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
125:     break;
126:   case STMATMODE_SHELL:
127:     STMatShellCreate(st,&st->mat);
128:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
129:     break;
130:   default:
131:     if (st->sigma != 0.0) {
132:       MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
133:       if (st->B) {
134:         MatAXPY(st->mat,-st->sigma,st->B,st->str);
135:       } else {
136:         MatShift(st->mat,-st->sigma);
137:       }
138:       KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
139:     } else {
140:       st->mat = PETSC_NULL;
141:       KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
142:     }
143:   }

145:   KSPSetUp(st->ksp);
146:   return(0);
147: }

151: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
152: {
154:   MatStructure   flg;


158:   /* Nothing to be done if STSetUp has not been called yet */
159:   if (!st->setupcalled) return(0);
160: 
161:   /* Check if the new KSP matrix has the same zero structure */
162:   if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
163:     flg = DIFFERENT_NONZERO_PATTERN;
164:   } else {
165:     flg = SAME_NONZERO_PATTERN;
166:   }

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

214: PetscErrorCode STCreate_Sinvert(ST st)
215: {
217:   st->data                 = 0;

219:   st->ops->apply           = STApply_Sinvert;
220:   st->ops->getbilinearform = STGetBilinearForm_Default;
221:   st->ops->applytrans      = STApplyTranspose_Sinvert;
222:   st->ops->postsolve       = STPostSolve_Sinvert;
223:   st->ops->backtr          = STBackTransform_Sinvert;
224:   st->ops->setup           = STSetUp_Sinvert;
225:   st->ops->setshift        = STSetShift_Sinvert;
226:   st->ops->view            = STView_Default;
227: 
228:   st->checknullspace      = STCheckNullSpace_Default;

230:   return(0);
231: }