Actual source code: cayley.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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>          /*I "slepcst.h" I*/

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

 21: PetscErrorCode STApply_Cayley(ST st,Vec x,Vec y)
 22: {

 26:   /* standard eigenproblem: y = (A - sI)^-1 (A + tI)x */
 27:   /* generalized eigenproblem: y = (A - sB)^-1 (A + tB)x */
 28:   MatMult(st->T[0],x,st->work[0]);
 29:   STMatSolve(st,st->work[0],y);
 30:   return(0);
 31: }

 33: PetscErrorCode STApplyTranspose_Cayley(ST st,Vec x,Vec y)
 34: {

 38:   /* standard eigenproblem: y =  (A + tI)^T (A - sI)^-T x */
 39:   /* generalized eigenproblem: y = (A + tB)^T (A - sB)^-T x */
 40:   STMatSolveTranspose(st,x,st->work[0]);
 41:   MatMultTranspose(st->T[0],st->work[0],y);
 42:   return(0);
 43: }

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

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

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

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

 72: static PetscErrorCode MatMultTranspose_Cayley(Mat B,Vec x,Vec y)
 73: {
 75:   ST             st;
 76:   ST_CAYLEY      *ctx;
 77:   PetscScalar    nu;

 80:   MatShellGetContext(B,(void**)&st);
 81:   ctx = (ST_CAYLEY*)st->data;
 82:   nu = ctx->nu;

 84:   if (st->shift_matrix == ST_MATMODE_INPLACE) { nu = nu + st->sigma; };
 85:   nu = PetscConj(nu);

 87:   if (st->nmat>1) {
 88:     /* generalized eigenproblem: y = (A + tB)x */
 89:     MatMultTranspose(st->A[0],x,y);
 90:     MatMultTranspose(st->A[1],x,st->work[1]);
 91:     VecAXPY(y,nu,st->work[1]);
 92:   } else {
 93:     /* standard eigenproblem: y = (A + tI)x */
 94:     MatMultTranspose(st->A[0],x,y);
 95:     VecAXPY(y,nu,x);
 96:   }
 97:   return(0);
 98: }

100: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
101: {

105:   STSetUp(st);
106:   *B = st->T[0];
107:   PetscObjectReference((PetscObject)*B);
108:   return(0);
109: }

111: PetscErrorCode STBackTransform_Cayley(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
112: {
113:   ST_CAYLEY   *ctx = (ST_CAYLEY*)st->data;
114:   PetscInt    j;
115: #if !defined(PETSC_USE_COMPLEX)
116:   PetscScalar t,i,r;
117: #endif

120: #if !defined(PETSC_USE_COMPLEX)
121:   for (j=0;j<n;j++) {
122:     if (eigi[j] == 0.0) eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
123:     else {
124:       r = eigr[j];
125:       i = eigi[j];
126:       r = st->sigma * (r * r + i * i - r) + ctx->nu * (r - 1);
127:       i = - st->sigma * i - ctx->nu * i;
128:       t = i * i + r * (r - 2.0) + 1.0;
129:       eigr[j] = r / t;
130:       eigi[j] = i / t;
131:     }
132:   }
133: #else
134:   for (j=0;j<n;j++) {
135:     eigr[j] = (ctx->nu + eigr[j] * st->sigma) / (eigr[j] - 1.0);
136:   }
137: #endif
138:   return(0);
139: }

141: PetscErrorCode STPostSolve_Cayley(ST st)
142: {

146:   if (st->shift_matrix == ST_MATMODE_INPLACE) {
147:     if (st->nmat>1) {
148:       MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
149:     } else {
150:       MatShift(st->A[0],st->sigma);
151:     }
152:     st->Astate[0] = ((PetscObject)st->A[0])->state;
153:     st->state = ST_STATE_INITIAL;
154:   }
155:   return(0);
156: }

158: PetscErrorCode STSetUp_Cayley(ST st)
159: {
161:   PetscInt       n,m;
162:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

165:   STSetWorkVecs(st,2);

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

170:   if (!ctx->nu_set) ctx->nu = st->sigma;
171:   if (ctx->nu == 0.0 && st->sigma == 0.0) SETERRQ(PetscObjectComm((PetscObject)st),1,"Values of shift and antishift cannot be zero simultaneously");

173:   /* T[0] = A+nu*B */
174:   if (st->shift_matrix==ST_MATMODE_INPLACE) {
175:     MatGetLocalSize(st->A[0],&n,&m);
176:     MatCreateShell(PetscObjectComm((PetscObject)st),n,m,PETSC_DETERMINE,PETSC_DETERMINE,st,&st->T[0]);
177:     MatShellSetOperation(st->T[0],MATOP_MULT,(void(*)(void))MatMult_Cayley);
178:     MatShellSetOperation(st->T[0],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_Cayley);
179:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->T[0]);
180:   } else {
181:     STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[0]);
182:   }

