Actual source code: cayley.c

  1: /*
  2:       Implements the Cayley spectral transform.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.
  9:       
 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc-private/stimpl.h>          /*I "slepcst.h" I*/

 26: typedef struct {
 27:   PetscScalar nu;
 28:   PetscBool   nu_set;
 29:   Vec         w2;
 30: } ST_CAYLEY;

 34: PetscErrorCode STApply_Cayley(ST st,Vec x,Vec y)
 35: {
 37:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
 38:   PetscScalar    nu = ctx->nu;
 39: 
 41:   if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };

 43:   if (st->B) {
 44:     /* generalized eigenproblem: y = (A - sB)^-1 (A + tB)x */
 45:     MatMult(st->A,x,st->w);
 46:     MatMult(st->B,x,ctx->w2);
 47:     VecAXPY(st->w,nu,ctx->w2);
 48:     STAssociatedKSPSolve(st,st->w,y);
 49:   }
 50:   else {
 51:     /* standard eigenproblem: y = (A - sI)^-1 (A + tI)x */
 52:     MatMult(st->A,x,st->w);
 53:     VecAXPY(st->w,nu,x);
 54:     STAssociatedKSPSolve(st,st->w,y);
 55:   }
 56:   return(0);
 57: }

 61: PetscErrorCode STApplyTranspose_Cayley(ST st,Vec x,Vec y)
 62: {
 64:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
 65:   PetscScalar    nu = ctx->nu;
 66: 
 68:   if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };

 70:   if (st->B) {
 71:     /* generalized eigenproblem: y = (A + tB)^T (A - sB)^-T x */
 72:     STAssociatedKSPSolveTranspose(st,x,st->w);
 73:     MatMultTranspose(st->A,st->w,y);
 74:     MatMultTranspose(st->B,st->w,ctx->w2);
 75:     VecAXPY(y,nu,ctx->w2);
 76:   }
 77:   else {
 78:     /* standard eigenproblem: y =  (A + tI)^T (A - sI)^-T x */
 79:     STAssociatedKSPSolveTranspose(st,x,st->w);
 80:     MatMultTranspose(st->A,st->w,y);
 81:     VecAXPY(y,nu,st->w);
 82:   }
 83:   return(0);
 84: }

 88: PetscErrorCode STBilinearMatMult_Cayley(Mat B,Vec x,Vec y)
 89: {
 91:   ST             st;
 92:   ST_CAYLEY      *ctx;
 93:   PetscScalar    nu;
 94: 
 96:   MatShellGetContext(B,(void**)&st);
 97:   ctx = (ST_CAYLEY*)st->data;
 98:   nu = ctx->nu;
 99: 
100:   if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };

102:   if (st->B) {
103:     /* generalized eigenproblem: y = (A + tB)x */
104:     MatMult(st->A,x,y);
105:     MatMult(st->B,x,ctx->w2);
106:     VecAXPY(y,nu,ctx->w2);
107:   }
108:   else {
109:     /* standard eigenproblem: y = (A + tI)x */
110:     MatMult(st->A,x,y);
111:     VecAXPY(y,nu,x);
112:   }
113:   return(0);
114: }

118: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
119: {
121:   PetscInt       n,m;

124:   MatGetLocalSize(st->B,&n,&m);
125:   MatCreateShell(((PetscObject)st)->comm,n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,B);
126:   MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))STBilinearMatMult_Cayley);
127:   return(0);
128: }

132: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
133: {
134:   ST_CAYLEY   *ctx = (ST_CAYLEY*)st->data;
135:   PetscInt    j;
136: #if !defined(PETSC_USE_COMPLEX)
137:   PetscScalar t,i,r;
138: #endif

141: #if !defined(PETSC_USE_COMPLEX)
142:   for (j=0;j<n;j++) {
143:     if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
144:     else {
145:       r = eigr[j];
146:       i = eigi[j];
147:       r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
148:       i = - st->sigma * i - ctx->nu * i;
149:       t = i * i + r * (r - 2.0) + 1.0;
150:       eigr[j] = r / t;
151:       eigi[j] = i / t;
152:     }
153:   }
154: #else
155:   for (j=0;j<n;j++) {
156:     eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
157:   }
158: #endif
159:   return(0);
160: }

