Actual source code: setup.c
1: /*
2: EPS routines related to problem setup.
3: */
4: #include src/eps/epsimpl.h
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.
13: Collective on EPS
15: Input Parameter:
16: . eps - eigenproblem solver context
18: Level: advanced
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.
25: This function sets a random initial vector if none has been provided.
27: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
28: @*/
29: PetscErrorCode EPSSetUp(EPS eps)
30: {
32: int i;
33: Vec v0,w0;
34: Mat A,B;
35: PetscInt N;
36:
40: if (eps->setupcalled) return(0);
42: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
44: /* Set default solver type */
45: if (!eps->type_name) {
46: EPSSetType(eps,EPSKRYLOVSCHUR);
47: }
49: /* Set default eta for orthogonalization */
50: if (eps->orthog_eta == PETSC_DEFAULT)
51: eps->orthog_eta = 0.7071;
52:
53: STGetOperators(eps->OP,&A,&B);
54: /* Set default problem type */
55: if (!eps->problem_type) {
56: if (B==PETSC_NULL) {
57: EPSSetProblemType(eps,EPS_NHEP);
58: }
59: else {
60: EPSSetProblemType(eps,EPS_GNHEP);
61: }
62: } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
63: SETERRQ(0,"Warning: Inconsistent EPS state");
64: }
65:
66: /* Create random initial vectors if not set */
67: /* right */
68: EPSGetInitialVector(eps,&v0);
69: if (!v0) {
70: MatGetVecs(A,&v0,PETSC_NULL);
71: SlepcVecSetRandom(v0);
72: eps->vec_initial = v0;
73: }
74: /* left */
75: EPSGetLeftInitialVector(eps,&w0);
76: if (!w0) {
77: MatGetVecs(A,PETSC_NULL,&w0);
78: SlepcVecSetRandom(w0);
79: eps->vec_initial_left = w0;
80: }
82: VecGetSize(eps->vec_initial,&N);
83: if (eps->nev > N) eps->nev = N;
84: if (eps->ncv > N) eps->ncv = N;
85: if (!eps->tol) eps->tol = 1.e-7;
87: (*eps->ops->setup)(eps);
88: STSetUp(eps->OP);
90: /* DSV is equal to the columns of DS followed by the ones in V */
91: PetscFree(eps->DSV);
92: PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);
93: for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
94: for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
95:
96: if (eps->nds>0) {
97: if (!eps->ds_ortho) {
98: /* orthonormalize vectors in DS if necessary */
99: EPSQRDecomposition(eps,eps->DS,0,eps->nds,PETSC_NULL,0);
100: }
101: EPSOrthogonalize(eps,eps->nds,PETSC_NULL,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL);
102: }
104: STCheckNullSpace(eps->OP,eps->nds,eps->DS);
105:
106: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
107: eps->setupcalled = 1;
108: return(0);
109: }
113: /*@
114: EPSSetInitialVector - Sets the initial vector from which the
115: eigensolver starts to iterate.
117: Collective on EPS and Vec
119: Input Parameters:
120: + eps - the eigensolver context
121: - vec - the vector
123: Level: intermediate
125: .seealso: EPSGetInitialVector(), EPSSetLeftInitialVector()
127: @*/
128: PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
129: {
131:
136: if (eps->vec_initial) {
137: VecDestroy(eps->vec_initial);
138: }
139: eps->vec_initial = vec;
140: PetscObjectReference((PetscObject)eps->vec_initial);
141: return(0);
142: }
146: /*@
147: EPSGetInitialVector - Gets the initial vector associated with the
148: eigensolver; if the vector was not set it will return a 0 pointer or
149: a vector randomly generated by EPSSetUp().
151: Not collective, but vector is shared by all processors that share the EPS
153: Input Parameter:
154: . eps - the eigensolver context
156: Output Parameter:
157: . vec - the vector
159: Level: intermediate
161: .seealso: EPSSetInitialVector(), EPSGetLeftInitialVector()
163: @*/
164: PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
165: {
168: *vec = eps->vec_initial;
169: return(0);
170: }
174: /*@
175: EPSSetLeftInitialVector - Sets the initial vector from which the eigensolver
176: starts to iterate, corresponding to the left recurrence (two-sided solvers).
178: Collective on EPS and Vec
180: Input Parameters:
181: + eps - the eigensolver context
182: - vec - the vector
184: Level: intermediate
186: .seealso: EPSGetLeftInitialVector(), EPSSetInitialVector()
188: @*/
189: PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
190: {
192:
197: if (eps->vec_initial_left) {
198: VecDestroy(eps->vec_initial_left);
199: }
200: eps->vec_initial_left = vec;
201: PetscObjectReference((PetscObject)eps->vec_initial_left);
202: return(0);
203: }
207: /*@
208: EPSGetLeftInitialVector - Gets the left initial vector associated with the
209: eigensolver; if the vector was not set it will return a 0 pointer or
210: a vector randomly generated by EPSSetUp().
212: Not collective, but vector is shared by all processors that share the EPS
214: Input Parameter:
215: . eps - the eigensolver context
217: Output Parameter:
218: . vec - the vector
220: Level: intermediate
222: .seealso: EPSSetLeftInitialVector(), EPSGetLeftInitialVector()
224: @*/
225: PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
226: {
229: *vec = eps->vec_initial_left;
230: return(0);
231: }
235: /*@
236: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
238: Collective on EPS and Mat
240: Input Parameters:
241: + eps - the eigenproblem solver context
242: . A - the matrix associated with the eigensystem
243: - B - the second matrix in the case of generalized eigenproblems
245: Notes:
246: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
248: Level: beginner
250: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
251: @*/
252: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
253: {
255: PetscInt m,n;
264: /* Check for square matrices */
265: MatGetSize(A,&m,&n);
266: if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
267: if (B) {
268: MatGetSize(B,&m,&n);
269: if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
270: }
272: STSetOperators(eps->OP,A,B);
273: eps->setupcalled = 0; /* so that next solve call will call setup */
275: /* Destroy randomly generated initial vectors */
276: if (eps->vec_initial) {
277: VecDestroy(eps->vec_initial);
278: eps->vec_initial = PETSC_NULL;
279: }
280: if (eps->vec_initial_left) {
281: VecDestroy(eps->vec_initial_left);
282: eps->vec_initial_left = PETSC_NULL;
283: }
285: return(0);
286: }
290: /*@
291: EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.
293: Not Collective
295: Input Parameter:
296: + eps - the eigenproblem solver context
297: . n - number of vectors to add
298: . ds - set of basis vectors of the deflation space
299: - ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal
301: Notes:
302: When a deflation space is given, the eigensolver seeks the eigensolution
303: in the restriction of the problem to the orthogonal complement of this
304: space. This can be used for instance in the case that an invariant
305: subspace is known beforehand (such as the nullspace of the matrix).
307: The basis vectors can be provided all at once or incrementally with
308: several calls to EPSAttachDeflationSpace().
310: Use a value of PETSC_TRUE for parameter ortho if all the vectors passed
311: in are known to be mutually orthonormal.
313: Level: intermediate
315: .seealso: EPSRemoveDeflationSpace()
316: @*/
317: PetscErrorCode EPSAttachDeflationSpace(EPS eps,int n,Vec *ds,PetscTruth ortho)
318: {
320: int i;
321: Vec *tvec;
322:
325: tvec = eps->DS;
326: if (n+eps->nds > 0) {
327: PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);
328: }
329: if (eps->nds > 0) {
330: for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
331: PetscFree(tvec);
332: }
333: for (i=0; i<n; i++) {
334: VecDuplicate(ds[i],&eps->DS[i + eps->nds]);
335: VecCopy(ds[i],eps->DS[i + eps->nds]);
336: }
337: eps->nds += n;
338: if (!ortho) eps->ds_ortho = PETSC_FALSE;
339: eps->setupcalled = 0;
340: return(0);
341: }
345: /*@
346: EPSRemoveDeflationSpace - Removes the deflation space.
348: Not Collective
350: Input Parameter:
351: . eps - the eigenproblem solver context
353: Level: intermediate
355: .seealso: EPSAttachDeflationSpace()
356: @*/
357: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
358: {
360:
363: if (eps->nds > 0) {
364: VecDestroyVecs(eps->DS, eps->nds);
365: }
366: eps->ds_ortho = PETSC_TRUE;
367: eps->setupcalled = 0;
368: return(0);
369: }