Actual source code: setup.c
1: /*
2: EPS routines related to problem setup.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc. See the README file for conditions of use
9: and additional information.
10: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
11: */
13: #include src/eps/epsimpl.h
17: /*@
18: EPSSetUp - Sets up all the internal data structures necessary for the
19: execution of the eigensolver. Then calls STSetUp() for any set-up
20: operations associated to the ST object.
22: Collective on EPS
24: Input Parameter:
25: . eps - eigenproblem solver context
27: Level: advanced
29: Notes:
30: This function need not be called explicitly in most cases, since EPSSolve()
31: calls it. It can be useful when one wants to measure the set-up time
32: separately from the solve time.
34: This function sets a random initial vector if none has been provided.
36: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
37: @*/
38: PetscErrorCode EPSSetUp(EPS eps)
39: {
41: int i;
42: Vec v0,w0;
43: Mat A,B;
44: PetscInt N;
45:
49: if (eps->setupcalled) return(0);
51: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
53: /* Set default solver type */
54: if (!eps->type_name) {
55: EPSSetType(eps,EPSKRYLOVSCHUR);
56: }
57:
58: STGetOperators(eps->OP,&A,&B);
59: /* Set default problem type */
60: if (!eps->problem_type) {
61: if (B==PETSC_NULL) {
62: EPSSetProblemType(eps,EPS_NHEP);
63: }
64: else {
65: EPSSetProblemType(eps,EPS_GNHEP);
66: }
67: } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
68: SETERRQ(0,"Warning: Inconsistent EPS state");
69: }
70:
71: if (eps->ispositive) {
72: STGetBilinearForm(eps->OP,&B);
73: IPSetBilinearForm(eps->ip,B,IPINNER_HERMITIAN);
74: MatDestroy(B);
75: } else {
76: IPSetBilinearForm(eps->ip,PETSC_NULL,IPINNER_HERMITIAN);
77: }
78:
79: /* Create random initial vectors if not set */
80: /* right */
81: EPSGetInitialVector(eps,&v0);
82: if (!v0) {
83: MatGetVecs(A,&v0,PETSC_NULL);
84: SlepcVecSetRandom(v0);
85: eps->vec_initial = v0;
86: }
87: /* left */
88: EPSGetLeftInitialVector(eps,&w0);
89: if (!w0) {
90: MatGetVecs(A,PETSC_NULL,&w0);
91: SlepcVecSetRandom(w0);
92: eps->vec_initial_left = w0;
93: }
95: VecGetSize(eps->vec_initial,&N);
96: if (eps->nev > N) eps->nev = N;
97: if (eps->ncv > N) eps->ncv = N;
99: (*eps->ops->setup)(eps);
100: STSetUp(eps->OP);
102: /* DSV is equal to the columns of DS followed by the ones in V */
103: PetscFree(eps->DSV);
104: PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);
105: for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
106: for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
107:
108: if (eps->nds>0) {
109: if (!eps->ds_ortho) {
110: /* orthonormalize vectors in DS if necessary */
111: IPQRDecomposition(eps->ip,eps->DS,0,eps->nds,PETSC_NULL,0,PETSC_NULL);
112: }
113: IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
114: }
116: STCheckNullSpace(eps->OP,eps->nds,eps->DS);
117:
118: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
119: eps->setupcalled = 1;
120: return(0);
121: }
125: /*@
126: EPSSetInitialVector - Sets the initial vector from which the
127: eigensolver starts to iterate.
129: Collective on EPS and Vec
131: Input Parameters:
132: + eps - the eigensolver context
133: - vec - the vector
135: Level: intermediate
137: .seealso: EPSGetInitialVector(), EPSSetLeftInitialVector()
139: @*/
140: PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
141: {
143:
148: PetscObjectReference((PetscObject)vec);
149: if (eps->vec_initial) {
150: VecDestroy(eps->vec_initial);
151: }
152: eps->vec_initial = vec;
153: return(0);
154: }
158: /*@
159: EPSGetInitialVector - Gets the initial vector associated with the
160: eigensolver; if the vector was not set it will return a 0 pointer or
161: a vector randomly generated by EPSSetUp().
163: Not collective, but vector is shared by all processors that share the EPS
165: Input Parameter:
166: . eps - the eigensolver context
168: Output Parameter:
169: . vec - the vector
171: Level: intermediate
173: .seealso: EPSSetInitialVector(), EPSGetLeftInitialVector()
175: @*/
176: PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
177: {
181: *vec = eps->vec_initial;
182: return(0);
183: }
187: /*@
188: EPSSetLeftInitialVector - Sets the initial vector from which the eigensolver
189: starts to iterate, corresponding to the left recurrence (two-sided solvers).
191: Collective on EPS and Vec
193: Input Parameters:
194: + eps - the eigensolver context
195: - vec - the vector
197: Level: intermediate
199: .seealso: EPSGetLeftInitialVector(), EPSSetInitialVector()
201: @*/
202: PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
203: {
205:
210: PetscObjectReference((PetscObject)vec);
211: if (eps->vec_initial_left) {
212: VecDestroy(eps->vec_initial_left);
213: }
214: eps->vec_initial_left = vec;
215: return(0);
216: }
220: /*@
221: EPSGetLeftInitialVector - Gets the left initial vector associated with the
222: eigensolver; if the vector was not set it will return a 0 pointer or
223: a vector randomly generated by EPSSetUp().
225: Not collective, but vector is shared by all processors that share the EPS
227: Input Parameter:
228: . eps - the eigensolver context
230: Output Parameter:
231: . vec - the vector
233: Level: intermediate
235: .seealso: EPSSetLeftInitialVector(), EPSGetLeftInitialVector()
237: @*/
238: PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
239: {
243: *vec = eps->vec_initial_left;
244: return(0);
245: }
249: /*@
250: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
252: Collective on EPS and Mat
254: Input Parameters:
255: + eps - the eigenproblem solver context
256: . A - the matrix associated with the eigensystem
257: - B - the second matrix in the case of generalized eigenproblems
259: Notes:
260: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
262: Level: beginner
264: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
265: @*/
266: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
267: {
269: PetscInt m,n;
278: /* Check for square matrices */
279: MatGetSize(A,&m,&n);
280: if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
281: if (B) {
282: MatGetSize(B,&m,&n);
283: if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
284: }
286: STSetOperators(eps->OP,A,B);
287: eps->setupcalled = 0; /* so that next solve call will call setup */
289: /* Destroy randomly generated initial vectors */
290: if (eps->vec_initial) {
291: VecDestroy(eps->vec_initial);
292: eps->vec_initial = PETSC_NULL;
293: }
294: if (eps->vec_initial_left) {
295: VecDestroy(eps->vec_initial_left);
296: eps->vec_initial_left = PETSC_NULL;
297: }
299: return(0);
300: }
304: /*@
305: EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.
307: Not Collective
309: Input Parameter:
310: + eps - the eigenproblem solver context
311: . n - number of vectors to add
312: . ds - set of basis vectors of the deflation space
313: - ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal
315: Notes:
316: When a deflation space is given, the eigensolver seeks the eigensolution
317: in the restriction of the problem to the orthogonal complement of this
318: space. This can be used for instance in the case that an invariant
319: subspace is known beforehand (such as the nullspace of the matrix).
321: The basis vectors can be provided all at once or incrementally with
322: several calls to EPSAttachDeflationSpace().
324: Use a value of PETSC_TRUE for parameter ortho if all the vectors passed
325: in are known to be mutually orthonormal.
327: Level: intermediate
329: .seealso: EPSRemoveDeflationSpace()
330: @*/
331: PetscErrorCode EPSAttachDeflationSpace(EPS eps,int n,Vec *ds,PetscTruth ortho)
332: {
334: int i;
335: Vec *tvec;
336:
339: tvec = eps->DS;
340: if (n+eps->nds > 0) {
341: PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);
342: }
343: if (eps->nds > 0) {
344: for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
345: PetscFree(tvec);
346: }
347: for (i=0; i<n; i++) {
348: VecDuplicate(ds[i],&eps->DS[i + eps->nds]);
349: VecCopy(ds[i],eps->DS[i + eps->nds]);
350: }
351: eps->nds += n;
352: if (!ortho) eps->ds_ortho = PETSC_FALSE;
353: eps->setupcalled = 0;
354: return(0);
355: }
359: /*@
360: EPSRemoveDeflationSpace - Removes the deflation space.
362: Not Collective
364: Input Parameter:
365: . eps - the eigenproblem solver context
367: Level: intermediate
369: .seealso: EPSAttachDeflationSpace()
370: @*/
371: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
372: {
374:
377: if (eps->nds > 0) {
378: VecDestroyVecs(eps->DS, eps->nds);
379: }
380: eps->ds_ortho = PETSC_TRUE;
381: eps->setupcalled = 0;
382: return(0);
383: }