Actual source code: cayley.c

slepc-3.13.4 2020-09-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Implements the Cayley spectral transform
 12: */

 14:  #include <slepc/private/stimpl.h>

 16: typedef struct {
 17:   PetscScalar nu;
 18:   PetscBool   nu_set;
 19: } ST_CAYLEY;

 21: static PetscErrorCode MatMult_Cayley(Mat B,Vec x,Vec y)
 22: {
 24:   ST             st;
 25:   ST_CAYLEY      *ctx;
 26:   PetscScalar    nu;

 29:   MatShellGetContext(B,(void**)&st);
 30:   ctx = (ST_CAYLEY*)st->data;
 31:   nu = ctx->nu;

 33:   if (st->matmode == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };

 35:   if (st->nmat>1) {
 36:     /* generalized eigenproblem: y = (A + tB)x */
 37:     MatMult(st->A[0],x,y);
 38:     MatMult(st->A[1],x,st->work[1]);
 39:     VecAXPY(y,nu,st->work[1]);
 40:   } else {
 41:     /* standard eigenproblem: y = (A + tI)x */
 42:     MatMult(st->A[0],x,y);
 43:     VecAXPY(y,nu,x);
 44:   }
 45:   return(0);
 46: }

 48: static PetscErrorCode MatMultTranspose_Cayley(Mat B,Vec x,Vec y)
 49: {
 51:   ST             st;
 52:   ST_CAYLEY      *ctx;
 53:   PetscScalar    nu;

 56:   MatShellGetContext(B,(void**)&st);
 57:   ctx = (ST_CAYLEY*)st->data;
 58:   nu = ctx->nu;

 60:   if (st->matmode == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
 61:   nu = PetscConj(nu);

 63:   if (st->nmat>1) {
 64:     /* generalized eigenproblem: y = (A + tB)x */
 65:     MatMultTranspose(st->A[0],x,y);
 66:     MatMultTranspose(st->A[1],x,st->work[1]);
 67:     VecAXPY(y,nu,st->work[1]);
 68:   } else {
 69:     /* standard eigenproblem: y = (A + tI)x */
 70:     MatMultTranspose(st->A[0],x,y);
 71:     VecAXPY(y,nu,x);
 72:   }
 73:   return(0);
 74: }

 76: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
 77: {

 81:   STSetUp(st);
 82:   *B = st->T[0];
 83:   PetscObjectReference((PetscObject)*B);
 84:   return(0);
 85: }

 87: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 88: {
 89:   ST_CAYLEY   *ctx = (ST_CAYLEY*)st->data;
 90:   PetscInt    j;
 91: #if !defined(PETSC_USE_COMPLEX)
 92:   PetscScalar t,i,r;
 93: #endif

 96: #if !defined(PETSC_USE_COMPLEX)
 97:   for (j=0;j<n;j++) {
 98:     if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
 99:     else {
100:       r = eigr[j];
101:       i = eigi[j];
102:       r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
103:       i = - st->sigma * i - ctx->nu * i;
104:       t = i * i + r * (r - 2.0) + 1.0;
105:       eigr[j] = r / t;
106:       eigi[j] = i / t;
107:     }
108:   }
109: #else
110:   for (j=0;j<n;j++) {
111:     eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
112:   }
113: #endif
114:   return(0);
115: }

117: PetscErrorCode STPostSolve_Cayley(ST st)
118: {

122:   if (st->matmode == ST_MATMODE_INPLACE) {
123:     if (st->nmat>1) {
124:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
125:     } else {
126:       MatShift(st->A[0],st->sigma);
127:     }
128:     st->Astate[0] = ((PetscObject)st->A[0])->state;
129:     st->state   = ST_STATE_INITIAL;
130:     st->opready = PETSC_FALSE;
131:   }
132:   return(0);
133: }

135: /*
136:    Operator (cayley):
137:                Op                  P         M
138:    if nmat=1:  (A-sI)^-1 (A+tI)    A-sI      A+tI
139:    if nmat=2:  (A-sB)^-1 (A+tB)    A-sB      A+tI
140: */
141: PetscErrorCode STComputeOperator_Cayley(ST st)
142: {
144:   PetscInt       n,m;
145:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

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

151:   if (!ctx->nu_set) ctx->nu = st->sigma;
152:   if (ctx->nu == 0.0 && st->sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)st),1,"Values of shift and antishift cannot be zero simultaneously");
153:   if (ctx->nu == -st->sigma) SETERRQ(PetscObjectComm((PetscObject)st),1,"It is not allowed to set the antishift equal to minus the shift (the target)");

155:   /* T[0] = A+nu*B */
156:   if (st->matmode==ST_MATMODE_INPLACE) {
157:     MatGetLocalSize(st->A[0],&n,&m);
158:     MatCreateShell(PetscObjectComm((PetscObject)st),n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,&st->T[0]);
159:     MatShellSetOperation(st->T[0],MATOP_MULT,(void(*)(void))MatMult_Cayley);
160:     MatShellSetOperation(st->T[0],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Cayley);
161:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->T[0]);
162:   } else {
163:     STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[0]);
164:   }
165:   st->M = st->T[0];

