Actual source code: fold.c
1: /*
2: Folding spectral transformation, applies (A + sigma I)^2 as operator, or
3: inv(B)(A + sigma I)^2 for generalized problems
4: */
5: #include src/st/stimpl.h
7: typedef struct {
8: PetscTruth left;
9: Vec w2;
10: } ST_FOLD;
14: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
15: {
17: ST_FOLD *ctx = (ST_FOLD *) st->data;
20: if (st->B) {
21: /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
22: MatMult(st->A,x,st->w);
23: STAssociatedKSPSolve(st,st->w,ctx->w2);
24: if (st->sigma != 0.0) {
25: VecAXPY(ctx->w2,-st->sigma,x);
26: }
27: MatMult(st->A,ctx->w2,st->w);
28: STAssociatedKSPSolve(st,st->w,y);
29: if (st->sigma != 0.0) {
30: VecAXPY(y,-st->sigma,ctx->w2);
31: }
32: } else {
33: /* standard eigenproblem: y = (A + sI)^2 x */
34: MatMult(st->A,x,st->w);
35: if (st->sigma != 0.0) {
36: VecAXPY(st->w,-st->sigma,x);
37: }
38: MatMult(st->A,st->w,y);
39: if (st->sigma != 0.0) {
40: VecAXPY(y,-st->sigma,st->w);
41: }
42: }
43: return(0);
44: }
48: PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
49: {
51: ST_FOLD *ctx = (ST_FOLD *) st->data;
54: if (st->B) {
55: /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
56: STAssociatedKSPSolveTranspose(st,x,st->w);
57: MatMult(st->A,st->w,ctx->w2);
58: if (st->sigma != 0.0) {
59: VecAXPY(ctx->w2,-st->sigma,x);
60: }
61: STAssociatedKSPSolveTranspose(st,ctx->w2,st->w);
62: MatMult(st->A,st->w,y);
63: if (st->sigma != 0.0) {
64: VecAXPY(y,-st->sigma,ctx->w2);
65: }
66: } else {
67: /* standard eigenproblem: y = (A^T + sI)^2 x */
68: MatMultTranspose(st->A,x,st->w);
69: if (st->sigma != 0.0) {
70: VecAXPY(st->w,-st->sigma,x);
71: }
72: MatMultTranspose(st->A,st->w,y);
73: if (st->sigma != 0.0) {
74: VecAXPY(y,-st->sigma,st->w);
75: }
76: }
77: return(0);
78: }
82: PetscErrorCode STApplyB_Fold(ST st,Vec x,Vec y)
83: {
87: VecCopy( x, y );
88: return(0);
89: }
93: PetscErrorCode STBackTransform_Fold(ST st,PetscScalar *eigr,PetscScalar *eigi)
94: {
95: ST_FOLD *ctx = (ST_FOLD *) st->data;
99: #if !defined(PETSC_USE_COMPLEX)
100: if (*eigi == 0) {
101: #endif
102: if (ctx->left) *eigr = st->sigma - PetscSqrtScalar(*eigr);
103: else *eigr = st->sigma + PetscSqrtScalar(*eigr);
104: #if !defined(PETSC_USE_COMPLEX)
105: } else {
106: PetscScalar r,x,y;
107: r = PetscSqrtScalar(*eigr * *eigr + *eigi * *eigi);
108: x = PetscSqrtScalar((r + *eigr) / 2);
109: y = PetscSqrtScalar((r - *eigr) / 2);
110: if (*eigi < 0) y = - y;
111: if (ctx->left) {
112: *eigr = st->sigma - x;
113: *eigi = - y;
114: } else {
115: *eigr = st->sigma + x;
116: *eigi = y;
117: }
118: }
119: #endif
120: return(0);
121: }
125: PetscErrorCode STSetUp_Fold(ST st)
126: {
128: ST_FOLD *ctx = (ST_FOLD *) st->data;
131: if (st->B) {
132: KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);
133: KSPSetUp(st->ksp);
134: if (ctx->w2) { VecDestroy(ctx->w2); }
135: MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
136: }
137: return(0);
138: }
143: PetscErrorCode STFoldSetLeftSide_Fold(ST st,PetscTruth left)
144: {
145: ST_FOLD *ctx = (ST_FOLD *) st->data;
148: ctx->left = left;
149: return(0);
150: }
155: /*@
156: STFoldSetLeftSide - Sets a flag to compute eigenvalues on the left side of shift.
158: Collective on ST
160: Input Parameters:
161: + st - the spectral transformation context
162: - left - if true compute eigenvalues on the left side
164: Options Database Key:
165: . -st_fold_leftside - Sets the value of the flag
167: Level: intermediate
169: .seealso: STSetShift()
170: @*/
171: PetscErrorCode STFoldSetLeftSide(ST st,PetscTruth left)
172: {
173: PetscErrorCode ierr, (*f)(ST,PetscTruth);
177: PetscObjectQueryFunction((PetscObject)st,"STFoldSetLeftSide_C",(void (**)(void))&f);
178: if (f) {
179: (*f)(st,left);
180: }
181: return(0);
182: }
186: PetscErrorCode STView_Fold(ST st,PetscViewer viewer)
187: {
189: ST_FOLD *ctx = (ST_FOLD *) st->data;
192: if (ctx->left) {
193: PetscViewerASCIIPrintf(viewer," computing eigenvalues on left side of shift\n");
194: }
195: if (st->B) {
196: STView_Default(st,viewer);
197: }
198: return(0);
199: }
203: PetscErrorCode STSetFromOptions_Fold(ST st)
204: {
206: ST_FOLD *ctx = (ST_FOLD *) st->data;
209: PetscOptionsHead("ST Fold Options");
210: PetscOptionsTruth("-st_fold_leftside","Compute eigenvalues on left side of shift","STFoldSetLeftSide",ctx->left,&ctx->left,PETSC_NULL);
211: PetscOptionsTail();
212: return(0);
213: }
217: PetscErrorCode STDestroy_Fold(ST st)
218: {
220: ST_FOLD *ctx = (ST_FOLD *) st->data;
223: if (ctx->w2) { VecDestroy(ctx->w2); }
224: PetscFree(ctx);
225: return(0);
226: }
231: PetscErrorCode STCreate_Fold(ST st)
232: {
234: ST_FOLD *ctx;
238: PetscNew(ST_FOLD,&ctx);
239: PetscLogObjectMemory(st,sizeof(ST_FOLD));
240: st->data = (void *) ctx;
242: st->ops->apply = STApply_Fold;
243: st->ops->applyB = STApplyB_Fold;
244: st->ops->applynoB = STApply_Fold;
245: st->ops->applytrans = STApplyTranspose_Fold;
246: st->ops->backtr = STBackTransform_Fold;
247: st->ops->setup = STSetUp_Fold;
248: st->ops->view = STView_Fold;
249: st->ops->setfromoptions = STSetFromOptions_Fold;
250: st->ops->destroy = STDestroy_Fold;
251: st->checknullspace = 0;
252:
253: ctx->left = PETSC_FALSE;
254:
255: PetscObjectComposeFunctionDynamic((PetscObject)st,"STFoldSetLeftSide_C","STFoldSetLeftSide_Fold",
256: STFoldSetLeftSide_Fold);
258: return(0);
259: }