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

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

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

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

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

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

102: PetscErrorCode STBackTransform_Fold(ST st,PetscScalar *eigr,PetscScalar *eigi)
103: {
104:   ST_FOLD *ctx = (ST_FOLD *) st->data;
108: #if !defined(PETSC_USE_COMPLEX)
109:   if (*eigi == 0) {
110: #endif
111:     if (ctx->left) *eigr = st->sigma - PetscSqrtScalar(*eigr);
112:     else *eigr = st->sigma + PetscSqrtScalar(*eigr);
113: #if !defined(PETSC_USE_COMPLEX)
114:   } else {
115:     PetscScalar r,x,y;
116:     r = PetscSqrtScalar(*eigr * *eigr + *eigi * *eigi);
117:     x = PetscSqrtScalar((r + *eigr) / 2);
118:     y = PetscSqrtScalar((r - *eigr) / 2);
119:     if (*eigi < 0) y = - y;
120:     if (ctx->left) {
121:       *eigr = st->sigma - x;
122:       *eigi = - y;
123:     } else {
124:       *eigr = st->sigma + x;
125:       *eigi = y;
126:     }
127:   }
128: #endif
129:   return(0);
130: }

134: PetscErrorCode STSetUp_Fold(ST st)
135: {
137:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

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

152: PetscErrorCode STFoldSetLeftSide_Fold(ST st,PetscTruth left)
153: {
154:   ST_FOLD *ctx = (ST_FOLD *) st->data;

157:   ctx->left = left;
158:   return(0);
159: }

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

167:    Collective on ST

169:    Input Parameters:
170: +  st  - the spectral transformation context
171: -  left - if true compute eigenvalues on the left side 

173:    Options Database Key:
174: .  -st_fold_leftside - Sets the value of the flag

176:    Level: intermediate

178: .seealso: STSetShift()
179: @*/
180: PetscErrorCode STFoldSetLeftSide(ST st,PetscTruth left)
181: {
182:   PetscErrorCode ierr, (*f)(ST,PetscTruth);

186:   PetscObjectQueryFunction((PetscObject)st,"STFoldSetLeftSide_C",(void (**)(void))&f);
187:   if (f) {
188:     (*f)(st,left);
189:   }
190:   return(0);
191: }

195: PetscErrorCode STView_Fold(ST st,PetscViewer viewer)
196: {
198:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

201:   if (ctx->left) {
202:     PetscViewerASCIIPrintf(viewer,"  computing eigenvalues on left side of shift\n");
203:   }
204:   if (st->B) {
205:     STView_Default(st,viewer);
206:   }
207:   return(0);
208: }

212: PetscErrorCode STSetFromOptions_Fold(ST st)
213: {
215:   ST_FOLD      *ctx = (ST_FOLD *) st->data;

218:   PetscOptionsHead("ST Fold Options");
219:   PetscOptionsTruth("-st_fold_leftside","Compute eigenvalues on left side of shift","STFoldSetLeftSide",ctx->left,&ctx->left,PETSC_NULL);
220:   PetscOptionsTail();
221:   return(0);
222: }

226: PetscErrorCode STDestroy_Fold(ST st)
227: {
229:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

232:   if (ctx->w2) { VecDestroy(ctx->w2); }
233:   PetscFree(ctx);
234:   return(0);
235: }

240: PetscErrorCode STCreate_Fold(ST st)
241: {
243:   ST_FOLD        *ctx;


247:   PetscNew(ST_FOLD,&ctx);
248:   PetscLogObjectMemory(st,sizeof(ST_FOLD));
249:   st->data                  = (void *) ctx;

251:   st->ops->apply           = STApply_Fold;
252:   st->ops->getbilinearform = STGetBilinearForm_Default;
253:   st->ops->applytrans      = STApplyTranspose_Fold;
254:   st->ops->backtr           = STBackTransform_Fold;
255:   st->ops->setup           = STSetUp_Fold;
256:   st->ops->view            = STView_Fold;
257:   st->ops->setfromoptions  = STSetFromOptions_Fold;
258:   st->ops->destroy           = STDestroy_Fold;
259:   st->checknullspace           = 0;
260: 
261:   ctx->left            = PETSC_FALSE;
262: 
263:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STFoldSetLeftSide_C","STFoldSetLeftSide_Fold",
264:                     STFoldSetLeftSide_Fold);

266:   return(0);
267: }