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
545 dsic.upv.es!jroman 1
/*
2
      EPS routines related to problem setup.
3
*/
527 dsic.upv.es!antodo 4
#include "src/eps/epsimpl.h"   /*I "slepceps.h" I*/
5
 
6
#undef __FUNCT__  
7
#define __FUNCT__ "EPSSetUp"
8
/*@
9
   EPSSetUp - Sets up all the internal data structures necessary for the
10
   execution of the eigensolver. Then calls STSetUp() for any set-up
11
   operations associated to the ST object.
12
 
13
   Collective on EPS
14
 
15
   Input Parameter:
16
.  eps   - eigenproblem solver context
17
 
18
   Level: advanced
19
 
20
   Notes:
21
   This function need not be called explicitly in most cases, since EPSSolve()
22
   calls it. It can be useful when one wants to measure the set-up time
23
   separately from the solve time.
24
 
25
   This function sets a random initial vector if none has been provided.
26
 
27
.seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
28
@*/
29
PetscErrorCode EPSSetUp(EPS eps)
30
{
31
  PetscErrorCode ierr;
32
  int            i;  
780 dsic.upv.es!jroman 33
  Vec            v0,w0;  
527 dsic.upv.es!antodo 34
  Mat            A,B;
35
 
36
  PetscFunctionBegin;
37
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
38
 
39
  if (eps->setupcalled) PetscFunctionReturn(0);
40
 
41
  ierr = PetscLogEventBegin(EPS_SetUp,eps,0,0,0);CHKERRQ(ierr);
42
 
43
  /* Set default solver type */
44
  if (!eps->type_name) {
45
    ierr = EPSSetType(eps,EPSARNOLDI);CHKERRQ(ierr);
46
  }
47
 
691 dsic.upv.es!antodo 48
  /* Set default eta for orthogonalization */
49
  if (eps->orthog_eta == PETSC_DEFAULT) {
50
    if (eps->orthog_type == EPS_MGS_ORTH) eps->orthog_eta = 0.7071;
51
    else eps->orthog_eta = 0.5;
52
  }
53
 
527 dsic.upv.es!antodo 54
  ierr = STGetOperators(eps->OP,&A,&B);CHKERRQ(ierr);
55
  /* Set default problem type */
56
  if (!eps->problem_type) {
57
    if (B==PETSC_NULL) {
58
      ierr = EPSSetProblemType(eps,EPS_NHEP);CHKERRQ(ierr);
59
    }
60
    else {
61
      ierr = EPSSetProblemType(eps,EPS_GNHEP);CHKERRQ(ierr);
62
    }
63
  } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
64
    SETERRQ(0,"Warning: Inconsistent EPS state");
65
  }
66
 
780 dsic.upv.es!jroman 67
  /* Create random initial vectors if not set */
68
  /* right */
527 dsic.upv.es!antodo 69
  ierr = EPSGetInitialVector(eps,&v0);CHKERRQ(ierr);
987 slepc 70
  if (!v0) {
527 dsic.upv.es!antodo 71
    ierr = MatGetVecs(A,&v0,PETSC_NULL);CHKERRQ(ierr);
72
    ierr = SlepcVecSetRandom(v0);CHKERRQ(ierr);
73
    eps->vec_initial = v0;
74
  }
780 dsic.upv.es!jroman 75
  /* left */
76
  ierr = EPSGetLeftInitialVector(eps,&w0);CHKERRQ(ierr);
987 slepc 77
  if (!w0) {
780 dsic.upv.es!jroman 78
    ierr = MatGetVecs(A,PETSC_NULL,&w0);CHKERRQ(ierr);
79
    ierr = SlepcVecSetRandom(w0);CHKERRQ(ierr);
80
    eps->vec_initial_left = w0;
81
  }
527 dsic.upv.es!antodo 82
 
83
  ierr = (*eps->ops->setup)(eps);CHKERRQ(ierr);
543 dsic.upv.es!jroman 84
  ierr = STSetUp(eps->OP); CHKERRQ(ierr);
527 dsic.upv.es!antodo 85
 
86
  /* DSV is equal to the columns of DS followed by the ones in V */
1040 slepc 87
  ierr = PetscFree(eps->DSV);CHKERRQ(ierr);
527 dsic.upv.es!antodo 88
  ierr = PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);CHKERRQ(ierr);    
89
  for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
90
  for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
91
 
92
  if (eps->nds>0) {
93
    if (!eps->ds_ortho) {
94
      /* orthonormalize vectors in DS if necessary */
95
      ierr = EPSQRDecomposition(eps,eps->DS,0,eps->nds,PETSC_NULL,0);CHKERRQ(ierr);
96
    }
97
    ierr = EPSOrthogonalize(eps,eps->nds,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
98
  }
99
 
100
  ierr = STCheckNullSpace(eps->OP,eps->nds,eps->DS);CHKERRQ(ierr);
101
 
102
  ierr = PetscLogEventEnd(EPS_SetUp,eps,0,0,0);CHKERRQ(ierr);
103
  eps->setupcalled = 1;
104
  PetscFunctionReturn(0);
105
}
106
 
