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

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

 13:  #include src/st/stimpl.h

 15: typedef struct {
 16:   PetscScalar tau;
 17:   PetscTruth  tau_set;
 18:   Vec         w2;
 19: } ST_CAYLEY;

 23: PetscErrorCode STApply_Cayley(ST st,Vec x,Vec y)
 24: {
 26:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;
 27:   PetscScalar    tau = ctx->tau;
 28: 
 30:   if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + st->sigma; };

 32:   if (st->B) {
 33:     /* generalized eigenproblem: y = (A - sB)^-1 (A + tB)x */
 34:     MatMult(st->A,x,st->w);
 35:     MatMult(st->B,x,ctx->w2);
 36:     VecAXPY(st->w,tau,ctx->w2);
 37:     STAssociatedKSPSolve(st,st->w,y);
 38:   }
 39:   else {
 40:     /* standard eigenproblem: y = (A - sI)^-1 (A + tI)x */
 41:     MatMult(st->A,x,st->w);
 42:     VecAXPY(st->w,tau,x);
 43:     STAssociatedKSPSolve(st,st->w,y);
 44:   }
 45:   return(0);
 46: }

 50: PetscErrorCode STApplyTranspose_Cayley(ST st,Vec x,Vec y)
 51: {
 53:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;
 54:   PetscScalar    tau = ctx->tau;
 55: 
 57:   if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + st->sigma; };

 59:   if (st->B) {
 60:     /* generalized eigenproblem: y = (A + tB)^T (A - sB)^-T x */
 61:     STAssociatedKSPSolve(st,x,st->w);
 62:     MatMult(st->A,st->w,y);
 63:     MatMult(st->B,st->w,ctx->w2);
 64:     VecAXPY(y,tau,ctx->w2);
 65:   }
 66:   else {
 67:     /* standard eigenproblem: y =  (A + tI)^T (A - sI)^-T x */
 68:     STAssociatedKSPSolve(st,x,st->w);
 69:     MatMult(st->A,st->w,y);
 70:     VecAXPY(y,tau,st->w);
 71:   }
 72:   return(0);
 73: }

 77: PetscErrorCode STBilinearMatMult_Cayley(Mat B,Vec x,Vec y)
 78: {
 80:   ST             st;
 81:   ST_CAYLEY      *ctx;
 82:   PetscScalar    tau;
 83: 
 85:   MatShellGetContext(B,(void**)&st);
 86:   ctx = (ST_CAYLEY *) st->data;
 87:   tau = ctx->tau;
 88: 
 89:   if (st->shift_matrix == STMATMODE_INPLACE) { tau = tau + st->sigma; };

 91:   if (st->B) {
 92:     /* generalized eigenproblem: y = (A + tB)x */
 93:     MatMult(st->A,x,y);
 94:     MatMult(st->B,x,ctx->w2);
 95:     VecAXPY(y,tau,ctx->w2);
 96:   }
 97:   else {
 98:     /* standard eigenproblem: y = (A + tI)x */
 99:     MatMult(st->A,x,y);
100:     VecAXPY(y,tau,x);
101:   }
102:   return(0);
103: }

107: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
108: {
110:   PetscInt       n,m;

113:   MatGetLocalSize(st->B,&n,&m);
114:   MatCreateShell(st->comm,n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,B);
115:   MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))STBilinearMatMult_Cayley);
116:   return(0);
117: }

121: PetscErrorCode STBackTransform_Cayley(ST st,PetscScalar *eigr,PetscScalar *eigi)
122: {
123:   ST_CAYLEY   *ctx = (ST_CAYLEY *) st->data;
124: #ifndef PETSC_USE_COMPLEX
125:   PetscScalar t,i,r;
129:   if (*eigi == 0.0) *eigr = (ctx->tau + *eigr * st->sigma) / (*eigr - 1.0);
130:   else {
131:     r = *eigr;
132:     i = *eigi;
133:     r = st->sigma * (r * r + i * i - r) + ctx->tau * (r - 1);
134:     i = - st->sigma * i - ctx->tau * i;
135:     t = i * i + r * (r - 2.0) + 1.0;
136:     *eigr = r / t;
137:     *eigi = i / t;
138:   }
139: #else
142:   *eigr = (ctx->tau + *eigr * st->sigma) / (*eigr - 1.0);
143: #endif
144:   return(0);
145: }