164: PetscErrorCode STPostSolve_Cayley(ST st)
165: {

169:   if (st->shift_matrix == ST_MATMODE_INPLACE) {
170:     if (st->B) {
171:       MatAXPY(st->A,st->sigma,st->B,st->str);
172:     } else {
173:       MatShift(st->A,st->sigma);
174:     }
175:     st->setupcalled = 0;
176:   }
177:   return(0);
178: }

182: PetscErrorCode STSetUp_Cayley(ST st)
183: {
185:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

188:   MatDestroy(&st->mat);

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

193:   if (!ctx->nu_set) { ctx->nu = st->sigma; }
194:   if (ctx->nu == 0.0 && st->sigma == 0.0) SETERRQ(((PetscObject)st)->comm,1,"Values of shift and antishift cannot be zero simultaneously");

196:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
197:   switch (st->shift_matrix) {
198:   case ST_MATMODE_INPLACE:
199:     st->mat = PETSC_NULL;
200:     if (st->sigma != 0.0) {
201:       if (st->B) {
202:         MatAXPY(st->A,-st->sigma,st->B,st->str);
203:       } else {
204:         MatShift(st->A,-st->sigma);
205:       }
206:     }
207:     /* TODO: should keep the Hermitian flag of st->A and restore at the end */
208:     STMatSetHermitian(st,st->A);
209:     KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
210:     break;
211:   case ST_MATMODE_SHELL:
212:     STMatShellCreate(st,&st->mat);
213:     STMatSetHermitian(st,st->mat);
214:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
215:     break;
216:   default:
217:     MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
218:     if (st->sigma != 0.0) {
219:       if (st->B) {
220:         MatAXPY(st->mat,-st->sigma,st->B,st->str);
221:       } else {
222:         MatShift(st->mat,-st->sigma);
223:       }
224:     }
225:     STMatSetHermitian(st,st->mat);
226:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
227:   }
228:   if (st->B) {
229:     VecDestroy(&ctx->w2);
230:     MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
231:   }
232:   KSPSetUp(st->ksp);
233:   return(0);
234: }

238: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
239: {
241:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
242:   MatStructure   flg;

245:   if (!ctx->nu_set) { ctx->nu = newshift; }
246:   if (ctx->nu == 0.0 && newshift == 0.0) SETERRQ(((PetscObject)st)->comm,1,"Values of shift and antishift cannot be zero simultaneously");

248:   /* Nothing to be done if STSetUp has not been called yet */
249:   if (!st->setupcalled) return(0);

251:   /* Check if the new KSP matrix has the same zero structure */
252:   if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
253:     flg = DIFFERENT_NONZERO_PATTERN;
254:   } else {
255:     flg = SAME_NONZERO_PATTERN;
256:   }

258:   switch (st->shift_matrix) {
259:   case ST_MATMODE_INPLACE:
260:     /* Undo previous operations */
261:     if (st->sigma != 0.0) {
262:       if (st->B) {
263:         MatAXPY(st->A,st->sigma,st->B,st->str);
264:       } else {
265:         MatShift(st->A,st->sigma);
266:       }
267:     }
268:     /* Apply new shift */
269:     if (newshift != 0.0) {
270:       if (st->B) {
271:         MatAXPY(st->A,-newshift,st->B,st->str);
272:       } else {
273:         MatShift(st->A,-newshift);
274:       }
275:     }
276:     STMatSetHermitian(st,st->A);
277:     KSPSetOperators(st->ksp,st->A,st->A,flg);
278:     break;
279:   case ST_MATMODE_SHELL:
280:     STMatSetHermitian(st,st->mat);
281:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
282:     break;
283:   default:
284:     MatCopy(st->A,st->mat,DIFFERENT_NONZERO_PATTERN);
285:     if (newshift != 0.0) {
286:       if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
287:       else { MatShift(st->mat,-newshift); }
288:     }
289:     STMatSetHermitian(st,st->mat);
290:     KSPSetOperators(st->ksp,st->mat,st->mat,flg);
291:   }
292:   st->sigma = newshift;
293:   KSPSetUp(st->ksp);
294:   return(0);
295: }

