Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
778 dsic.upv.es!antodo 1
/*
804 dsic.upv.es!antodo 2
    Folding spectral transformation, applies (A + sigma I)^2 as operator, or
3
    inv(B)(A + sigma I)^2 for generalized problems
778 dsic.upv.es!antodo 4
*/
5
#include "src/st/stimpl.h"          /*I "slepcst.h" I*/
6
 
7
typedef struct {
8
  PetscTruth  left;
9
  Vec         w2;
10
} ST_FOLD;
11
 
12
#undef __FUNCT__  
13
#define __FUNCT__ "STApply_Fold"
14
PetscErrorCode STApply_Fold(ST st,Vec x,Vec y)
15
{
16
  PetscErrorCode ierr;
17
  ST_FOLD        *ctx = (ST_FOLD *) st->data;
18
 
19
  PetscFunctionBegin;
20
  if (st->B) {
804 dsic.upv.es!antodo 21
    /* generalized eigenproblem: y = (B^-1 A + sI)^2 x */
828 dsic.upv.es!antodo 22
    ierr = MatMult(st->A,x,st->w);CHKERRQ(ierr);
23
    ierr = STAssociatedKSPSolve(st,st->w,ctx->w2);CHKERRQ(ierr);
778 dsic.upv.es!antodo 24
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 25
      ierr = VecAXPY(ctx->w2,-st->sigma,x);CHKERRQ(ierr);
778 dsic.upv.es!antodo 26
    }
828 dsic.upv.es!antodo 27
    ierr = MatMult(st->A,ctx->w2,st->w);CHKERRQ(ierr);
28
    ierr = STAssociatedKSPSolve(st,st->w,y);CHKERRQ(ierr);
29
    if (st->sigma != 0.0) {
30
      ierr = VecAXPY(y,-st->sigma,ctx->w2);CHKERRQ(ierr);
31
    }
778 dsic.upv.es!antodo 32
  } else {
804 dsic.upv.es!antodo 33
    /* standard eigenproblem: y = (A + sI)^2 x */
778 dsic.upv.es!antodo 34
    ierr = MatMult(st->A,x,st->w);CHKERRQ(ierr);
35
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 36
      ierr = VecAXPY(st->w,-st->sigma,x);CHKERRQ(ierr);
778 dsic.upv.es!antodo 37
    }
38
    ierr = MatMult(st->A,st->w,y);CHKERRQ(ierr);
39
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 40
      ierr = VecAXPY(y,-st->sigma,st->w);CHKERRQ(ierr);
778 dsic.upv.es!antodo 41
    }
42
  }
43
  PetscFunctionReturn(0);
44
}
45
 
46
#undef __FUNCT__  
813 dsic.upv.es!antodo 47
#define __FUNCT__ "STApplyTranspose_Fold"
48
PetscErrorCode STApplyTranspose_Fold(ST st,Vec x,Vec y)
49
{
50
  PetscErrorCode ierr;
51
  ST_FOLD        *ctx = (ST_FOLD *) st->data;
52
 
53
  PetscFunctionBegin;
54
  if (st->B) {
55
    /* generalized eigenproblem: y = (A^T B^-T + sI)^2 x */
828 dsic.upv.es!antodo 56
    ierr = STAssociatedKSPSolveTranspose(st,x,st->w);CHKERRQ(ierr);
57
    ierr = MatMult(st->A,st->w,ctx->w2);CHKERRQ(ierr);
813 dsic.upv.es!antodo 58
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 59
      ierr = VecAXPY(ctx->w2,-st->sigma,x);CHKERRQ(ierr);
813 dsic.upv.es!antodo 60
    }
828 dsic.upv.es!antodo 61
    ierr = STAssociatedKSPSolveTranspose(st,ctx->w2,st->w);CHKERRQ(ierr);
62
    ierr = MatMult(st->A,st->w,y);CHKERRQ(ierr);
63
    if (st->sigma != 0.0) {
64
      ierr = VecAXPY(y,-st->sigma,ctx->w2);CHKERRQ(ierr);
65
    }
813 dsic.upv.es!antodo 66
  } else {
67
    /* standard eigenproblem: y = (A^T + sI)^2 x */
68
    ierr = MatMultTranspose(st->A,x,st->w);CHKERRQ(ierr);
69
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 70
      ierr = VecAXPY(st->w,-st->sigma,x);CHKERRQ(ierr);
813 dsic.upv.es!antodo 71
    }