149: PetscErrorCode STPostSolve_Cayley(ST st)
150: {

154:   if (st->shift_matrix == STMATMODE_INPLACE) {
155:     if (st->B) {
156:       MatAXPY(st->A,st->sigma,st->B,st->str);
157:     } else {
158:       MatShift(st->A,st->sigma);
159:     }
160:     st->setupcalled = 0;
161:   }
162:   return(0);
163: }

167: PetscErrorCode STSetUp_Cayley(ST st)
168: {
170:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;


174:   if (st->mat) { MatDestroy(st->mat); }
175:   if (!ctx->tau_set) { ctx->tau = st->sigma; }
176:   if (ctx->tau == 0.0 &&  st->sigma == 0.0) {
177:     SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
178:   }

180:   switch (st->shift_matrix) {
181:   case STMATMODE_INPLACE:
182:     st->mat = PETSC_NULL;
183:     if (st->sigma != 0.0) {
184:       if (st->B) {
185:         MatAXPY(st->A,-st->sigma,st->B,st->str);
186:       } else {
187:         MatShift(st->A,-st->sigma);
188:       }
189:     }
190:     KSPSetOperators(st->ksp,st->A,st->A,DIFFERENT_NONZERO_PATTERN);
191:     break;
192:   case STMATMODE_SHELL:
193:     STMatShellCreate(st,&st->mat);
194:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
195:     break;
196:   default:
197:     MatDuplicate(st->A,MAT_COPY_VALUES,&st->mat);
198:     if (st->sigma != 0.0) {
199:       if (st->B) {
200:         MatAXPY(st->mat,-st->sigma,st->B,st->str);
201:       } else {
202:         MatShift(st->mat,-st->sigma);
203:       }
204:     }
205:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
206:   }
207:   if (st->B) {
208:    if (ctx->w2) { VecDestroy(ctx->w2); }
209:    MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
210:   }
211:   KSPSetUp(st->ksp);
212:   return(0);
213: }

217: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
218: {
220:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;
221:   MatStructure   flg;

224:   if (!ctx->tau_set) { ctx->tau = newshift; }
225:   if (ctx->tau == 0.0 &&  newshift == 0.0) {
226:     SETERRQ(1,"Values of shift and antishift cannot be zero simultaneously");
227:   }

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

232:   /* Check if the new KSP matrix has the same zero structure */
233:   if (st->B && st->str == DIFFERENT_NONZERO_PATTERN && (st->sigma == 0.0 || newshift == 0.0)) {
234:     flg = DIFFERENT_NONZERO_PATTERN;
235:   } else {
236:     flg = SAME_NONZERO_PATTERN;
237:   }

239:   switch (st->shift_matrix) {
240:   case STMATMODE_INPLACE:
241:     /* Undo previous operations */
242:     if (st->sigma != 0.0) {
243:       if (st->B) {
244:         MatAXPY(st->A,st->sigma,st->B,st->str);
245:       } else {
246:         MatShift(st->A,st->sigma);
247:       }
248:     }
249:     /* Apply new shift */
250:     if (newshift != 0.0) {
251:       if (st->B) {
252:         MatAXPY(st->A,-newshift,st->B,st->str);
253:       } else {
254:         MatShift(st->A,-newshift);
255:       }
256:     }
257:     KSPSetOperators(st->ksp,st->A,st->A,flg);
258:     break;
259:   case STMATMODE_SHELL:
260:     KSPSetOperators(st->ksp,st->mat,st->mat,DIFFERENT_NONZERO_PATTERN);
261:     break;
262:   default:
263:     MatCopy(st->A, st->mat,SUBSET_NONZERO_PATTERN);
264:     if (newshift != 0.0) {
265:       if (st->B) { MatAXPY(st->mat,-newshift,st->B,st->str); }
266:       else { MatShift(st->mat,-newshift); }
267:     }
268:     KSPSetOperators(st->ksp,st->mat,st->mat,flg);
269:   }
270:   st->sigma = newshift;
271:   KSPSetUp(st->ksp);
272:   return(0);
273: }

