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

  8:       This file is part of SLEPc. See the README file for conditions of use
  9:       and additional information.
 10:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 11: */

 13:  #include src/st/stimpl.h

 17: PetscErrorCode STApply_Sinvert(ST st,Vec x,Vec y)
 18: {

 22:   if (st->B) {
 23:     /* generalized eigenproblem: y = (A - sB)^-1 B x */
 24:     MatMult(st->B,x,st->w);
 25:     STAssociatedKSPSolve(st,st->w,y);
 26:   }
 27:   else {
 28:     /* standard eigenproblem: y = (A - sI)^-1 x */
 29:     STAssociatedKSPSolve(st,x,y);
 30:   }
 31:   return(0);
 32: }

 36: PetscErrorCode STApplyTranspose_Sinvert(ST st,Vec x,Vec y)
 37: {

 41:   if (st->B) {
 42:     /* generalized eigenproblem: y = B^T (A - sB)^-T x */
 43:     STAssociatedKSPSolveTranspose(st,x,st->w);
 44:     MatMultTranspose(st->B,st->w,y);
 45:   }
 46:   else {
 47:     /* standard eigenproblem: y = (A - sI)^-T x */
 48:     STAssociatedKSPSolveTranspose(st,x,y);
 49:   }
 50:   return(0);
 51: }

 55: PetscErrorCode STBackTransform_Sinvert(ST st,PetscScalar *eigr,PetscScalar *eigi)
 56: {
 57: #ifndef PETSC_USE_COMPLEX
 58:   PetscScalar t;
 62:   if (*eigi == 0) *eigr = 1.0 / *eigr + st->sigma;
 63:   else {
 64:     t = *eigr * *eigr + *eigi * *eigi;
 65:     *eigr = *eigr / t + st->sigma;
 66:     *eigi = - *eigi / t;
 67:   }
 68: #else
 71:   *eigr = 1.0 / *eigr + st->sigma;
 72: #endif
 73:   return(0);
 74: }

 78: PetscErrorCode STPostSolve_Sinvert(ST st)
 79: {

 83:   if (st->shift_matrix == STMATMODE_INPLACE) {
 84:     if( st->B ) {
 85:       MatAXPY(st->A,st->sigma,st->B,st->str);
 86:     } else {
 87:       MatShift(st->A,st->sigma);
 88:     }
 89:     st->setupcalled = 0;
 90:   }
 91:   return(0);
 92: }

 96: PetscErrorCode STSetUp_Sinvert(ST st)
 97: {

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

103:   switch (st->shift_matrix) {
104:   case STMATMODE_INPLACE:
105:     st->mat = PETSC_NULL;
106:     if (st->sigma != 0.0) {
107:       if (st->B) {
108:         MatAXPY(st->A,-st->sigma,st->B,st->str);
109:       } else {
110:         MatShift(st->A,-st->sigma);
111:       }
112:     }
113:     KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
114:     break;
115:   case STMATMODE_SHELL:
116:     STMatShellCreate(st,&st->mat);
117:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
118:     break;
119:   default:
120:     MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
121:     if (st->sigma != 0.0) {
122:       if (st->B) {
123:         MatAXPY(st->mat,-st->sigma,st->B,st->str);
124:       } else {
125:         MatShift(st->mat,-st->sigma);
126:       }
127:       KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
128:     } else {
129:       st->mat = PETSC_NULL;
130:       KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
131:     }
132:   }

134:   KSPSetUp(st->ksp);
135:   return(0);
136: }

140: PetscErrorCode STSetShift_Sinvert(ST st,PetscScalar newshift)
141: {
143:   MatStructure   flg;


147:   /* Nothing to be done if STSetUp has not been called yet */
148:   if (!st->setupcalled) return(0);
149: 
150:   /* Check if the new KSP matrix has the same zero structure */
151:   if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
152:     flg = DIFFERENT_NONZERO_PATTERN;
153:   } else {
154:     flg = SAME_NONZERO_PATTERN;
155:   }

157:   switch (st->shift_matrix) {
158:   case STMATMODE_INPLACE:
159:     /* Undo previous operations */
160:     if (st->sigma != 0.0) {
161:       if (st->B) {
162:         MatAXPY(st->A,st->sigma,st->B,st->str);
163:       } else {
164:         MatShift(st->A,st->sigma);
165:       }
166:     }
167:     /* Apply new shift */
168:     if (newshift != 0.0) {
169:       if (st->B) {
170:         MatAXPY(st->A,-newshift,st->B,st->str);
171:       } else {
172:         MatShift(st->A,-newshift);
173:       }
174:     }
175:     KSPSetOperators(st->ksp,st->A,st->A,flg);
176:     break;
177:   case STMATMODE_SHELL:
178:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
179:     break;
180:   default:
181:     if (st->mat) {
182:       MatCopy(st->A,st->mat,SUBSET_NONZERO_PATTERN);
183:     } else {
184:       MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
185:     }
186:     if (newshift != 0.0) {
187:       if (st->B) {
188:         MatAXPY(st->mat,-newshift,st->B,st->str);
189:       } else {
190:         MatShift(st->mat,-newshift);
191:       }
192:     }
193:     KSPSetOperators(st->ksp,st->mat,st->mat,flg);
194:   }
195:   st->sigma = newshift;
196:   KSPSetUp(st->ksp);
197:   return(0);
198: }

203: PetscErrorCode STCreate_Sinvert(ST st)
204: {
206:   st->data                 = 0;

208:   st->ops->apply           = STApply_Sinvert;
209:   st->ops->getbilinearform = STGetBilinearForm_Default;
210:   st->ops->applytrans      = STApplyTranspose_Sinvert;
211:   st->ops->postsolve       = STPostSolve_Sinvert;
212:   st->ops->backtr          = STBackTransform_Sinvert;
213:   st->ops->setup           = STSetUp_Sinvert;
214:   st->ops->setshift        = STSetShift_Sinvert;
215:   st->ops->view            = STView_Default;
216: 
217:   st->checknullspace      = STCheckNullSpace_Default;

219:   return(0);
220: }