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

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

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

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

149: PetscErrorCode STView_Fold(ST st,PetscViewer viewer)
150: {

154:   if (st->B) {
155:     STView_Default(st,viewer);
156:   }
157:   return(0);
158: }

162: PetscErrorCode STSetFromOptions_Fold(ST st)
163: {
165:   PC             pc;
166:   const PCType   pctype;
167:   const KSPType  ksptype;


171:   KSPGetPC(st->ksp,&pc);
172:   KSPGetType(st->ksp,&ksptype);
173:   PCGetType(pc,&pctype);
174:   if (!pctype && !ksptype) {
175:     if (st->shift_matrix == ST_MATMODE_SHELL) {
176:       /* in shell mode use GMRES with Jacobi as the default */
177:       KSPSetType(st->ksp,KSPGMRES);
178:       PCSetType(pc,PCJACOBI);
179:     } else {
180:       /* use direct solver as default */
181:       KSPSetType(st->ksp,KSPPREONLY);
182:       PCSetType(pc,PCREDUNDANT);
183:     }
184:   }

186:   return(0);
187: }

191: PetscErrorCode STDestroy_Fold(ST st)
192: {
194:   ST_FOLD        *ctx = (ST_FOLD *) st->data;

197:   if (ctx->w2) { VecDestroy(ctx->w2); }
198:   PetscFree(ctx);
199:   return(0);
200: }

205: PetscErrorCode STCreate_Fold(ST st)
206: {
208:   ST_FOLD        *ctx;


212:   PetscNew(ST_FOLD,&ctx);
213:   PetscLogObjectMemory(st,sizeof(ST_FOLD));
214:   st->data                  = (void *) ctx;

216:   st->ops->apply           = STApply_Fold;
217:   st->ops->getbilinearform = STGetBilinearForm_Default;
218:   st->ops->applytrans      = STApplyTranspose_Fold;
219:   st->ops->backtr           = STBackTransform_Fold;
220:   st->ops->setup           = STSetUp_Fold;
221:   st->ops->view            = STView_Fold;
222:   st->ops->setfromoptions  = STSetFromOptions_Fold;
223:   st->ops->destroy           = STDestroy_Fold;
224:   st->checknullspace           = 0;
225: 
226:   return(0);
227: }