Actual source code: fold.c

  1: /*
  2:     Folding spectral transformation, applies (A + sigma I)^2 as operator, or 
  3:     inv(B)(A + sigma I)^2 for generalized problems

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2012, Universitat 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 <slepc-private/stimpl.h>          /*I "slepcst.h" I*/

 27: typedef struct {
 28:   Vec         w2;
 29: } ST_FOLD;

 33: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
 34: {
 36:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

 39:   if (st->B) {
 40:     /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
 41:     MatMult(st->A,x,st->w);
 42:     STAssociatedKSPSolve(st,st->w,ctx->w2);
 43:     if (st->sigma != 0.0) {
 44:       VecAXPY(ctx->w2,-st->sigma,x);
 45:     }
 46:     MatMult(st->A,ctx->w2,st->w);
 47:     STAssociatedKSPSolve(st,st->w,y);
 48:     if (st->sigma != 0.0) {
 49:       VecAXPY(y,-st->sigma,ctx->w2);
 50:     }
 51:   } else {
 52:     /* standard eigenproblem: y = (A + sI)^2 x */
 53:     MatMult(st->A,x,st->w);
 54:     if (st->sigma != 0.0) {
 55:       VecAXPY(st->w,-st->sigma,x);
 56:     }
 57:     MatMult(st->A,st->w,y);
 58:     if (st->sigma != 0.0) {
 59:       VecAXPY(y,-st->sigma,st->w);
 60:     }
 61:   }
 62:   return(0);
 63: }

 67: PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
 68: {
 70:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

 73:   if (st->B) {
 74:     /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
 75:     STAssociatedKSPSolveTranspose(st,x,st->w);
 76:     MatMult(st->A,st->w,ctx->w2);
 77:     if (st->sigma != 0.0) {
 78:       VecAXPY(ctx->w2,-st->sigma,x);
 79:     }
 80:     STAssociatedKSPSolveTranspose(st,ctx->w2,st->w);
 81:     MatMult(st->A,st->w,y);
 82:     if (st->sigma != 0.0) {
 83:       VecAXPY(y,-st->sigma,ctx->w2);
 84:     }
 85:   } else {
 86:     /* standard eigenproblem: y = (A^T + sI)^2 x */
 87:     MatMultTranspose(st->A,x,st->w);
 88:     if (st->sigma != 0.0) {
 89:       VecAXPY(st->w,-st->sigma,x);
 90:     }
 91:     MatMultTranspose(st->A,st->w,y);
 92:     if (st->sigma != 0.0) {
 93:       VecAXPY(y,-st->sigma,st->w);
 94:     }
 95:   }
 96:   return(0);
 97: }

101: PetscErrorCode STBackTransform_Fold(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
102: {
103:   PetscInt    j;
104: #if !defined(PETSC_USE_COMPLEX)
105:   PetscScalar r,x,y;
106: #endif

109:   for (j=0;j<n;j++) {
110: #if !defined(PETSC_USE_COMPLEX)
111:     if (eigi[j] == 0.0) {
112: #endif
113:       eigr[j] = st->sigma + PetscSqrtScalar(eigr[j]);
114: #if !defined(PETSC_USE_COMPLEX)
115:     } else {
116:       r = SlepcAbsEigenvalue(eigr[j],eigi[j]);
117:       x = PetscSqrtScalar((r + eigr[j]) / 2);
118:       y = PetscSqrtScalar((r - eigr[j]) / 2);
119:       if (eigi[j] < 0.0) y = - y;
120:       eigr[j] = st->sigma + x;
121:       eigi[j] = y;
122:     }
123: #endif
124:   }
125:   return(0);
126: }

130: PetscErrorCode STSetUp_Fold(ST st)
131: {
133:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

136:   /* if the user did not set the shift, use the target value */
137:   if (!st->sigma_set) st->sigma = st->defsigma;

139:   if (st->B) {
140:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
141:     KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);
142:     KSPSetUp(st->ksp);
143:     VecDestroy(&ctx->w2);
144:     MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
145:   }
146:   return(0);
147: }

151: PetscErrorCode STSetFromOptions_Fold(ST st)
152: {
154:   PC             pc;
155:   const PCType   pctype;
156:   const KSPType  ksptype;

159:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
160:   KSPGetPC(st->ksp,&pc);
161:   KSPGetType(st->ksp,&ksptype);
162:   PCGetType(pc,&pctype);
163:   if (!pctype && !ksptype) {
164:     if (st->shift_matrix == ST_MATMODE_SHELL) {
165:       /* in shell mode use GMRES with Jacobi as the default */
166:       KSPSetType(st->ksp,KSPGMRES);
167:       PCSetType(pc,PCJACOBI);
168:     } else {
169:       /* use direct solver as default */
170:       KSPSetType(st->ksp,KSPPREONLY);
171:       PCSetType(pc,PCREDUNDANT);
172:     }
173:   }
174:   return(0);
175: }

179: PetscErrorCode STReset_Fold(ST st)
180: {
182:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

185:   VecDestroy(&ctx->w2);
186:   return(0);
187: }

191: PetscErrorCode STDestroy_Fold(ST st)
192: {

196:   PetscFree(st->data);
197:   return(0);
198: }

200: EXTERN_C_BEGIN
203: PetscErrorCode STCreate_Fold(ST st)
204: {

208:   PetscNewLog(st,ST_FOLD,&st->data);
209:   st->ops->apply           = STApply_Fold;
210:   st->ops->getbilinearform = STGetBilinearForm_Default;
211:   st->ops->applytrans      = STApplyTranspose_Fold;
212:   st->ops->backtr          = STBackTransform_Fold;
213:   st->ops->setup           = STSetUp_Fold;
214:   st->ops->setfromoptions  = STSetFromOptions_Fold;
215:   st->ops->destroy         = STDestroy_Fold;
216:   st->ops->reset           = STReset_Fold;
217:   return(0);
218: }
219: EXTERN_C_END