277: PetscErrorCode STSetFromOptions_Cayley(ST st)
278: {
280:   PetscScalar    tau;
281:   PetscTruth     flg;
282:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;

285:   PetscOptionsHead("ST Cayley Options");
286:   PetscOptionsScalar("-st_antishift","Value of the antishift","STSetAntishift",ctx->tau,&tau,&flg);
287:   if (flg) {
288:     STCayleySetAntishift(st,tau);
289:   }
290:   PetscOptionsTail();
291:   return(0);
292: }

297: PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
298: {
299:   ST_CAYLEY *ctx = (ST_CAYLEY *) st->data;

302:   ctx->tau = newshift;
303:   ctx->tau_set = PETSC_TRUE;
304:   return(0);
305: }

310: /*@
311:    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
312:    spectral transformation.

314:    Collective on ST

316:    Input Parameters:
317: +  st  - the spectral transformation context
318: -  tau - the anti-shift

320:    Options Database Key:
321: .  -st_antishift - Sets the value of the anti-shift

323:    Level: intermediate

325:    Note:
326:    In the generalized Cayley transform, the operator can be expressed as
327:    OP = inv(A - sigma B)*(A + tau B). This function sets the value of tau.
328:    Use STSetShift() for setting sigma.

330: .seealso: STSetShift()
331: @*/
332: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar tau)
333: {
334:   PetscErrorCode ierr, (*f)(ST,PetscScalar);

338:   PetscObjectQueryFunction((PetscObject)st,"STCayleySetAntishift_C",(void (**)(void))&f);
339:   if (f) {
340:     (*f)(st,tau);
341:   }
342:   return(0);
343: }

347: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
348: {
350:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;

353: #if !defined(PETSC_USE_COMPLEX)
354:   PetscViewerASCIIPrintf(viewer,"  antishift: %g\n",ctx->tau);
355: #else
356:   PetscViewerASCIIPrintf(viewer,"  antishift: %g+%g i\n",PetscRealPart(ctx->tau),PetscImaginaryPart(ctx->tau));
357: #endif
358:   STView_Default(st,viewer);
359:   return(0);
360: }

364: PetscErrorCode STDestroy_Cayley(ST st)
365: {
367:   ST_CAYLEY      *ctx = (ST_CAYLEY *) st->data;

370:   if (ctx->w2) { VecDestroy(ctx->w2); }
371:   PetscFree(ctx);
372:   return(0);
373: }

378: PetscErrorCode STCreate_Cayley(ST st)
379: {
381:   ST_CAYLEY      *ctx;

384:   PetscNew(ST_CAYLEY,&ctx);
385:   PetscLogObjectMemory(st,sizeof(ST_CAYLEY));
386:   st->data                 = (void *) ctx;

388:   st->ops->apply           = STApply_Cayley;
389:   st->ops->getbilinearform = STGetBilinearForm_Cayley;
390:   st->ops->applytrans      = STApplyTranspose_Cayley;
391:   st->ops->postsolve       = STPostSolve_Cayley;
392:   st->ops->backtr          = STBackTransform_Cayley;
393:   st->ops->setfromoptions  = STSetFromOptions_Cayley;
394:   st->ops->setup           = STSetUp_Cayley;
395:   st->ops->setshift        = STSetShift_Cayley;
396:   st->ops->destroy         = STDestroy_Cayley;
397:   st->ops->view            = STView_Cayley;
398: 
399:   st->checknullspace      = STCheckNullSpace_Default;

401:   ctx->tau                = 0.0;
402:   ctx->tau_set            = PETSC_FALSE;

404:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STCayleySetAntishift_C","STCayleySetAntishift_Cayley",
405:                     STCayleySetAntishift_Cayley);

407:   return(0);
408: }