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
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <private/stimpl.h> /*I "slepcst.h" I*/
27: typedef struct {
28: Vec w2;
29: } ST_FOLD;
33: PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
34: {
36: ST_FOLD *ctx = (ST_FOLD*)st->data;
39: if (st->B) {
40: /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
41: MatMult(st->A,x,st->w);
42: STAssociatedKSPSolve(st,st->w,ctx->w2);
43: if (st->sigma != 0.0) {
44: VecAXPY(ctx->w2,-st->sigma,x);
45: }
46: MatMult(st->A,ctx->w2,st->w);
47: STAssociatedKSPSolve(st,st->w,y);
48: if (st->sigma != 0.0) {
49: VecAXPY(y,-st->sigma,ctx->w2);
50: }
51: } else {
52: /* standard eigenproblem: y = (A + sI)^2 x */
53: MatMult(st->A,x,st->w);
54: if (st->sigma != 0.0) {
55: VecAXPY(st->w,-st->sigma,x);
56: }
57: MatMult(st->A,st->w,y);
58: if (st->sigma != 0.0) {
59: VecAXPY(y,-st->sigma,st->w);
60: }
61: }
62: return(0);
63: }
67: PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
68: {
70: ST_FOLD *ctx = (ST_FOLD*)st->data;
73: if (st->B) {
74: /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
75: STAssociatedKSPSolveTranspose(st,x,st->w);
76: MatMult(st->A,st->w,ctx->w2);
77: if (st->sigma != 0.0) {
78: VecAXPY(ctx->w2,-st->sigma,x);
79: }
80: STAssociatedKSPSolveTranspose(st,ctx->w2,st->w);
81: MatMult(st->A,st->w,y);
82: if (st->sigma != 0.0) {
83: VecAXPY(y,-st->sigma,ctx->w2);
84: }
85: } else {
86: /* standard eigenproblem: y = (A^T + sI)^2 x */
87: MatMultTranspose(st->A,x,st->w);
88: if (st->sigma != 0.0) {
89: VecAXPY(st->w,-st->sigma,x);
90: }
91: MatMultTranspose(st->A,st->w,y);
92: if (st->sigma != 0.0) {
93: VecAXPY(y,-st->sigma,st->w);
94: }
95: }
96: return(0);
97: }
101: PetscErrorCode STBackTransform_Fold(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
102: {
103: PetscInt j;
106: for (j=0;j<n;j++) {
107: #if !defined(PETSC_USE_COMPLEX)
108: if (eigi[j] == 0) {
109: #endif
110: eigr[j] = st->sigma + PetscSqrtScalar(eigr[j]);
111: #if !defined(PETSC_USE_COMPLEX)
112: } else {
113: PetscScalar r,x,y;
114: r = PetscSqrtScalar(eigr[j] * eigr[j] + eigi[j] * eigi[j]);
115: x = PetscSqrtScalar((r + eigr[j]) / 2);
116: y = PetscSqrtScalar((r - eigr[j]) / 2);
117: if (eigi[j] < 0) y = - y;
118: eigr[j] = st->sigma + x;
119: eigi[j] = y;
120: }
121: #endif
122: }
123: return(0);
124: }
128: PetscErrorCode STSetUp_Fold(ST st)
129: {
131: ST_FOLD *ctx = (ST_FOLD*)st->data;
134: /* if the user did not set the shift, use the target value */
135: if (!st->sigma_set) st->sigma = st->defsigma;
137: if (st->B) {
138: if (!st->ksp) { STGetKSP(st,&st->ksp); }
139: KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);
140: KSPSetUp(st->ksp);
141: VecDestroy(&ctx->w2);
142: MatGetVecs(st->B,&ctx->w2,PETSC_NULL);
143: }
144: return(0);
145: }
149: PetscErrorCode STSetFromOptions_Fold(ST st)
150: {
152: PC pc;
153: const PCType pctype;
154: const KSPType ksptype;
157: if (!st->ksp) { STGetKSP(st,&st->ksp); }
158: KSPGetPC(st->ksp,&pc);
159: KSPGetType(st->ksp,&ksptype);
160: PCGetType(pc,&pctype);
161: if (!pctype && !ksptype) {
162: if (st->shift_matrix == ST_MATMODE_SHELL) {
163: /* in shell mode use GMRES with Jacobi as the default */
164: KSPSetType(st->ksp,KSPGMRES);
165: PCSetType(pc,PCJACOBI);
166: } else {
167: /* use direct solver as default */
168: KSPSetType(st->ksp,KSPPREONLY);
169: PCSetType(pc,PCREDUNDANT);
170: }
171: }
172: return(0);
173: }
177: PetscErrorCode STReset_Fold(ST st)
178: {
180: ST_FOLD *ctx = (ST_FOLD*)st->data;
183: VecDestroy(&ctx->w2);
184: return(0);
185: }
189: PetscErrorCode STDestroy_Fold(ST st)
190: {
194: PetscFree(st->data);
195: return(0);
196: }
198: EXTERN_C_BEGIN
201: PetscErrorCode STCreate_Fold(ST st)
202: {
206: PetscNewLog(st,ST_FOLD,&st->data);
207: st->ops->apply = STApply_Fold;
208: st->ops->getbilinearform = STGetBilinearForm_Default;
209: st->ops->applytrans = STApplyTranspose_Fold;
210: st->ops->backtr = STBackTransform_Fold;
211: st->ops->setup = STSetUp_Fold;
212: st->ops->setfromoptions = STSetFromOptions_Fold;
213: st->ops->destroy = STDestroy_Fold;
214: st->ops->reset = STReset_Fold;
215: return(0);
216: }
217: EXTERN_C_END