72
    ierr = MatMultTranspose(st->A,st->w,y);CHKERRQ(ierr);
73
    if (st->sigma != 0.0) {
828 dsic.upv.es!antodo 74
      ierr = VecAXPY(y,-st->sigma,st->w);CHKERRQ(ierr);
813 dsic.upv.es!antodo 75
    }
76
  }
77
  PetscFunctionReturn(0);
78
}
79
 
80
#undef __FUNCT__  
778 dsic.upv.es!antodo 81
#define __FUNCT__ "STBackTransform_Fold"
82
PetscErrorCode STBackTransform_Fold(ST st,PetscScalar *eigr,PetscScalar *eigi)
83
{
84
  ST_FOLD *ctx = (ST_FOLD *) st->data;
85
  PetscFunctionBegin;
86
  PetscValidScalarPointer(eigr,2);
87
  PetscValidScalarPointer(eigi,3);
88
#if !defined(PETSC_USE_COMPLEX)
89
  if (*eigi == 0) {
90
#endif
91
    if (ctx->left) *eigr = st->sigma - PetscSqrtScalar(*eigr);
92
    else *eigr = st->sigma + PetscSqrtScalar(*eigr);
93
#if !defined(PETSC_USE_COMPLEX)
94
  } else {
95
    PetscScalar r,x,y;
96
    r = PetscSqrtScalar(*eigr * *eigr + *eigi * *eigi);
97
    x = PetscSqrtScalar((r + *eigr) / 2);
98
    y = PetscSqrtScalar((r - *eigr) / 2);
99
    if (*eigi < 0) y = - y;
100
    if (ctx->left) {
101
      *eigr = st->sigma - x;
102
      *eigi = - y;
103
    } else {
104
      *eigr = st->sigma + x;
105
      *eigi = y;
106
    }
107
  }
108
#endif
109
  PetscFunctionReturn(0);
110
}
111
 
