Actual source code: cayley.c

slepc-3.18.2 2023-01-26
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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: {
 23:   ST             st;
 24:   ST_CAYLEY      *ctx;
 25:   PetscScalar    nu;

 27:   MatShellGetContext(B,&st);
 28:   ctx = (ST_CAYLEY*)st->data;
 29:   nu = ctx->nu;

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

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

 46: static PetscErrorCode MatMultTranspose_Cayley(Mat B,Vec x,Vec y)
 47: {
 48:   ST             st;
 49:   ST_CAYLEY      *ctx;
 50:   PetscScalar    nu;

 52:   MatShellGetContext(B,&st);
 53:   ctx = (ST_CAYLEY*)st->data;
 54:   nu = ctx->nu;

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

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

 72: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
 73: {
 74:   STSetUp(st);
 75:   *B = st->T[0];
 76:   PetscObjectReference((PetscObject)*B);
 77:   return 0;
 78: }

 80: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
 81: {
 82:   ST_CAYLEY   *ctx = (ST_CAYLEY*)st->data;
 83:   PetscInt    j;
 84: #if !defined(PETSC_USE_COMPLEX)
 85:   PetscScalar t,i,r;
 86: #endif

 88: #if !defined(PETSC_USE_COMPLEX)
 89:   for (j=0;j<n;j++) {
 90:     if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
 91:     else {
 92:       r = eigr[j];
 93:       i = eigi[j];
 94:       r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
 95:       i = - st->sigma * i - ctx->nu * i;
 96:       t = i * i + r * (r - 2.0) + 1.0;
 97:       eigr[j] = r / t;
 98:       eigi[j] = i / t;
 99:     }
100:   }
101: #else
102:   for (j=0;j<n;j++) {
103:     eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
104:   }
105: #endif
106:   return 0;
107: }

109: PetscErrorCode STPostSolve_Cayley(ST st)
110: {
111:   if (st->matmode == ST_MATMODE_INPLACE) {
112:     if (st->nmat>1) MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
113:     else MatShift(st->A[0],st->sigma);
114:     st->Astate[0] = ((PetscObject)st->A[0])->state;
115:     st->state   = ST_STATE_INITIAL;
116:     st->opready = PETSC_FALSE;
117:   }
118:   return 0;
119: }

121: /*
122:    Operator (cayley):
123:                Op                  P         M
124:    if nmat=1:  (A-sI)^-1 (A+tI)    A-sI      A+tI
125:    if nmat=2:  (A-sB)^-1 (A+tB)    A-sB      A+tI
126: */
127: PetscErrorCode STComputeOperator_Cayley(ST st)
128: {
129:   PetscInt       n,m;
130:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

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

135:   if (!ctx->nu_set) ctx->nu = st->sigma;

139:   /* T[0] = A+nu*B */
140:   if (st->matmode==ST_MATMODE_INPLACE) {
141:     MatGetLocalSize(st->A[0],&n,&m);
142:     MatCreateShell(PetscObjectComm((PetscObject)st),n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,&st->T[0]);
143:     MatShellSetOperation(st->T[0],MATOP_MULT,(void(*)(void))MatMult_Cayley);
144:     MatShellSetOperation(st->T[0],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Cayley);
145:   } else STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[0]);
146:   st->M = st->T[0];

148:   /* T[1] = A-sigma*B */
149:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]);
150:   PetscObjectReference((PetscObject)st->T[1]);
151:   MatDestroy(&st->P);
152:   st->P = st->T[1];
153:   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
154:     STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat);
155:   }
156:   return 0;
157: }

159: PetscErrorCode STSetUp_Cayley(ST st)
160: {
162:   STSetWorkVecs(st,2);
163:   KSPSetUp(st->ksp);
164:   return 0;
165: }

