Actual source code: cayley.c
slepc-3.17.1 2022-04-11
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: PetscFunctionReturn(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: PetscFunctionReturn(0);
70: }
72: PetscErrorCode STGetBilinearForm_Cayley(ST st,Mat *B)
73: {
74: STSetUp(st);
75: *B = st->T[0];
76: PetscObjectReference((PetscObject)*B);
77: PetscFunctionReturn(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: PetscFunctionReturn(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: PetscFunctionReturn(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: PetscLogObjectParent((PetscObject)st,(PetscObject)st->T[0]);
146: } else STMatMAXPY_Private(st,ctx->nu,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[0]);
147: st->M = st->T[0];
149: /* T[1] = A-sigma*B */
150: STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PetscNot(st->state==ST_STATE_UPDATED),PETSC_FALSE,&st->T[1]);
151: PetscObjectReference((PetscObject)st->T[1]);
152: MatDestroy(&st->P);
153: st->P = st->T[1];
154: if (st->Psplit) { /* build custom preconditioner from the split matrices */
155: STMatMAXPY_Private(st,-st->sigma,0.0,0,NULL,PETSC_TRUE,PETSC_TRUE,&st->Pmat);
156: }
157: PetscFunctionReturn(0);
158: }
160: PetscErrorCode STSetUp_Cayley(ST st)
161: {
163: STSetWorkVecs(st,2);
164: KSPSetUp(st->ksp);
165: PetscFunctionReturn(0);
166: }
168: PetscErrorCode STSetShift_Cayley(ST st,PetscScalar newshift)
169: {
170: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
175: if (!ctx->nu_set) {
176: if (st->matmode!=ST_MATMODE_INPLACE) STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]);
177: ctx->nu = newshift;
178: }
179: STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[1]);
180: if (st->P!=st->T[1]) {
181: PetscObjectReference((PetscObject)st->T[1]);
182: MatDestroy(&st->P);
183: st->P = st->T[1];
184: }
185: if (st->Psplit) { /* build custom preconditioner from the split matrices */
186: STMatMAXPY_Private(st,-newshift,-st->sigma,0,NULL,PETSC_FALSE,PETSC_TRUE,&st->Pmat);
187: }
188: ST_KSPSetOperators(st,st->P,st->Pmat?st->Pmat:st->P);
189: KSPSetUp(st->ksp);
190: PetscFunctionReturn(0);
191: }
193: PetscErrorCode STSetFromOptions_Cayley(PetscOptionItems *PetscOptionsObject,ST st)
194: {
195: PetscScalar nu;
196: PetscBool flg;
197: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
199: PetscOptionsHead(PetscOptionsObject,"ST Cayley Options");
201: PetscOptionsScalar("-st_cayley_antishift","Value of the antishift","STCayleySetAntishift",ctx->nu,&nu,&flg);
202: if (flg) STCayleySetAntishift(st,nu);
204: PetscOptionsTail();
205: PetscFunctionReturn(0);
206: }
208: static PetscErrorCode STCayleySetAntishift_Cayley(ST st,PetscScalar newshift)
209: {
210: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
212: if (ctx->nu != newshift) {
213: STCheckNotSeized(st,1);
214: if (st->state && st->matmode!=ST_MATMODE_INPLACE) STMatMAXPY_Private(st,newshift,ctx->nu,0,NULL,PETSC_FALSE,PETSC_FALSE,&st->T[0]);
215: ctx->nu = newshift;
216: }
217: ctx->nu_set = PETSC_TRUE;
218: PetscFunctionReturn(0);
219: }
221: /*@
222: STCayleySetAntishift - Sets the value of the anti-shift for the Cayley
223: spectral transformation.
225: Logically Collective on st
227: Input Parameters:
228: + st - the spectral transformation context
229: - nu - the anti-shift
231: Options Database Key:
232: . -st_cayley_antishift - Sets the value of the anti-shift
234: Level: intermediate
236: Note:
237: In the generalized Cayley transform, the operator can be expressed as
238: OP = inv(A - sigma B)*(A + nu B). This function sets the value of nu.
239: Use STSetShift() for setting sigma. The value nu=-sigma is not allowed.
241: .seealso: STSetShift(), STCayleyGetAntishift()
242: @*/
243: PetscErrorCode STCayleySetAntishift(ST st,PetscScalar nu)
244: {
247: PetscTryMethod(st,"STCayleySetAntishift_C",(ST,PetscScalar),(st,nu));
248: PetscFunctionReturn(0);
249: }
251: static PetscErrorCode STCayleyGetAntishift_Cayley(ST st,PetscScalar *nu)
252: {
253: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
255: *nu = ctx->nu;
256: PetscFunctionReturn(0);
257: }
259: /*@
260: STCayleyGetAntishift - Gets the value of the anti-shift used in the Cayley
261: spectral transformation.
263: Not Collective
265: Input Parameter:
266: . st - the spectral transformation context
268: Output Parameter:
269: . nu - the anti-shift
271: Level: intermediate
273: .seealso: STGetShift(), STCayleySetAntishift()
274: @*/
275: PetscErrorCode STCayleyGetAntishift(ST st,PetscScalar *nu)
276: {
279: PetscUseMethod(st,"STCayleyGetAntishift_C",(ST,PetscScalar*),(st,nu));
280: PetscFunctionReturn(0);
281: }
283: PetscErrorCode STView_Cayley(ST st,PetscViewer viewer)
284: {
285: char str[50];
286: ST_CAYLEY *ctx = (ST_CAYLEY*)st->data;
287: PetscBool isascii;
289: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
290: if (isascii) {
291: SlepcSNPrintfScalar(str,sizeof(str),ctx->nu,PETSC_FALSE);
292: PetscViewerASCIIPrintf(viewer," antishift: %s\n",str);
293: }
294: PetscFunctionReturn(0);
295: }
297: PetscErrorCode STDestroy_Cayley(ST st)
298: {
299: PetscFree(st->data);
300: PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",NULL);
301: PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",NULL);
302: PetscFunctionReturn(0);
303: }
305: SLEPC_EXTERN PetscErrorCode STCreate_Cayley(ST st)
306: {
307: ST_CAYLEY *ctx;
309: PetscNewLog(st,&ctx);
310: st->data = (void*)ctx;
312: st->usesksp = PETSC_TRUE;
314: st->ops->apply = STApply_Generic;
315: st->ops->applytrans = STApplyTranspose_Generic;
316: st->ops->backtransform = STBackTransform_Cayley;
317: st->ops->setshift = STSetShift_Cayley;
318: st->ops->getbilinearform = STGetBilinearForm_Cayley;
319: st->ops->setup = STSetUp_Cayley;
320: st->ops->computeoperator = STComputeOperator_Cayley;
321: st->ops->setfromoptions = STSetFromOptions_Cayley;
322: st->ops->postsolve = STPostSolve_Cayley;
323: st->ops->destroy = STDestroy_Cayley;
324: st->ops->view = STView_Cayley;
325: st->ops->checknullspace = STCheckNullSpace_Default;
326: st->ops->setdefaultksp = STSetDefaultKSP_Default;
328: PetscObjectComposeFunction((PetscObject)st,"STCayleySetAntishift_C",STCayleySetAntishift_Cayley);
329: PetscObjectComposeFunction((PetscObject)st,"STCayleyGetAntishift_C",STCayleyGetAntishift_Cayley);
330: PetscFunctionReturn(0);
331: }