112
#undef __FUNCT__  
113
#define __FUNCT__ "STSetUp_Fold"
114
PetscErrorCode STSetUp_Fold(ST st)
115
{
116
  PetscErrorCode ierr;
117
  ST_FOLD        *ctx = (ST_FOLD *) st->data;
118
 
119
  PetscFunctionBegin;
120
  if (st->B) {
121
    ierr = KSPSetOperators(st->ksp,st->B,st->B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
122
    ierr = KSPSetUp(st->ksp);CHKERRQ(ierr);
123
    if (ctx->w2) { ierr = VecDestroy(ctx->w2);CHKERRQ(ierr); }
124
    ierr = MatGetVecs(st->B,&ctx->w2,PETSC_NULL);CHKERRQ(ierr);
125
  }
126
  PetscFunctionReturn(0);
127
}
128
 
129
EXTERN_C_BEGIN
130
#undef __FUNCT__  
131
#define __FUNCT__ "STFoldSetLeftSide_Fold"
132
PetscErrorCode STFoldSetLeftSide_Fold(ST st,PetscTruth left)
133
{
134
  ST_FOLD *ctx = (ST_FOLD *) st->data;
135
 
136
  PetscFunctionBegin;
137
  ctx->left = left;
138
  PetscFunctionReturn(0);
139
}
140
EXTERN_C_END
141
 
142
#undef __FUNCT__  
143
#define __FUNCT__ "STFoldSetLeftSide"
144
/*@
145
   STFoldSetLeftSide - Sets a flag to compute eigenvalues on the left side of shift.
146
 
147
   Collective on ST
148
 
149
   Input Parameters:
150
+  st  - the spectral transformation context
151
-  left - if true compute eigenvalues on the left side
152
 
153
   Options Database Key:
154
.  -st_fold_leftside - Sets the value of the flag
155
 
156
   Level: intermediate
157
 
158
.seealso: STSetShift()
159
@*/
160
PetscErrorCode STFoldSetLeftSide(ST st,PetscTruth left)
161
{
162
  PetscErrorCode ierr, (*f)(ST,PetscTruth);
163
 
164
  PetscFunctionBegin;
165
  PetscValidHeaderSpecific(st,ST_COOKIE,1);
166
  ierr = PetscObjectQueryFunction((PetscObject)st,"STFoldSetLeftSide_C",(void (**)(void))&f);CHKERRQ(ierr);
167
  if (f) {
168
    ierr = (*f)(st,left);CHKERRQ(ierr);
169
  }
170
  PetscFunctionReturn(0);
171
}
172
 
173
#undef __FUNCT__  
174
#define __FUNCT__ "STView_Fold"
175
PetscErrorCode STView_Fold(ST st,PetscViewer viewer)
176
{
177
  PetscErrorCode ierr;
178
  ST_FOLD        *ctx = (ST_FOLD *) st->data;
179
 
180
  PetscFunctionBegin;
181
  if (ctx->left) {
182
    ierr = PetscViewerASCIIPrintf(viewer,"  computing eigenvalues on left side of shift\n");CHKERRQ(ierr);
183
  }  
184
  if (st->B) {
185
    ierr = STView_Default(st,viewer);CHKERRQ(ierr);
186
  }
187
  PetscFunctionReturn(0);
188
}
189
 
190
#undef __FUNCT__  
191
#define __FUNCT__ "STSetFromOptions_Fold"
192
PetscErrorCode STSetFromOptions_Fold(ST st)
193
{
194
  PetscErrorCode ierr;
195
  ST_FOLD      *ctx = (ST_FOLD *) st->data;
196
 
197
  PetscFunctionBegin;
198
  ierr = PetscOptionsHead("ST Fold Options");CHKERRQ(ierr);
830 dsic.upv.es!antodo 199
  ierr = PetscOptionsTruth("-st_fold_leftside","Compute eigenvalues on left side of shift","STFoldSetLeftSide",ctx->left,&ctx->left,PETSC_NULL); CHKERRQ(ierr);
778 dsic.upv.es!antodo 200
  ierr = PetscOptionsTail();CHKERRQ(ierr);
201
  PetscFunctionReturn(0);
202
}
203
 
204
#undef __FUNCT__  
205
#define __FUNCT__ "STDestroy_Fold"
206
PetscErrorCode STDestroy_Fold(ST st)
207
{
208
  PetscErrorCode ierr;
209
  ST_FOLD        *ctx = (ST_FOLD *) st->data;
210
 
211
  PetscFunctionBegin;
212
  if (ctx->w2) { ierr = VecDestroy(ctx->w2);CHKERRQ(ierr); }
213
  ierr = PetscFree(ctx);CHKERRQ(ierr);
214
  PetscFunctionReturn(0);
215
}
216
 
217
EXTERN_C_BEGIN
218
#undef __FUNCT__  
219
#define __FUNCT__ "STCreate_Fold"
220
PetscErrorCode STCreate_Fold(ST st)
221
{
222
  PetscErrorCode ierr;
223
  ST_FOLD        *ctx;
224
 
225
  PetscFunctionBegin;
226
 
227
  ierr = PetscNew(ST_FOLD,&ctx); CHKERRQ(ierr);
228
  PetscLogObjectMemory(st,sizeof(ST_FOLD));
229
  st->data                = (void *) ctx;
230
 
1358 slepc 231
  st->ops->apply           = STApply_Fold;
232
  st->ops->getbilinearform = STGetBilinearForm_Default;
233
  st->ops->applytrans      = STApplyTranspose_Fold;
234
  st->ops->backtr          = STBackTransform_Fold;
235
  st->ops->setup           = STSetUp_Fold;
236
  st->ops->view            = STView_Fold;
237
  st->ops->setfromoptions  = STSetFromOptions_Fold;
238
  st->ops->destroy         = STDestroy_Fold;
239
  st->checknullspace       = 0;
778 dsic.upv.es!antodo 240
 
241
  ctx->left            = PETSC_FALSE;
242
 
243
  ierr = PetscObjectComposeFunctionDynamic((PetscObject)st,"STFoldSetLeftSide_C","STFoldSetLeftSide_Fold",
244
                    STFoldSetLeftSide_Fold);CHKERRQ(ierr);
245
 
246
  PetscFunctionReturn(0);
247
}
248
EXTERN_C_END
249