299: PetscErrorCode STSetFromOptions_Cayley(ST st)
300: {
302:   PetscScalar    nu;
303:   PetscBool      flg;
304:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
305:   PC             pc;
306:   const PCType   pctype;
307:   const KSPType  ksptype;

310:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
311:   KSPGetPC(st->ksp,&pc);
312:   KSPGetType(st->ksp,&ksptype);
313:   PCGetType(pc,&pctype);
314:   if (!pctype && !ksptype) {
315:     if (st->shift_matrix == ST_MATMODE_SHELL) {
316:       /* in shell mode use GMRES with Jacobi as the default */
317:       KSPSetType(st->ksp,KSPGMRES);
318:       PCSetType(pc,PCJACOBI);
319:     } else {
320:       /* use direct solver as default */
321:       KSPSetType(st->ksp,KSPPREONLY);
322:       PCSetType(pc,PCREDUNDANT);
323:     }
324:   }

326:   PetscOptionsHead("ST Cayley Options");
327:   PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
328:   if (flg) {
329:     STCayleySetAntishift(st,nu);
330:   }
331:   PetscOptionsTail();
332:   return(0);
333: }

335: EXTERN_C_BEGIN
338: PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
339: {
340:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

343:   ctx->nu     = newshift;
344:   ctx->nu_set = PETSC_TRUE;
345:   return(0);
346: }
347: EXTERN_C_END

351: /*@
352:    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
353:    spectral transformation.

355:    Logically Collective on ST

357:    Input Parameters:
358: +  st  - the spectral transformation context
359: -  nu  - the anti-shift

361:    Options Database Key:
362: .  -st_cayley_antishift - Sets the value of the anti-shift

364:    Level: intermediate

366:    Note:
367:    In the generalized Cayley transform, the operator can be expressed as
368:    OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
369:    Use STSetShift() for setting sigma.

371: .seealso: STSetShift(), STCayleyGetAntishift()
372: @*/
373: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
374: {

380:   PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
381:   return(0);
382: }
383: EXTERN_C_BEGIN
386: PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
387: {
388:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

391:   *nu = ctx->nu;
392:   return(0);
393: }
394: EXTERN_C_END

398: /*@
399:    STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
400:    spectral transformation.

402:    Not Collective

404:    Input Parameter:
405: .  st  - the spectral transformation context

407:    Output Parameter:
408: .  nu  - the anti-shift

410:    Level: intermediate

412: .seealso: STGetShift(), STCayleySetAntishift()
413: @*/
414: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
415: {

421:   PetscTryMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
422:   return(0);
423: }

427: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
428: {
430:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

433: #if !defined(PETSC_USE_COMPLEX)
434:   PetscViewerASCIIPrintf(viewer,"  Cayley: antishift: %G\n",ctx->nu);
435: #else
436:   PetscViewerASCIIPrintf(viewer,"  Cayley: antishift: %G+%G i\n",PetscRealPart(ctx->nu),PetscImaginaryPart(ctx->nu));
437: #endif
438:   return(0);
439: }

443: PetscErrorCode STReset_Cayley(ST st)
444: {
446:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

449:   VecDestroy(&ctx->w2);
450:   return(0);
451: }

455: PetscErrorCode STDestroy_Cayley(ST st)
456: {

460:   PetscFree(st->data);
461:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","",PETSC_NULL);
462:   return(0);
463: }

465: EXTERN_C_BEGIN
468: PetscErrorCode STCreate_Cayley(ST st)
469: {

473:   PetscNewLog(st,ST_CAYLEY,&st->data);
474:   st->ops->apply           = STApply_Cayley;
475:   st->ops->getbilinearform = STGetBilinearForm_Cayley;
476:   st->ops->applytrans      = STApplyTranspose_Cayley;
477:   st->ops->postsolve       = STPostSolve_Cayley;
478:   st->ops->backtr          = STBackTransform_Cayley;
479:   st->ops->setfromoptions  = STSetFromOptions_Cayley;
480:   st->ops->setup           = STSetUp_Cayley;
481:   st->ops->setshift        = STSetShift_Cayley;
482:   st->ops->destroy         = STDestroy_Cayley;
483:   st->ops->reset           = STReset_Cayley;
484:   st->ops->view            = STView_Cayley;
485:   st->ops->checknullspace  = STCheckNullSpace_Default;
486:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","STCayleySetAntishift_Cayley",STCayleySetAntishift_Cayley);
487:   return(0);
488: }
489: EXTERN_C_END