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> /*I "slepcip.h" I*/
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: /* call specific solver setup */
136: (*eps->ops->setup)(eps);
138: /* set tolerance if not yet set */
139: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
141: /* Build balancing matrix if required */
142: if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
143: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
144: if (!eps->D) {
145: VecDuplicate(eps->V[0],&eps->D);
146: } else {
147: VecSet(eps->D,1.0);
148: }
149: EPSBuildBalance_Krylov(eps);
150: STSetBalanceMatrix(eps->OP,eps->D);
151: }
153: /* Setup ST */
154: STSetUp(eps->OP);
155:
156: PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&flg);
157: if (flg && eps->problem_type == EPS_PGNHEP)
158: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
160: PetscTypeCompare((PetscObject)eps->OP,STFOLD,&flg);
161: if (flg && !eps->ishermitian)
162: SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Fold spectral transformation requires a Hermitian problem");
164: if (eps->nds>0) {
165: if (!eps->ds_ortho) {
166: /* allocate memory and copy deflation basis vectors into DS */
167: VecDuplicateVecs(eps->t,eps->nds,&newDS);
168: for (i=0;i<eps->nds;i++) {
169: VecCopy(eps->DS[i],newDS[i]);
170: VecDestroy(&eps->DS[i]);
171: }
172: PetscFree(eps->DS);
173: eps->DS = newDS;
174: /* orthonormalize vectors in DS */
175: k = 0;
176: for (i=0;i<eps->nds;i++) {
177: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->DS,eps->DS[k],PETSC_NULL,&norm,&lindep);
178: if (norm==0.0 || lindep) {
179: PetscInfo(eps,"Linearly dependent deflation vector found, removing...\n");
180: } else {
181: VecScale(eps->DS[k],1.0/norm);
182: k++;
183: }
184: }
185: for (i=k;i<eps->nds;i++) { VecDestroy(&eps->DS[i]); }
186: eps->nds = k;
187: eps->ds_ortho = PETSC_TRUE;
188: }
189: }
190: STCheckNullSpace(eps->OP,eps->nds,eps->DS);
192: /* process initial vectors */
193: if (eps->nini<0) {
194: eps->nini = -eps->nini;
195: if (eps->nini>eps->ncv) SETERRQ(((PetscObject)eps)->comm,1,"The number of initial vectors is larger than ncv");
196: k = 0;
197: for (i=0;i<eps->nini;i++) {
198: VecCopy(eps->IS[i],eps->V[k]);
199: VecDestroy(&eps->IS[i]);
200: IPOrthogonalize(eps->ip,eps->nds,eps->DS,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&lindep);
201: if (norm==0.0 || lindep) {
202: PetscInfo(eps,"Linearly dependent initial vector found, removing...\n");
203: } else {
204: VecScale(eps->V[k],1.0/norm);
205: k++;
206: }
207: }
208: eps->nini = k;
209: PetscFree(eps->IS);
210: }
211: if (eps->ninil<0) {
212: if (!eps->leftvecs) {
213: PetscInfo(eps,"Ignoring initial left vectors\n");
214: } else {
215: eps->ninil = -eps->ninil;
216: if (eps->ninil>eps->ncv) SETERRQ(((PetscObject)eps)->comm,1,"The number of initial left vectors is larger than ncv");
217: k = 0;
218: for (i=0;i<eps->ninil;i++) {
219: VecCopy(eps->ISL[i],eps->W[k]);
220: VecDestroy(&eps->ISL[i]);
221: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->W,eps->W[k],PETSC_NULL,&norm,&lindep);
222: if (norm==0.0 || lindep) {
223: PetscInfo(eps,"Linearly dependent initial left vector found, removing...\n");
224: } else {
225: VecScale(eps->W[k],1.0/norm);
226: k++;
227: }
228: }
229: eps->ninil = k;
230: PetscFree(eps->ISL);
231: }
232: }
234: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
235: eps->setupcalled = 1;
236: return(0);
237: }
241: /*@
242: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
244: Collective on EPS and Mat
246: Input Parameters:
247: + eps - the eigenproblem solver context
248: . A - the matrix associated with the eigensystem
249: - B - the second matrix in the case of generalized eigenproblems
251: Notes:
252: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
254: It must be called after EPSSetUp(). If it is called again after EPSSetUp() then
255: the EPS object is reset.
257: Level: beginner
259: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
260: @*/
261: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
262: {
264: PetscInt m,n,m0;
273: /* Check for square matrices */
274: MatGetSize(A,&m,&n);
275: if (m!=n) SETERRQ(((PetscObject)eps)->comm,1,"A is a non-square matrix");
276: if (B) {
277: MatGetSize(B,&m0,&n);
278: if (m0!=n) SETERRQ(((PetscObject)eps)->comm,1,"B is a non-square matrix");
279: if (m!=m0) SETERRQ(((PetscObject)eps)->comm,1,"Dimensions of A and B do not match");
280: }
282: if (eps->setupcalled) { EPSReset(eps); }
283: if (!eps->OP) { EPSGetST(eps,&eps->OP); }
284: STSetOperators(eps->OP,A,B);
285: return(0);
286: }
290: /*@
291: EPSGetOperators - Gets the matrices associated with the eigensystem.
293: Collective on EPS and Mat
295: Input Parameter:
296: . eps - the EPS context
298: Output Parameters:
299: + A - the matrix associated with the eigensystem
300: - B - the second matrix in the case of generalized eigenproblems
302: Level: intermediate
304: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
305: @*/
306: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
307: {
309: ST st;
315: EPSGetST(eps,&st);
316: STGetOperators(st,A,B);
317: return(0);
318: }
322: /*@
323: EPSSetDeflationSpace - Specify a basis of vectors that constitute
324: the deflation space.
326: Collective on EPS and Vec
328: Input Parameter:
329: + eps - the eigenproblem solver context
330: . n - number of vectors
331: - ds - set of basis vectors of the deflation space
333: Notes:
334: When a deflation space is given, the eigensolver seeks the eigensolution
335: in the restriction of the problem to the orthogonal complement of this
336: space. This can be used for instance in the case that an invariant
337: subspace is known beforehand (such as the nullspace of the matrix).
339: Basis vectors set by a previous call to EPSSetDeflationSpace() are
340: replaced.
342: The vectors do not need to be mutually orthonormal, since they are explicitly
343: orthonormalized internally.
345: These vectors persist from one EPSSolve() call to the other, use
346: EPSRemoveDeflationSpace() to eliminate them.
348: Level: intermediate
350: .seealso: EPSRemoveDeflationSpace()
351: @*/
352: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *ds)
353: {
355: PetscInt i;
356:
360: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
362: /* free previous vectors */
363: EPSRemoveDeflationSpace(eps);
365: /* get references of passed vectors */
366: if (n>0) {
367: PetscMalloc(n*sizeof(Vec),&eps->DS);
368: for (i=0;i<n;i++) {
369: PetscObjectReference((PetscObject)ds[i]);
370: eps->DS[i] = ds[i];
371: }
372: eps->setupcalled = 0;
373: eps->ds_ortho = PETSC_FALSE;
374: }
376: eps->nds = n;
377: return(0);
378: }
382: /*@
383: EPSRemoveDeflationSpace - Removes the deflation space.
385: Collective on EPS
387: Input Parameter:
388: . eps - the eigenproblem solver context
390: Level: intermediate
392: .seealso: EPSSetDeflationSpace()
393: @*/
394: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
395: {
397:
400: VecDestroyVecs(eps->nds,&eps->DS);
401: eps->nds = 0;
402: eps->setupcalled = 0;
403: eps->ds_ortho = PETSC_FALSE;
404: return(0);
405: }
409: /*@
410: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
411: space, that is, the subspace from which the solver starts to iterate.
413: Collective on EPS and Vec
415: Input Parameter:
416: + eps - the eigenproblem solver context
417: . n - number of vectors
418: - is - set of basis vectors of the initial space
420: Notes:
421: Some solvers start to iterate on a single vector (initial vector). In that case,
422: the other vectors are ignored.
424: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
425: EPSSolve() call to the other, so the initial space should be set every time.
427: The vectors do not need to be mutually orthonormal, since they are explicitly
428: orthonormalized internally.
430: Common usage of this function is when the user can provide a rough approximation
431: of the wanted eigenspace. Then, convergence may be faster.
433: Level: intermediate
435: .seealso: EPSSetInitialSpaceLeft(), EPSSetDeflationSpace()
436: @*/
437: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
438: {
440: PetscInt i;
441:
445: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
447: /* free previous non-processed vectors */
448: if (eps->nini<0) {
449: for (i=0;i<-eps->nini;i++) {
450: VecDestroy(&eps->IS[i]);
451: }
452: PetscFree(eps->IS);
453: }
455: /* get references of passed vectors */
456: if (n>0) {
457: PetscMalloc(n*sizeof(Vec),&eps->IS);
458: for (i=0;i<n;i++) {
459: PetscObjectReference((PetscObject)is[i]);
460: eps->IS[i] = is[i];
461: }
462: eps->setupcalled = 0;
463: }
465: eps->nini = -n;
466: return(0);
467: }
471: /*@
472: EPSSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
473: left space, that is, the subspace from which the solver starts to iterate for
474: building the left subspace (in methods that work with two subspaces).
476: Collective on EPS and Vec
478: Input Parameter:
479: + eps - the eigenproblem solver context
480: . n - number of vectors
481: - is - set of basis vectors of the initial left space
483: Notes:
484: Some solvers start to iterate on a single vector (initial left vector). In that case,
485: the other vectors are ignored.
487: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
488: EPSSolve() call to the other, so the initial left space should be set every time.
490: The vectors do not need to be mutually orthonormal, since they are explicitly
491: orthonormalized internally.
493: Common usage of this function is when the user can provide a rough approximation
494: of the wanted left eigenspace. Then, convergence may be faster.
496: Level: intermediate
498: .seealso: EPSSetInitialSpace(), EPSSetDeflationSpace()
499: @*/
500: PetscErrorCode EPSSetInitialSpaceLeft(EPS eps,PetscInt n,Vec *is)
501: {
503: PetscInt i;
504:
508: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
510: /* free previous non-processed vectors */
511: if (eps->ninil<0) {
512: for (i=0;i<-eps->ninil;i++) {
513: VecDestroy(&eps->ISL[i]);
514: }
515: PetscFree(eps->ISL);
516: }
518: /* get references of passed vectors */
519: if (n>0) {
520: PetscMalloc(n*sizeof(Vec),&eps->ISL);
521: for (i=0;i<n;i++) {
522: PetscObjectReference((PetscObject)is[i]);
523: eps->ISL[i] = is[i];
524: }
525: eps->setupcalled = 0;
526: }
528: eps->ninil = -n;
529: return(0);
530: }