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-2011, Universitat 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> /*I "slepceps.h" I*/
25: #include <private/ipimpl.h>
29: /*@
30: EPSSetUp - Sets up all the internal data structures necessary for the
31: execution of the eigensolver. Then calls STSetUp() for any set-up
32: operations associated to the ST object.
34: Collective on EPS
36: Input Parameter:
37: . eps - eigenproblem solver context
39: Notes:
40: This function need not be called explicitly in most cases, since EPSSolve()
41: calls it. It can be useful when one wants to measure the set-up time
42: separately from the solve time.
44: Level: advanced
46: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
47: @*/
48: PetscErrorCode EPSSetUp(EPS eps)
49: {
51: Mat A,B;
52: PetscInt i,k;
53: PetscBool flg,lindep;
54: Vec *newDS;
55: PetscReal norm;
56: #if defined(PETSC_USE_COMPLEX)
57: PetscScalar sigma;
58: #endif
59:
62: if (eps->setupcalled) return(0);
63: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
65: /* Set default solver type (EPSSetFromOptions was not called) */
66: if (!((PetscObject)eps)->type_name) {
67: EPSSetType(eps,EPSKRYLOVSCHUR);
68: }
69: if (!eps->OP) { EPSGetST(eps,&eps->OP); }
70: if (!((PetscObject)eps->OP)->type_name) {
71: PetscTypeCompareAny((PetscObject)eps,&flg,EPSGD,EPSJD,"");
72: STSetType(eps->OP,flg?STPRECOND:STSHIFT);
73: }
74: if (!eps->ip) { EPSGetIP(eps,&eps->ip); }
75: if (!((PetscObject)eps->ip)->type_name) {
76: IPSetDefaultType_Private(eps->ip);
77: }
78: if (!((PetscObject)eps->rand)->type_name) {
79: PetscRandomSetFromOptions(eps->rand);
80: }
81:
82: /* Set problem dimensions */
83: STGetOperators(eps->OP,&A,&B);
84: if (!A) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
85: MatGetSize(A,&eps->n,PETSC_NULL);
86: MatGetLocalSize(A,&eps->nloc,PETSC_NULL);
87: VecDestroy(&eps->t);
88: SlepcMatGetVecsTemplate(A,&eps->t,PETSC_NULL);
90: /* Set default problem type */
91: if (!eps->problem_type) {
92: if (B==PETSC_NULL) {
93: EPSSetProblemType(eps,EPS_NHEP);
94: } else {
95: EPSSetProblemType(eps,EPS_GNHEP);
96: }
97: } else if (!B && eps->isgeneralized) {
98: PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
99: eps->isgeneralized = PETSC_FALSE;
100: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
101: } else if (B && !eps->isgeneralized) {
102: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");
103: }
104: #if defined(PETSC_USE_COMPLEX)
105: STGetShift(eps->OP,&sigma);
106: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0)
107: SETERRQ(((PetscObject)eps)->comm,1,"Hermitian problems are not compatible with complex shifts");
108: #endif
109: if (eps->ishermitian && eps->leftvecs)
110: SETERRQ(((PetscObject)eps)->comm,1,"Requesting left eigenvectors not allowed in Hermitian problems");
111:
112: if (eps->ispositive) {
113: STGetBilinearForm(eps->OP,&B);
114: IPSetMatrix(eps->ip,B);
115: MatDestroy(&B);
116: } else {
117: IPSetMatrix(eps->ip,PETSC_NULL);
118: }
119:
120: if (eps->nev > eps->n) eps->nev = eps->n;
121: if (eps->ncv > eps->n) eps->ncv = eps->n;
123: /* initialization of matrix norms */
124: if (eps->nrma == PETSC_DETERMINE) {
125: MatHasOperation(A,MATOP_NORM,&flg);
126: if (flg) { MatNorm(A,NORM_INFINITY,&eps->nrma); }
127: else eps->nrma = 1.0;
128: }
129: if (eps->nrmb == PETSC_DETERMINE) {
130: MatHasOperation(B,MATOP_NORM,&flg);
131: if (flg) { MatNorm(B,NORM_INFINITY,&eps->nrmb); }
132: else eps->nrmb = 1.0;
133: }
135: if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
137: /* call specific solver setup */
138: (*eps->ops->setup)(eps);
140: /* set tolerance if not yet set */
141: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
143: /* Build balancing matrix if required */
144: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
145: if (!eps->D) {
146: VecDuplicate(eps->V[0],&eps->D);
147: } else {
148: VecSet(eps->D,1.0);
149: }
150: EPSBuildBalance_Krylov(eps);
151: STSetBalanceMatrix(eps->OP,eps->D);
152: }
154: /* Setup ST */
155: STSetUp(eps->OP);
156:
157: PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&flg);
158: if (flg && eps->problem_type == EPS_PGNHEP)
159: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
161: PetscTypeCompare((PetscObject)eps->OP,STFOLD,&flg);
162: if (flg && !eps->ishermitian)
163: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Fold spectral transformation requires a Hermitian problem");
165: if (eps->nds>0) {
166: if (!eps->ds_ortho) {
167: /* allocate memory and copy deflation basis vectors into DS */
168: VecDuplicateVecs(eps->t,eps->nds,&newDS);
169: for (i=0;i<eps->nds;i++) {
170: VecCopy(eps->DS[i],newDS[i]);
171: VecDestroy(&eps->DS[i]);
172: }
173: PetscFree(eps->DS);
174: eps->DS = newDS;
175: /* orthonormalize vectors in DS */
176: k = 0;
177: for (i=0;i<eps->nds;i++) {
178: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->DS,eps->DS[k],PETSC_NULL,&norm,&lindep);
179: if (norm==0.0 || lindep) {
180: PetscInfo(eps,"Linearly dependent deflation vector found, removing...\n");
181: } else {
182: VecScale(eps->DS[k],1.0/norm);
183: k++;
184: }
185: }
186: for (i=k;i<eps->nds;i++) { VecDestroy(&eps->DS[i]); }
187: eps->nds = k;
188: eps->ds_ortho = PETSC_TRUE;
189: }
190: }
191: STCheckNullSpace(eps->OP,eps->nds,eps->DS);
193: /* process initial vectors */
194: if (eps->nini<0) {
195: eps->nini = -eps->nini;
196: if (eps->nini>eps->ncv) SETERRQ(((PetscObject)eps)->comm,1,"The number of initial vectors is larger than ncv");
197: k = 0;
198: for (i=0;i<eps->nini;i++) {
199: VecCopy(eps->IS[i],eps->V[k]);
200: VecDestroy(&eps->IS[i]);
201: IPOrthogonalize(eps->ip,eps->nds,eps->DS,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&lindep);
202: if (norm==0.0 || lindep) {
203: PetscInfo(eps,"Linearly dependent initial vector found, removing...\n");
204: } else {
205: VecScale(eps->V[k],1.0/norm);
206: k++;
207: }
208: }
209: eps->nini = k;
210: PetscFree(eps->IS);
211: }
212: if (eps->ninil<0) {
213: if (!eps->leftvecs) {
214: PetscInfo(eps,"Ignoring initial left vectors\n");
215: } else {
216: eps->ninil = -eps->ninil;
217: if (eps->ninil>eps->ncv) SETERRQ(((PetscObject)eps)->comm,1,"The number of initial left vectors is larger than ncv");
218: k = 0;
219: for (i=0;i<eps->ninil;i++) {
220: VecCopy(eps->ISL[i],eps->W[k]);
221: VecDestroy(&eps->ISL[i]);
222: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->W,eps->W[k],PETSC_NULL,&norm,&lindep);
223: if (norm==0.0 || lindep) {
224: PetscInfo(eps,"Linearly dependent initial left vector found, removing...\n");
225: } else {
226: VecScale(eps->W[k],1.0/norm);
227: k++;
228: }
229: }
230: eps->ninil = k;
231: PetscFree(eps->ISL);
232: }
233: }
235: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
236: eps->setupcalled = 1;
237: return(0);
238: }
242: /*@
243: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
245: Collective on EPS and Mat
247: Input Parameters:
248: + eps - the eigenproblem solver context
249: . A - the matrix associated with the eigensystem
250: - B - the second matrix in the case of generalized eigenproblems
252: Notes:
253: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
255: It must be called after EPSSetUp(). If it is called again after EPSSetUp() then
256: the EPS object is reset.
258: Level: beginner
260: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
261: @*/
262: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
263: {
265: PetscInt m,n,m0;
274: /* Check for square matrices */
275: MatGetSize(A,&m,&n);
276: if (m!=n) SETERRQ(((PetscObject)eps)->comm,1,"A is a non-square matrix");
277: if (B) {
278: MatGetSize(B,&m0,&n);
279: if (m0!=n) SETERRQ(((PetscObject)eps)->comm,1,"B is a non-square matrix");
280: if (m!=m0) SETERRQ(((PetscObject)eps)->comm,1,"Dimensions of A and B do not match");
281: }
283: if (eps->setupcalled) { EPSReset(eps); }
284: if (!eps->OP) { EPSGetST(eps,&eps->OP); }
285: STSetOperators(eps->OP,A,B);
286: return(0);
287: }
291: /*@
292: EPSGetOperators - Gets the matrices associated with the eigensystem.
294: Collective on EPS and Mat
296: Input Parameter:
297: . eps - the EPS context
299: Output Parameters:
300: + A - the matrix associated with the eigensystem
301: - B - the second matrix in the case of generalized eigenproblems
303: Level: intermediate
305: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
306: @*/
307: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
308: {
310: ST st;
316: EPSGetST(eps,&st);
317: STGetOperators(st,A,B);
318: return(0);
319: }
323: /*@
324: EPSSetDeflationSpace - Specify a basis of vectors that constitute
325: the deflation space.
327: Collective on EPS and Vec
329: Input Parameter:
330: + eps - the eigenproblem solver context
331: . n - number of vectors
332: - ds - set of basis vectors of the deflation space
334: Notes:
335: When a deflation space is given, the eigensolver seeks the eigensolution
336: in the restriction of the problem to the orthogonal complement of this
337: space. This can be used for instance in the case that an invariant
338: subspace is known beforehand (such as the nullspace of the matrix).
340: Basis vectors set by a previous call to EPSSetDeflationSpace() are
341: replaced.
343: The vectors do not need to be mutually orthonormal, since they are explicitly
344: orthonormalized internally.
346: These vectors persist from one EPSSolve() call to the other, use
347: EPSRemoveDeflationSpace() to eliminate them.
349: Level: intermediate
351: .seealso: EPSRemoveDeflationSpace()
352: @*/
353: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *ds)
354: {
356: PetscInt i;
357:
361: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
363: /* free previous vectors */
364: EPSRemoveDeflationSpace(eps);
366: /* get references of passed vectors */
367: if (n>0) {
368: PetscMalloc(n*sizeof(Vec),&eps->DS);
369: for (i=0;i<n;i++) {
370: PetscObjectReference((PetscObject)ds[i]);
371: eps->DS[i] = ds[i];
372: }
373: eps->setupcalled = 0;
374: eps->ds_ortho = PETSC_FALSE;
375: }
377: eps->nds = n;
378: return(0);
379: }
383: /*@
384: EPSRemoveDeflationSpace - Removes the deflation space.
386: Collective on EPS
388: Input Parameter:
389: . eps - the eigenproblem solver context
391: Level: intermediate
393: .seealso: EPSSetDeflationSpace()
394: @*/
395: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
396: {
398:
401: VecDestroyVecs(eps->nds,&eps->DS);
402: eps->nds = 0;
403: eps->setupcalled = 0;
404: eps->ds_ortho = PETSC_FALSE;
405: return(0);
406: }
410: /*@
411: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
412: space, that is, the subspace from which the solver starts to iterate.
414: Collective on EPS and Vec
416: Input Parameter:
417: + eps - the eigenproblem solver context
418: . n - number of vectors
419: - is - set of basis vectors of the initial space
421: Notes:
422: Some solvers start to iterate on a single vector (initial vector). In that case,
423: the other vectors are ignored.
425: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
426: EPSSolve() call to the other, so the initial space should be set every time.
428: The vectors do not need to be mutually orthonormal, since they are explicitly
429: orthonormalized internally.
431: Common usage of this function is when the user can provide a rough approximation
432: of the wanted eigenspace. Then, convergence may be faster.
434: Level: intermediate
436: .seealso: EPSSetInitialSpaceLeft(), EPSSetDeflationSpace()
437: @*/
438: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
439: {
441: PetscInt i;
442:
446: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
448: /* free previous non-processed vectors */
449: if (eps->nini<0) {
450: for (i=0;i<-eps->nini;i++) {
451: VecDestroy(&eps->IS[i]);
452: }
453: PetscFree(eps->IS);
454: }
456: /* get references of passed vectors */
457: if (n>0) {
458: PetscMalloc(n*sizeof(Vec),&eps->IS);
459: for (i=0;i<n;i++) {
460: PetscObjectReference((PetscObject)is[i]);
461: eps->IS[i] = is[i];
462: }
463: eps->setupcalled = 0;
464: }
466: eps->nini = -n;
467: return(0);
468: }
472: /*@
473: EPSSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
474: left space, that is, the subspace from which the solver starts to iterate for
475: building the left subspace (in methods that work with two subspaces).
477: Collective on EPS and Vec
479: Input Parameter:
480: + eps - the eigenproblem solver context
481: . n - number of vectors
482: - is - set of basis vectors of the initial left space
484: Notes:
485: Some solvers start to iterate on a single vector (initial left vector). In that case,
486: the other vectors are ignored.
488: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
489: EPSSolve() call to the other, so the initial left space should be set every time.
491: The vectors do not need to be mutually orthonormal, since they are explicitly
492: orthonormalized internally.
494: Common usage of this function is when the user can provide a rough approximation
495: of the wanted left eigenspace. Then, convergence may be faster.
497: Level: intermediate
499: .seealso: EPSSetInitialSpace(), EPSSetDeflationSpace()
500: @*/
501: PetscErrorCode EPSSetInitialSpaceLeft(EPS eps,PetscInt n,Vec *is)
502: {
504: PetscInt i;
505:
509: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
511: /* free previous non-processed vectors */
512: if (eps->ninil<0) {
513: for (i=0;i<-eps->ninil;i++) {
514: VecDestroy(&eps->ISL[i]);
515: }
516: PetscFree(eps->ISL);
517: }
519: /* get references of passed vectors */
520: if (n>0) {
521: PetscMalloc(n*sizeof(Vec),&eps->ISL);
522: for (i=0;i<n;i++) {
523: PetscObjectReference((PetscObject)is[i]);
524: eps->ISL[i] = is[i];
525: }
526: eps->setupcalled = 0;
527: }
529: eps->ninil = -n;
530: return(0);
531: }