167: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
168: {
169:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;


174:   if (!ctx->nu_set) {
175:     if (st->matmode!=ST_MATMODE_INPLACE) STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]);
176:     ctx->nu = newshift;
177:   }
178:   STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[1]);
179:   if (st->P!=st->T[1]) {
180:     PetscObjectReference((PetscObject)st->T[1]);
181:     MatDestroy(&st->P);
182:     st->P = st->T[1];
183:   }
184:   if (st->Psplit) {  /* build custom preconditioner from the split matrices */
185:     STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_TRUE,&st->Pmat);
186:   }
187:   ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
188:   KSPSetUp(st->ksp);
189:   return 0;
190: }

192: PetscErrorCode STSetFromOptions_Cayley(ST st,PetscOptionItems *PetscOptionsObject)
193: {
194:   PetscScalar    nu;
195:   PetscBool      flg;
196:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

198:   PetscOptionsHeadBegin(PetscOptionsObject,"ST Cayley Options");

200:     PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
201:     if (flg) STCayleySetAntishift(st,nu);

203:   PetscOptionsHeadEnd();
204:   return 0;
205: }

207: static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
208: {
209:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

211:   if (ctx->nu != newshift) {
212:     STCheckNotSeized(st,1);
213:     if (st->state && st->matmode!=ST_MATMODE_INPLACE) STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]);
214:     ctx->nu = newshift;
215:   }
216:   ctx->nu_set = PETSC_TRUE;
217:   return 0;
218: }

220: /*@
221:    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
222:    spectral transformation.

224:    Logically Collective on st

226:    Input Parameters:
227: +  st  - the spectral transformation context
228: -  nu  - the anti-shift

230:    Options Database Key:
231: .  -st_cayley_antishift - Sets the value of the anti-shift

233:    Level: intermediate

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

240: .seealso: STSetShift(), STCayleyGetAntishift()
241: @*/
242: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
243: {
246:   PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
247:   return 0;
248: }

250: static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
251: {
252:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

254:   *nu = ctx->nu;
255:   return 0;
256: }

258: /*@
259:    STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
260:    spectral transformation.

262:    Not Collective

264:    Input Parameter:
265: .  st  - the spectral transformation context

267:    Output Parameter:
268: .  nu  - the anti-shift

270:    Level: intermediate

272: .seealso: STGetShift(), STCayleySetAntishift()
273: @*/
274: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
275: {
278:   PetscUseMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
279:   return 0;
280: }

282: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
283: {
284:   char           str[50];
285:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
286:   PetscBool      isascii;

288:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
289:   if (isascii) {
290:     SlepcSNPrintfScalar(str,sizeof(str),ctx->nu,PETSC_FALSE);
291:     PetscViewerASCIIPrintf(viewer,"  antishift: %s\n",str);
292:   }
293:   return 0;
294: }

296: PetscErrorCode STDestroy_Cayley(ST st)
297: {
298:   PetscFree(st->data);
299:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL);
300:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL);
301:   return 0;
302: }

304: SLEPC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
305: {
306:   ST_CAYLEY      *ctx;

308:   PetscNew(&ctx);
309:   st->data = (void*)ctx;

311:   st->usesksp = PETSC_TRUE;

313:   st->ops->apply           = STApply_Generic;
314:   st->ops->applytrans      = STApplyTranspose_Generic;
315:   st->ops->backtransform   = STBackTransform_Cayley;
316:   st->ops->setshift        = STSetShift_Cayley;
317:   st->ops->getbilinearform = STGetBilinearForm_Cayley;
318:   st->ops->setup           = STSetUp_Cayley;
319:   st->ops->computeoperator = STComputeOperator_Cayley;
320:   st->ops->setfromoptions  = STSetFromOptions_Cayley;
321:   st->ops->postsolve       = STPostSolve_Cayley;
322:   st->ops->destroy         = STDestroy_Cayley;
323:   st->ops->view            = STView_Cayley;
324:   st->ops->checknullspace  = STCheckNullSpace_Default;
325:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;

327:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley);
328:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley);
329:   return 0;
330: }