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-2013, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 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: {
 35:   PetscErrorCode     ierr;
 36:   ST_FOLD            *ctx = (ST_FOLD*)st->data;
 37:   PetscInt           its;
 38:   KSPConvergedReason reason;

 41:   if (st->nmat>1) {
 42:     /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
 43:     MatMult(st->A[0],x,st->w);
 44:     KSPSolve(st->ksp,st->w,ctx->w2);
 45:     KSPGetConvergedReason(st->ksp,&reason);
 46:     if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 47:     KSPGetIterationNumber(st->ksp,&its);
 48:     st->lineariterations += its;
 49:     if (st->sigma != 0.0) {
 50:       VecAXPY(ctx->w2,-st->sigma,x);
 51:     }
 52:     MatMult(st->A[0],ctx->w2,st->w);
 53:     KSPSolve(st->ksp,st->w,y);
 54:     KSPGetConvergedReason(st->ksp,&reason);
 55:     if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 56:     KSPGetIterationNumber(st->ksp,&its);
 57:     st->lineariterations += its;
 58:     if (st->sigma != 0.0) {
 59:       VecAXPY(y,-st->sigma,ctx->w2);
 60:     }
 61:   } else {
 62:     /* standard eigenproblem: y = (A + sI)^2 x */
 63:     MatMult(st->A[0],x,st->w);
 64:     if (st->sigma != 0.0) {
 65:       VecAXPY(st->w,-st->sigma,x);
 66:     }
 67:     MatMult(st->A[0],st->w,y);
 68:     if (st->sigma != 0.0) {
 69:       VecAXPY(y,-st->sigma,st->w);
 70:     }
 71:   }
 72:   return(0);
 73: }

 77: PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
 78: {
 79:   PetscErrorCode     ierr;
 80:   ST_FOLD            *ctx = (ST_FOLD*)st->data;
 81:   PetscInt           its;
 82:   KSPConvergedReason reason;

 85:   if (st->nmat>1) {
 86:     /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
 87:     KSPSolveTranspose(st->ksp,x,st->w);
 88:     KSPGetConvergedReason(st->ksp,&reason);
 89:     if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 90:     KSPGetIterationNumber(st->ksp,&its);
 91:     st->lineariterations += its;
 92:     MatMult(st->A[0],st->w,ctx->w2);
 93:     if (st->sigma != 0.0) {
 94:       VecAXPY(ctx->w2,-st->sigma,x);
 95:     }
 96:     KSPSolveTranspose(st->ksp,ctx->w2,st->w);
 97:     KSPGetConvergedReason(st->ksp,&reason);
 98:     if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 99:     KSPGetIterationNumber(st->ksp,&its);
100:     st->lineariterations += its;
101:     MatMult(st->A[0],st->w,y);
102:     if (st->sigma != 0.0) {
103:       VecAXPY(y,-st->sigma,ctx->w2);
104:     }
105:   } else {
106:     /* standard eigenproblem: y = (A^T + sI)^2 x */
107:     MatMultTranspose(st->A[0],x,st->w);
108:     if (st->sigma != 0.0) {
109:       VecAXPY(st->w,-st->sigma,x);
110:     }
111:     MatMultTranspose(st->A[0],st->w,y);
112:     if (st->sigma != 0.0) {
113:       VecAXPY(y,-st->sigma,st->w);
114:     }
115:   }
116:   return(0);
117: }

121: PetscErrorCode STBackTransform_Fold(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
122: {
123:   PetscInt    j;
124: #if !defined(PETSC_USE_COMPLEX)
125:   PetscScalar r,x,y;
126: #endif

129:   for (j=0;j<n;j++) {
130: #if !defined(PETSC_USE_COMPLEX)
131:     if (eigi[j] == 0.0) {
132: #endif
133:       eigr[j] = st->sigma + PetscSqrtScalar(eigr[j]);
134: #if !defined(PETSC_USE_COMPLEX)
135:     } else {
136:       r = SlepcAbsEigenvalue(eigr[j],eigi[j]);
137:       x = PetscSqrtScalar((r + eigr[j]) / 2);
138:       y = PetscSqrtScalar((r - eigr[j]) / 2);
139:       if (eigi[j] < 0.0) y = - y;
140:       eigr[j] = st->sigma + x;
141:       eigi[j] = y;
142:     }
143: #endif
144:   }
145:   return(0);
146: }

150: PetscErrorCode STSetUp_Fold(ST st)
151: {
153:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

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

159:   if (st->nmat>1) {
160:     if (!st->ksp) { STGetKSP(st,&st->ksp); }
161:     KSPSetOperators(st->ksp,st->A[1],st->A[1],DIFFERENT_NONZERO_PATTERN);
162:     KSPSetUp(st->ksp);
163:     VecDestroy(&ctx->w2);
164:     MatGetVecs(st->A[1],&ctx->w2,NULL);
165:     PetscLogObjectParent(st,ctx->w2);
166:   }
167:   return(0);
168: }

172: PetscErrorCode STSetFromOptions_Fold(ST st)
173: {
175:   PC             pc;
176:   PCType         pctype;
177:   KSPType        ksptype;

180:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
181:   KSPGetPC(st->ksp,&pc);
182:   KSPGetType(st->ksp,&ksptype);
183:   PCGetType(pc,&pctype);
184:   if (!pctype && !ksptype) {
185:     if (st->shift_matrix == ST_MATMODE_SHELL) {
186:       /* in shell mode use GMRES with Jacobi as the default */
187:       KSPSetType(st->ksp,KSPGMRES);
188:       PCSetType(pc,PCJACOBI);
189:     } else {
190:       /* use direct solver as default */
191:       KSPSetType(st->ksp,KSPPREONLY);
192:       PCSetType(pc,PCREDUNDANT);
193:     }
194:   }
195:   return(0);
196: }

200: PetscErrorCode STReset_Fold(ST st)
201: {
203:   ST_FOLD        *ctx = (ST_FOLD*)st->data;

206:   VecDestroy(&ctx->w2);
207:   return(0);
208: }

212: PetscErrorCode STDestroy_Fold(ST st)
213: {

217:   PetscFree(st->data);
218:   return(0);
219: }

223: PETSC_EXTERN PetscErrorCode STCreate_Fold(ST st)
224: {

228:   PetscNewLog(st,ST_FOLD,&st->data);
229:   st->ops->apply           = STApply_Fold;
230:   st->ops->getbilinearform = STGetBilinearForm_Default;
231:   st->ops->applytrans      = STApplyTranspose_Fold;
232:   st->ops->backtransform   = STBackTransform_Fold;
233:   st->ops->setup           = STSetUp_Fold;
234:   st->ops->setfromoptions  = STSetFromOptions_Fold;
235:   st->ops->destroy         = STDestroy_Fold;
236:   st->ops->reset           = STReset_Fold;
237:   return(0);
238: }