184:   /* T[1] = A-sigma*B */
185:   STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),&st->T[1]);
186:   PetscObjectReference((PetscObject)st->T[1]);
187:   MatDestroy(&st->P);
188:   st->P = st->T[1];
189:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
190:   STCheckFactorPackage(st);
191:   KSPSetOperators(st->ksp,st->P,st->P);
192:   KSPSetUp(st->ksp);
193:   return(0);
194: }

196: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
197: {
199:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

202:   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");

204:   if (!ctx->nu_set) {
205:     if (st->shift_matrix!=ST_MATMODE_INPLACE) {
206:       STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
207:     }
208:     ctx->nu = newshift;
209:   }
210:   STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,&st->T[1]);
211:   if (st->P!=st->T[1]) {
212:     PetscObjectReference((PetscObject)st->T[1]);
213:     MatDestroy(&st->P);
214:     st->P = st->T[1];
215:   }
216:   KSPSetOperators(st->ksp,st->P,st->P);
217:   KSPSetUp(st->ksp);
218:   return(0);
219: }

221: PetscErrorCode STSetFromOptions_Cayley(PetscOptionItems *PetscOptionsObject,ST st)
222: {
224:   PetscScalar    nu;
225:   PetscBool      flg;
226:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;

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

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

234:   PetscOptionsTail();
235:   return(0);
236: }

238: static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
239: {
241:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

244:   if (st->state && st->shift_matrix!=ST_MATMODE_INPLACE) {
245:     STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,&st->T[0]);
246:   }
247:   ctx->nu     = newshift;
248:   ctx->nu_set = PETSC_TRUE;
249:   return(0);
250: }

252: /*@
253:    STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
254:    spectral transformation.

256:    Logically Collective on ST

258:    Input Parameters:
259: +  st  - the spectral transformation context
260: -  nu  - the anti-shift

262:    Options Database Key:
263: .  -st_cayley_antishift - Sets the value of the anti-shift

265:    Level: intermediate

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

272: .seealso: STSetShift(), STCayleyGetAntishift()
273: @*/
274: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
275: {

281:   PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
282:   return(0);
283: }

285: static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
286: {
287:   ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;

290:   *nu = ctx->nu;
291:   return(0);
292: }

294: /*@
295:    STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
296:    spectral transformation.

298:    Not Collective

300:    Input Parameter:
301: .  st  - the spectral transformation context

303:    Output Parameter:
304: .  nu  - the anti-shift

306:    Level: intermediate

308: .seealso: STGetShift(), STCayleySetAntishift()
309: @*/
310: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
311: {

317:   PetscUseMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
318:   return(0);
319: }

321: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
322: {
324:   char           str[50];
325:   ST_CAYLEY      *ctx = (ST_CAYLEY*)st->data;
326:   PetscBool      isascii;

329:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
330:   if (isascii) {
331:     SlepcSNPrintfScalar(str,50,ctx->nu,PETSC_FALSE);
332:     PetscViewerASCIIPrintf(viewer,"  antishift: %s\n",str);
333:   }
334:   return(0);
335: }

337: PetscErrorCode STDestroy_Cayley(ST st)
338: {

342:   PetscFree(st->data);
343:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL);
344:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL);
345:   return(0);
346: }

348: SLEPC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
349: {
351:   ST_CAYLEY      *ctx;

354:   PetscNewLog(st,&ctx);
355:   st->data = (void*)ctx;

357:   st->ops->apply           = STApply_Cayley;
358:   st->ops->getbilinearform = STGetBilinearForm_Cayley;
359:   st->ops->applytrans      = STApplyTranspose_Cayley;
360:   st->ops->postsolve       = STPostSolve_Cayley;
361:   st->ops->backtransform   = STBackTransform_Cayley;
362:   st->ops->setfromoptions  = STSetFromOptions_Cayley;
363:   st->ops->setup           = STSetUp_Cayley;
364:   st->ops->setshift        = STSetShift_Cayley;
365:   st->ops->destroy         = STDestroy_Cayley;
366:   st->ops->view            = STView_Cayley;
367:   st->ops->checknullspace  = STCheckNullSpace_Default;
368:   st->ops->setdefaultksp   = STSetDefaultKSP_Default;
369:   PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley);
370:   PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley);
371:   return(0);
372: }