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

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

 14:  #include src/st/stimpl.h

 16: typedef struct {
 17:   PetscTruth  left;
 18:   Vec         w2;
 19: } ST_FOLD;

 23: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
 24: {
 26:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

 29:   if (st->B) {
 30:     /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
 31:     MatMult(st->A,x,st->w);
 32:     STAssociatedKSPSolve(st,st->w,ctx->w2);
 33:     if (st->sigma != 0.0) {
 34:       VecAXPY(ctx->w2,-st->sigma,x);
 35:     }
 36:     MatMult(st->A,ctx->w2,st->w);
 37:     STAssociatedKSPSolve(st,st->w,y);
 38:     if (st->sigma != 0.0) {
 39:       VecAXPY(y,-st->sigma,ctx->w2);
 40:     }
 41:   } else {
 42:     /* standard eigenproblem: y = (A + sI)^2 x */
 43:     MatMult(st->A,x,st->w);
 44:     if (st->sigma != 0.0) {
 45:       VecAXPY(st->w,-st->sigma,x);
 46:     }
 47:     MatMult(st->A,st->w,y);
 48:     if (st->sigma != 0.0) {
 49:       VecAXPY(y,-st->sigma,st->w);
 50:     }
 51:   }
 52:   return(0);
 53: }

 57: PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
 58: {
 60:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

 63:   if (st->B) {
 64:     /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
 65:     STAssociatedKSPSolveTranspose(st,x,st->w);
 66:     MatMult(st->A,st->w,ctx->w2);
 67:     if (st->sigma != 0.0) {
 68:       VecAXPY(ctx->w2,-st->sigma,x);
 69:     }
 70:     STAssociatedKSPSolveTranspose(st,ctx->w2,st->w);
 71:     MatMult(st->A,st->w,y);
 72:     if (st->sigma != 0.0) {
 73:       VecAXPY(y,-st->sigma,ctx->w2);
 74:     }
 75:   } else {
 76:     /* standard eigenproblem: y = (A^T + sI)^2 x */
 77:     MatMultTranspose(st->A,x,st->w);
 78:     if (st->sigma != 0.0) {
 79:       VecAXPY(st->w,-st->sigma,x);
 80:     }
 81:     MatMultTranspose(st->A,st->w,y);
 82:     if (st->sigma != 0.0) {
 83:       VecAXPY(y,-st->sigma,st->w);
 84:     }
 85:   }
 86:   return(0);
 87: }

 91: PetscErrorCode STBackTransform_Fold(ST st,PetscScalar *eigr,PetscScalar *eigi)
 92: {
 93:   ST_FOLD *ctx = (ST_FOLD *) st->data;
 97: #if !defined(PETSC_USE_COMPLEX)
 98:   if (*eigi == 0) {
 99: #endif
100:     if (ctx->left) *eigr = st->sigma - PetscSqrtScalar(*eigr);
101:     else *eigr = st->sigma + PetscSqrtScalar(*eigr);
102: #if !defined(PETSC_USE_COMPLEX)
103:   } else {
104:     PetscScalar r,x,y;
105:     r = PetscSqrtScalar(*eigr * *eigr + *eigi * *eigi);
106:     x = PetscSqrtScalar((r + *eigr) / 2);
107:     y = PetscSqrtScalar((r - *eigr) / 2);
108:     if (*eigi < 0) y = - y;
109:     if (ctx->left) {
110:       *eigr = st->sigma - x;
111:       *eigi = - y;
112:     } else {
113:       *eigr = st->sigma + x;
114:       *eigi = y;
115:     }
116:   }
117: #endif
118:   return(0);
119: }

123: PetscErrorCode STSetUp_Fold(ST st)
124: {
126:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

129:   if (st->B) {
130:     KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);
131:     KSPSetUp(st->ksp);
132:     if (ctx->w2) { VecDestroy(ctx->w2); }
133:     MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
134:   }
135:   return(0);
136: }

141: PetscErrorCode STFoldSetLeftSide_Fold(ST st,PetscTruth left)
142: {
143:   ST_FOLD *ctx = (ST_FOLD *) st->data;

146:   ctx->left = left;
147:   return(0);
148: }

153: /*@
154:    STFoldSetLeftSide - Sets a flag to compute eigenvalues on the left side of shift.

156:    Collective on ST

158:    Input Parameters:
159: +  st  - the spectral transformation context
160: -  left - if true compute eigenvalues on the left side 

162:    Options Database Key:
163: .  -st_fold_leftside - Sets the value of the flag

165:    Level: intermediate

167: .seealso: STSetShift()
168: @*/
169: PetscErrorCode STFoldSetLeftSide(ST st,PetscTruth left)
170: {
171:   PetscErrorCode ierr, (*f)(ST,PetscTruth);

175:   PetscObjectQueryFunction((PetscObject)st,"STFoldSetLeftSide_C",(void (**)(void))&f);
176:   if (f) {
177:     (*f)(st,left);
178:   }
179:   return(0);
180: }

184: PetscErrorCode STView_Fold(ST st,PetscViewer viewer)
185: {
187:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

190:   if (ctx->left) {
191:     PetscViewerASCIIPrintf(viewer,"  computing eigenvalues on left side of shift\n");
192:   }
193:   if (st->B) {
194:     STView_Default(st,viewer);
195:   }
196:   return(0);
197: }

201: PetscErrorCode STSetFromOptions_Fold(ST st)
202: {
204:   ST_FOLD      *ctx = (ST_FOLD *) st->data;

207:   PetscOptionsHead("ST Fold Options");
208:   PetscOptionsTruth("-st_fold_leftside","Compute eigenvalues on left side of shift","STFoldSetLeftSide",ctx->left,&ctx->left,PETSC_NULL);
209:   PetscOptionsTail();
210:   return(0);
211: }

215: PetscErrorCode STDestroy_Fold(ST st)
216: {
218:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

221:   if (ctx->w2) { VecDestroy(ctx->w2); }
222:   PetscFree(ctx);
223:   return(0);
224: }

229: PetscErrorCode STCreate_Fold(ST st)
230: {
232:   ST_FOLD        *ctx;


236:   PetscNew(ST_FOLD,&ctx);
237:   PetscLogObjectMemory(st,sizeof(ST_FOLD));
238:   st->data                  = (void *) ctx;

240:   st->ops->apply           = STApply_Fold;
241:   st->ops->getbilinearform = STGetBilinearForm_Default;
242:   st->ops->applytrans      = STApplyTranspose_Fold;
243:   st->ops->backtr           = STBackTransform_Fold;
244:   st->ops->setup           = STSetUp_Fold;
245:   st->ops->view            = STView_Fold;
246:   st->ops->setfromoptions  = STSetFromOptions_Fold;
247:   st->ops->destroy           = STDestroy_Fold;
248:   st->checknullspace           = 0;
249: 
250:   ctx->left            = PETSC_FALSE;
251: 
252:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STFoldSetLeftSide_C","STFoldSetLeftSide_Fold",
253:                     STFoldSetLeftSide_Fold);

255:   return(0);
256: }