107
#undef __FUNCT__  
108
#define __FUNCT__ "EPSSetInitialVector"
109
/*@
110
   EPSSetInitialVector - Sets the initial vector from which the
111
   eigensolver starts to iterate.
112
 
113
   Collective on EPS and Vec
114
 
115
   Input Parameters:
116
+  eps - the eigensolver context
117
-  vec - the vector
118
 
119
   Level: intermediate
120
 
780 dsic.upv.es!jroman 121
.seealso: EPSGetInitialVector(), EPSSetLeftInitialVector()
527 dsic.upv.es!antodo 122
 
123
@*/
124
PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
125
{
126
  PetscErrorCode ierr;
127
 
128
  PetscFunctionBegin;
129
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
130
  PetscValidHeaderSpecific(vec,VEC_COOKIE,2);
131
  PetscCheckSameComm(eps,1,vec,2);
132
  if (eps->vec_initial) {
133
    ierr = VecDestroy(eps->vec_initial); CHKERRQ(ierr);
134
  }
135
  eps->vec_initial = vec;
863 dsic.upv.es!jroman 136
  ierr = PetscObjectReference((PetscObject)eps->vec_initial);CHKERRQ(ierr);
527 dsic.upv.es!antodo 137
  PetscFunctionReturn(0);
138
}
139
 
140
#undef __FUNCT__  
141
#define __FUNCT__ "EPSGetInitialVector"
142
/*@
143
   EPSGetInitialVector - Gets the initial vector associated with the
144
   eigensolver; if the vector was not set it will return a 0 pointer or
145
   a vector randomly generated by EPSSetUp().
146
 
147
   Not collective, but vector is shared by all processors that share the EPS
148
 
149
   Input Parameter:
150
.  eps - the eigensolver context
151
 
152
   Output Parameter:
153
.  vec - the vector
154
 
155
   Level: intermediate
156
 
780 dsic.upv.es!jroman 157
.seealso: EPSSetInitialVector(), EPSGetLeftInitialVector()
527 dsic.upv.es!antodo 158
 
159
@*/
160
PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
161
{
162
  PetscFunctionBegin;
163
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
164
  *vec = eps->vec_initial;
165
  PetscFunctionReturn(0);
166
}
167
 
168
#undef __FUNCT__  
780 dsic.upv.es!jroman 169
#define __FUNCT__ "EPSSetLeftInitialVector"
170
/*@
794 dsic.upv.es!antodo 171
   EPSSetLeftInitialVector - Sets the initial vector from which the eigensolver
780 dsic.upv.es!jroman 172
   starts to iterate, corresponding to the left recurrence (two-sided solvers).
173
 
174
   Collective on EPS and Vec
175
 
176
   Input Parameters:
177
+  eps - the eigensolver context
178
-  vec - the vector
179
 
180
   Level: intermediate
181
 
182
.seealso: EPSGetLeftInitialVector(), EPSSetInitialVector()
183
 
184
@*/
185
PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
186
{
187
  PetscErrorCode ierr;
188
 
189
  PetscFunctionBegin;
190
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
191
  PetscValidHeaderSpecific(vec,VEC_COOKIE,2);
192
  PetscCheckSameComm(eps,1,vec,2);
193
  if (eps->vec_initial_left) {
194
    ierr = VecDestroy(eps->vec_initial_left); CHKERRQ(ierr);
195
  }
196
  eps->vec_initial_left = vec;
863 dsic.upv.es!jroman 197
  ierr = PetscObjectReference((PetscObject)eps->vec_initial_left);CHKERRQ(ierr);
780 dsic.upv.es!jroman 198
  PetscFunctionReturn(0);
199
}
200
 
201
#undef __FUNCT__  
202
#define __FUNCT__ "EPSGetLeftInitialVector"
203
/*@
794 dsic.upv.es!antodo 204
   EPSGetLeftInitialVector - Gets the left initial vector associated with the
780 dsic.upv.es!jroman 205
   eigensolver; if the vector was not set it will return a 0 pointer or
206
   a vector randomly generated by EPSSetUp().
207
 
208
   Not collective, but vector is shared by all processors that share the EPS
209
 
210
   Input Parameter:
211
.  eps - the eigensolver context
212
 
213
   Output Parameter:
214
.  vec - the vector
215
 
216
   Level: intermediate
217
 
218
.seealso: EPSSetLeftInitialVector(), EPSGetLeftInitialVector()
219
 
220
@*/
221
PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
222
{
223
  PetscFunctionBegin;
224
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
225
  *vec = eps->vec_initial_left;
226
  PetscFunctionReturn(0);
227
}
228
 
