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-2010, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include private/epsimpl.h
28: /*@
29: EPSSetUp - Sets up all the internal data structures necessary for the
30: execution of the eigensolver. Then calls STSetUp() for any set-up
31: operations associated to the ST object.
33: Collective on EPS
35: Input Parameter:
36: . eps - eigenproblem solver context
38: Notes:
39: This function need not be called explicitly in most cases, since EPSSolve()
40: calls it. It can be useful when one wants to measure the set-up time
41: separately from the solve time.
43: Level: advanced
45: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
46: @*/
47: PetscErrorCode EPSSetUp(EPS eps)
48: {
50: Vec vds;
51: Mat A,B;
52: PetscInt i,k;
53: PetscTruth flg,lindep;
54: PetscScalar *pDS;
55: PetscReal norm;
56: #if defined(PETSC_USE_COMPLEX)
57: PetscScalar sigma;
58: #endif
59:
63: if (eps->setupcalled) return(0);
65: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
67: /* Set default solver type */
68: if (!((PetscObject)eps)->type_name) {
69: EPSSetType(eps,EPSKRYLOVSCHUR);
70: }
71: if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
72:
73: /* Set problem dimensions */
74: STGetOperators(eps->OP,&A,&B);
75: MatGetSize(A,&eps->n,PETSC_NULL);
76: MatGetLocalSize(A,&eps->nloc,PETSC_NULL);
78: /* Set default problem type */
79: if (!eps->problem_type) {
80: if (B==PETSC_NULL) {
81: EPSSetProblemType(eps,EPS_NHEP);
82: }
83: else {
84: EPSSetProblemType(eps,EPS_GNHEP);
85: }
86: } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
87: SETERRQ(0,"Warning: Inconsistent EPS state");
88: }
89: #if defined(PETSC_USE_COMPLEX)
90: STGetShift(eps->OP,&sigma);
91: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) {
92: SETERRQ(1,"Hermitian problems are not compatible with complex shifts")
93: }
94: #endif
95: if (eps->ishermitian && eps->leftvecs) {
96: SETERRQ(1,"Requesting left eigenvectors not allowed in Hermitian problems")
97: }
98:
99: if (eps->ispositive) {
100: STGetBilinearForm(eps->OP,&B);
101: IPSetBilinearForm(eps->ip,B,IP_INNER_HERMITIAN);
102: MatDestroy(B);
103: } else {
104: IPSetBilinearForm(eps->ip,PETSC_NULL,IP_INNER_HERMITIAN);
105: }
106:
107: if (eps->nev > eps->n) eps->nev = eps->n;
108: if (eps->ncv > eps->n) eps->ncv = eps->n;
110: /* initialize the random number generator */
111: PetscRandomCreate(((PetscObject)eps)->comm,&eps->rand);
112: PetscRandomSetFromOptions(eps->rand);
114: /* initialization of matrix norms */
115: if (eps->nrma == PETSC_DETERMINE) {
116: MatHasOperation(A,MATOP_NORM,&flg);
117: if (flg) { MatNorm(A,NORM_INFINITY,&eps->nrma); }
118: else eps->nrma = 1.0;
119: }
120: if (eps->nrmb == PETSC_DETERMINE) {
121: MatHasOperation(B,MATOP_NORM,&flg);
122: if (flg) { MatNorm(B,NORM_INFINITY,&eps->nrmb); }
123: else eps->nrmb = 1.0;
124: }
126: /* call specific solver setup */
127: (*eps->ops->setup)(eps);
128: STSetUp(eps->OP);
129:
130: PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&flg);
131: if (flg && eps->problem_type == EPS_PGNHEP) {
132: SETERRQ(PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
133: }
135: PetscTypeCompare((PetscObject)eps->OP,STFOLD,&flg);
136: if (flg && !eps->ishermitian) {
137: SETERRQ(PETSC_ERR_SUP,"Fold spectral transformation requires a Hermitian problem");
138: }
140: if (eps->nds>0) {
141: if (!eps->ds_ortho) {
142: /* allocate memory and copy deflation basis vectors into DS */
143: PetscMalloc(eps->nds*eps->nloc*sizeof(PetscScalar),&pDS);
144: for (i=0;i<eps->nds;i++) {
145: VecCreateMPIWithArray(((PetscObject)eps)->comm,eps->nloc,PETSC_DECIDE,pDS+i*eps->nloc,&vds);
146: VecCopy(eps->DS[i],vds);
147: VecDestroy(eps->DS[i]);
148: eps->DS[i] = vds;
149: }
150: /* orthonormalize vectors in DS */
151: k = 0;
152: for (i=0;i<eps->nds;i++) {
153: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->DS,eps->DS[k],PETSC_NULL,&norm,&lindep);
154: if (norm==0.0 || lindep) PetscInfo(eps,"Linearly dependent deflation vector found, removing...\n");
155: else {
156: VecScale(eps->DS[k],1.0/norm);
157: k++;
158: }
159: }
160: eps->nds = k;
161: eps->ds_ortho = PETSC_TRUE;
162: }
163: }
164: STCheckNullSpace(eps->OP,eps->nds,eps->DS);
166: /* process initial vectors */
167: if (eps->nini<0) {
168: eps->nini = -eps->nini;
169: if (eps->nini>eps->ncv) SETERRQ(1,"The number of initial vectors is larger than ncv")
170: k = 0;
171: for (i=0;i<eps->nini;i++) {
172: VecCopy(eps->IS[i],eps->V[k]);
173: VecDestroy(eps->IS[i]);
174: IPOrthogonalize(eps->ip,eps->nds,eps->DS,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&lindep);
175: if (norm==0.0 || lindep) PetscInfo(eps,"Linearly dependent initial vector found, removing...\n");
176: else {
177: VecScale(eps->V[k],1.0/norm);
178: k++;
179: }
180: }
181: eps->nini = k;
182: PetscFree(eps->IS);
183: }
184: if (eps->ninil<0) {
185: if (!eps->leftvecs) PetscInfo(eps,"Ignoring initial left vectors\n");
186: else {
187: eps->ninil = -eps->ninil;
188: if (eps->ninil>eps->ncv) SETERRQ(1,"The number of initial left vectors is larger than ncv")
189: k = 0;
190: for (i=0;i<eps->ninil;i++) {
191: VecCopy(eps->ISL[i],eps->W[k]);
192: VecDestroy(eps->ISL[i]);
193: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->W,eps->W[k],PETSC_NULL,&norm,&lindep);
194: if (norm==0.0 || lindep) PetscInfo(eps,"Linearly dependent initial left vector found, removing...\n");
195: else {
196: VecScale(eps->W[k],1.0/norm);
197: k++;
198: }
199: }
200: eps->ninil = k;
201: PetscFree(eps->ISL);
202: }
203: }
205: /* Build balancing matrix if required */
206: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
207: if (!eps->D) {
208: VecDuplicate(eps->V[0],&eps->D);
209: }
210: else {
211: VecSet(eps->D,1.0);
212: }
213: EPSBuildBalance_Krylov(eps);
214: STSetBalanceMatrix(eps->OP,eps->D);
215: }
217: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
218: eps->setupcalled = 1;
219: return(0);
220: }
224: /*@
225: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
227: Collective on EPS and Mat
229: Input Parameters:
230: + eps - the eigenproblem solver context
231: . A - the matrix associated with the eigensystem
232: - B - the second matrix in the case of generalized eigenproblems
234: Notes:
235: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
237: Level: beginner
239: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
240: @*/
241: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
242: {
244: PetscInt m,n,m0;
253: /* Check for square matrices */
254: MatGetSize(A,&m,&n);
255: if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
256: if (B) {
257: MatGetSize(B,&m0,&n);
258: if (m0!=n) { SETERRQ(1,"B is a non-square matrix"); }
259: if (m!=m0) { SETERRQ(1,"Dimensions of A and B do not match"); }
260: }
262: STSetOperators(eps->OP,A,B);
263: eps->setupcalled = 0; /* so that next solve call will call setup */
265: if (eps->D) {
266: VecDestroy(eps->D);
267: eps->D = PETSC_NULL;
268: }
270: return(0);
271: }
275: /*@
276: EPSGetOperators - Gets the matrices associated with the eigensystem.
278: Collective on EPS and Mat
280: Input Parameter:
281: . eps - the EPS context
283: Output Parameters:
284: + A - the matrix associated with the eigensystem
285: - B - the second matrix in the case of generalized eigenproblems
287: Level: intermediate
289: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
290: @*/
291: PetscErrorCode EPSGetOperators(EPS eps, Mat *A, Mat *B)
292: {
294: ST st;
300: EPSGetST(eps,&st);
301: STGetOperators(st,A,B);
302: return(0);
303: }
307: /*@
308: EPSSetDeflationSpace - Specify a basis of vectors that constitute
309: the deflation space.
311: Collective on EPS and Vec
313: Input Parameter:
314: + eps - the eigenproblem solver context
315: . n - number of vectors
316: - ds - set of basis vectors of the deflation space
318: Notes:
319: When a deflation space is given, the eigensolver seeks the eigensolution
320: in the restriction of the problem to the orthogonal complement of this
321: space. This can be used for instance in the case that an invariant
322: subspace is known beforehand (such as the nullspace of the matrix).
324: Basis vectors set by a previous call to EPSSetDeflationSpace() are
325: replaced.
327: The vectors do not need to be mutually orthonormal, since they are explicitly
328: orthonormalized internally.
330: These vectors persist from one EPSSolve() call to the other, use
331: EPSRemoveDeflationSpace() to eliminate them.
333: Level: intermediate
335: .seealso: EPSRemoveDeflationSpace()
336: @*/
337: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *ds)
338: {
340: PetscInt i;
341:
344: if (n<=0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
346: /* free previous vectors */
347: EPSRemoveDeflationSpace(eps);
349: /* get references of passed vectors */
350: PetscMalloc(n*sizeof(Vec),&eps->DS);
351: for (i=0;i<n;i++) {
352: PetscObjectReference((PetscObject)ds[i]);
353: eps->DS[i] = ds[i];
354: }
356: eps->nds = n;
357: eps->setupcalled = 0;
358: eps->ds_ortho = PETSC_FALSE;
359: return(0);
360: }
364: /*@
365: EPSRemoveDeflationSpace - Removes the deflation space.
367: Collective on EPS
369: Input Parameter:
370: . eps - the eigenproblem solver context
372: Level: intermediate
374: .seealso: EPSSetDeflationSpace()
375: @*/
376: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
377: {
379: PetscInt i;
380: PetscScalar *pV;
381:
384: if (eps->nds > 0) {
385: VecGetArray(eps->DS[0],&pV);
386: VecRestoreArray(eps->DS[0],PETSC_NULL);
387: for (i=0;i<eps->nds;i++) {
388: VecDestroy(eps->DS[i]);
389: }
390: if (eps->setupcalled) { /* before EPSSetUp, DS are just references */
391: PetscFree(pV);
392: }
393: PetscFree(eps->DS);
394: }
395: eps->nds = 0;
396: eps->setupcalled = 0;
397: eps->ds_ortho = PETSC_FALSE;
398: return(0);
399: }
403: /*@
404: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
405: space, that is, the subspace from which the solver starts to iterate.
407: Collective on EPS and Vec
409: Input Parameter:
410: + eps - the eigenproblem solver context
411: . n - number of vectors
412: - is - set of basis vectors of the initial space
414: Notes:
415: Some solvers start to iterate on a single vector (initial vector). In that case,
416: the other vectors are ignored.
418: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
419: EPSSolve() call to the other, so the initial space should be set every time.
421: The vectors do not need to be mutually orthonormal, since they are explicitly
422: orthonormalized internally.
424: Common usage of this function is when the user can provide a rough approximation
425: of the wanted eigenspace. Then, convergence may be faster.
427: Level: intermediate
429: .seealso: EPSSetInitialSpaceLeft(), EPSSetDeflationSpace()
430: @*/
431: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
432: {
434: PetscInt i;
435:
438: if (n<0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
440: /* free previous non-processed vectors */
441: if (eps->nini<0) {
442: for (i=0;i<-eps->nini;i++) {
443: VecDestroy(eps->IS[i]);
444: }
445: PetscFree(eps->IS);
446: }
448: /* get references of passed vectors */
449: PetscMalloc(n*sizeof(Vec),&eps->IS);
450: for (i=0;i<n;i++) {
451: PetscObjectReference((PetscObject)is[i]);
452: eps->IS[i] = is[i];
453: }
455: eps->nini = -n;
456: eps->setupcalled = 0;
457: return(0);
458: }
462: /*@
463: EPSSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
464: left space, that is, the subspace from which the solver starts to iterate for
465: building the left subspace (in methods that work with two subspaces).
467: Collective on EPS and Vec
469: Input Parameter:
470: + eps - the eigenproblem solver context
471: . n - number of vectors
472: - is - set of basis vectors of the initial left space
474: Notes:
475: Some solvers start to iterate on a single vector (initial left vector). In that case,
476: the other vectors are ignored.
478: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
479: EPSSolve() call to the other, so the initial left space should be set every time.
481: The vectors do not need to be mutually orthonormal, since they are explicitly
482: orthonormalized internally.
484: Common usage of this function is when the user can provide a rough approximation
485: of the wanted left eigenspace. Then, convergence may be faster.
487: Level: intermediate
489: .seealso: EPSSetInitialSpace(), EPSSetDeflationSpace()
490: @*/
491: PetscErrorCode EPSSetInitialSpaceLeft(EPS eps,PetscInt n,Vec *is)
492: {
494: PetscInt i;
495:
498: if (n<0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
500: /* free previous non-processed vectors */
501: if (eps->ninil<0) {
502: for (i=0;i<-eps->ninil;i++) {
503: VecDestroy(eps->ISL[i]);
504: }
505: PetscFree(eps->ISL);
506: }
508: /* get references of passed vectors */
509: PetscMalloc(n*sizeof(Vec),&eps->ISL);
510: for (i=0;i<n;i++) {
511: PetscObjectReference((PetscObject)is[i]);
512: eps->ISL[i] = is[i];
513: }
515: eps->ninil = -n;
516: eps->setupcalled = 0;
517: return(0);
518: }