167:   /* T[1] = A-sigma*B */
168:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[1]);
169:   PetscObjectReference((PetscObject)st->T[1]);
170:   MatDestroy(&st->P);
171:   st->P = st->T[1];
172:   return(0);
173: }

175: PetscErrorCode STSetUp_Cayley(ST st)
176: {

180:   if (st->nmat>2) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"Cayley transform cannot be used in polynomial eigenproblems");
181:   STSetWorkVecs(st,2);
182:   KSPSetUp(st->ksp);
183:   return(0);
184: }

186: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
187: {
189:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

192:   if (newshift==0.0 && (!ctx->nu_set || (ctx->nu_set && ctx->nu==0.0))) SETERRQ(PetscObjectComm((PetscObject)st),1,"Values of shift and antishift cannot be zero simultaneously");
193:   if (ctx->nu == -newshift) SETERRQ(PetscObjectComm((PetscObject)st),1,"It is not allowed to set the shift equal to minus the antishift");

195:   if (!ctx->nu_set) {
196:     if (st->matmode!=ST_MATMODE_INPLACE) {
197:       STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
198:     }
199:     ctx->nu = newshift;
200:   }
201:   STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,&st->T[1]);
202:   if (st->P!=st->T[1]) {
203:     PetscObjectReference((PetscObject)st->T[1]);
204:     MatDestroy(&st->P);
205:     st->P = st->T[1];
206:   }
207:   KSPSetOperators(st->ksp,st->P,st->P);
208:   KSPSetUp(st->ksp);
209:   return(0);
210: }

212: PetscErrorCode STSetFromOptions_Cayley(PetscOptionItems *PetscOptionsObject,ST st)
213: {
215:   PetscScalar    nu;
216:   PetscBool      flg;
217:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

220:   PetscOptionsHead(PetscOptionsObject,"ST Cayley Options");

222:     PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
223:     if (flg) { STCayleySetAntishift(st,nu); }

225:   PetscOptionsTail();
226:   return(0);
227: }

229: static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
230: {
232:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

235:   if (ctx->nu != newshift) {
236:     STCheckNotSeized(st,1);
237:     if (st->state && st->matmode!=ST_MATMODE_INPLACE) {
238:       STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
239:     }
240:     ctx->nu = newshift;
241:   }
242:   ctx->nu_set = PETSC_TRUE;
243:   return(0);
244: }

246: /*@
247:    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
248:    spectral transformation.

250:    Logically Collective on st

252:    Input Parameters:
253: +  st  - the spectral transformation context
254: -  nu  - the anti-shift

256:    Options Database Key:
257: .  -st_cayley_antishift - Sets the value of the anti-shift

259:    Level: intermediate

261:    Note:
262:    In the generalized Cayley transform, the operator can be expressed as
263:    OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
264:    Use STSetShift() for setting sigma. The value nu=-sigma is not allowed.

266: .seealso: STSetShift(), STCayleyGetAntishift()
267: @*/
268: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
269: {

275:   PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
276:   return(0);
277: }

279: static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
280: {
281:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

284:   *nu = ctx->nu;
285:   return(0);
286: }

288: /*@
289:    STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
290:    spectral transformation.

292:    Not Collective

294:    Input Parameter:
295: .  st  - the spectral transformation context

297:    Output Parameter:
298: .  nu  - the anti-shift

300:    Level: intermediate

302: .seealso: STGetShift(), STCayleySetAntishift()
303: @*/
304: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
305: {

311:   PetscUseMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
312:   return(0);
313: }

315: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
316: {
318:   char           str[50];
319:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
320:   PetscBool      isascii;

323:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
324:   if (isascii) {
325:     SlepcSNPrintfScalar(str,50,ctx->nu,PETSC_FALSE);
326:     PetscViewerASCIIPrintf(viewer,"  antishift: %s\n",str);
327:   }
328:   return(0);
329: }

331: PetscErrorCode STDestroy_Cayley(ST st)
332: {

336:   PetscFree(st->data);
337:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL);
338:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL);
339:   return(0);
340: }

342: SLEPC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
343: {
345:   ST_CAYLEY      *ctx;

348:   PetscNewLog(st,&ctx);
349:   st->data = (void*)ctx;

351:   st->usesksp = PETSC_TRUE;

353:   st->ops->apply           = STApply_Generic;
354:   st->ops->applytrans      = STApplyTranspose_Generic;
355:   st->ops->backtransform   = STBackTransform_Cayley;
356:   st->ops->setshift        = STSetShift_Cayley;
357:   st->ops->getbilinearform = STGetBilinearForm_Cayley;
358:   st->ops->setup           = STSetUp_Cayley;
359:   st->ops->computeoperator = STComputeOperator_Cayley;
360:   st->ops->setfromoptions  = STSetFromOptions_Cayley;
361:   st->ops->postsolve       = STPostSolve_Cayley;
362:   st->ops->destroy         = STDestroy_Cayley;
363:   st->ops->view            = STView_Cayley;
364:   st->ops->checknullspace  = STCheckNullSpace_Default;
365:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;

367:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley);
368:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley);
369:   return(0);
370: }