229
#undef __FUNCT__  
527 dsic.upv.es!antodo 230
#define __FUNCT__ "EPSSetOperators"
231
/*@
232
   EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
233
 
234
   Collective on EPS and Mat
235
 
236
   Input Parameters:
237
+  eps - the eigenproblem solver context
238
.  A  - the matrix associated with the eigensystem
239
-  B  - the second matrix in the case of generalized eigenproblems
240
 
241
   Notes:
242
   To specify a standard eigenproblem, use PETSC_NULL for parameter B.
243
 
244
   Level: beginner
245
 
246
.seealso: EPSSolve(), EPSGetST(), STGetOperators()
247
@*/
248
PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
249
{
250
  PetscErrorCode ierr;
982 slepc 251
  PetscInt       m,n;
527 dsic.upv.es!antodo 252
 
253
  PetscFunctionBegin;
254
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
255
  PetscValidHeaderSpecific(A,MAT_COOKIE,2);
256
  if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
1014 slepc 257
  PetscCheckSameComm(eps,1,A,2);
258
  if (B) PetscCheckSameComm(eps,1,B,3);
527 dsic.upv.es!antodo 259
 
260
  /* Check for square matrices */
261
  ierr = MatGetSize(A,&m,&n);CHKERRQ(ierr);
262
  if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
263
  if (B) {
264
    ierr = MatGetSize(B,&m,&n);CHKERRQ(ierr);
265
    if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
266
  }
267
 
268
  ierr = STSetOperators(eps->OP,A,B);CHKERRQ(ierr);
269
  eps->setupcalled = 0;  /* so that next solve call will call setup */
780 dsic.upv.es!jroman 270
 
271
  /* Destroy randomly generated initial vectors */
987 slepc 272
  if (eps->vec_initial) {
527 dsic.upv.es!antodo 273
    ierr = VecDestroy(eps->vec_initial);CHKERRQ(ierr);
274
    eps->vec_initial = PETSC_NULL;
275
  }
987 slepc 276
  if (eps->vec_initial_left) {
780 dsic.upv.es!jroman 277
    ierr = VecDestroy(eps->vec_initial_left);CHKERRQ(ierr);
278
    eps->vec_initial_left = PETSC_NULL;
279
  }
527 dsic.upv.es!antodo 280
 
281
  PetscFunctionReturn(0);
282
}
283
 
284
#undef __FUNCT__  
285
#define __FUNCT__ "EPSAttachDeflationSpace"
286
/*@
287
   EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.
288
 
289
   Not Collective
290
 
291
   Input Parameter:
292
+  eps   - the eigenproblem solver context
293
.  n     - number of vectors to add
294
.  ds    - set of basis vectors of the deflation space
295
-  ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal
296
 
297
   Notes:
298
   When a deflation space is given, the eigensolver seeks the eigensolution
299
   in the restriction of the problem to the orthogonal complement of this
300
   space. This can be used for instance in the case that an invariant
301
   subspace is known beforehand (such as the nullspace of the matrix).
302
 
303
   The basis vectors can be provided all at once or incrementally with
304
   several calls to EPSAttachDeflationSpace().
305
 
306
   Use a value of PETSC_TRUE for parameter ortho if all the vectors passed
307
   in are known to be mutually orthonormal.
308
 
309
   Level: intermediate
310
 
311
.seealso: EPSRemoveDeflationSpace()
312
@*/
313
PetscErrorCode EPSAttachDeflationSpace(EPS eps,int n,Vec *ds,PetscTruth ortho)
314
{
315
  PetscErrorCode ierr;
316
  int            i;
317
  Vec            *tvec;
318
 
319
  PetscFunctionBegin;
320
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
321
  tvec = eps->DS;
322
  if (n+eps->nds > 0) {
323
     ierr = PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);CHKERRQ(ierr);
324
  }
325
  if (eps->nds > 0) {
326
    for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
327
    ierr = PetscFree(tvec);CHKERRQ(ierr);
328
  }
329
  for (i=0; i<n; i++) {
330
    ierr = VecDuplicate(ds[i],&eps->DS[i + eps->nds]);CHKERRQ(ierr);
331
    ierr = VecCopy(ds[i],eps->DS[i + eps->nds]);CHKERRQ(ierr);
332
  }
333
  eps->nds += n;
334
  if (!ortho) eps->ds_ortho = PETSC_FALSE;
335
  eps->setupcalled = 0;
336
  PetscFunctionReturn(0);
337
}
338
 
339
#undef __FUNCT__  
340
#define __FUNCT__ "EPSRemoveDeflationSpace"
341
/*@
342
   EPSRemoveDeflationSpace - Removes the deflation space.
343
 
344
   Not Collective
345
 
346
   Input Parameter:
347
.  eps   - the eigenproblem solver context
348
 
349
   Level: intermediate
350
 
351
.seealso: EPSAttachDeflationSpace()
352
@*/
353
PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
354
{
355
  PetscErrorCode ierr;
356
 
357
  PetscFunctionBegin;
358
  PetscValidHeaderSpecific(eps,EPS_COOKIE,1);
359
  if (eps->nds > 0) {
360
    ierr = VecDestroyVecs(eps->DS, eps->nds);CHKERRQ(ierr);
361
  }
362
  eps->ds_ortho = PETSC_TRUE;
363
  eps->setupcalled = 0;
364
  PetscFunctionReturn(0);
365
}