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-2013, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
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 <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
25: #include <slepc-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,nmat;
53: PetscBool flg,lindep;
54: Vec *newDS;
55: PetscReal norm;
56: #if defined(PETSC_USE_COMPLEX)
57: PetscScalar sigma;
58: #endif
62: if (eps->setupcalled) return(0);
63: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
65: /* reset the convergence flag from the previous solves */
66: eps->reason = EPS_CONVERGED_ITERATING;
68: /* Set default solver type (EPSSetFromOptions was not called) */
69: if (!((PetscObject)eps)->type_name) {
70: EPSSetType(eps,EPSKRYLOVSCHUR);
71: }
72: if (!eps->st) { EPSGetST(eps,&eps->st); }
73: if (!((PetscObject)eps->st)->type_name) {
74: PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSGD,EPSJD,EPSRQCG,EPSBLOPEX,EPSPRIMME,"");
75: STSetType(eps->st,flg?STPRECOND:STSHIFT);
76: }
77: if (!eps->ip) { EPSGetIP(eps,&eps->ip); }
78: if (!((PetscObject)eps->ip)->type_name) {
79: IPSetType_Default(eps->ip);
80: }
81: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
82: DSReset(eps->ds);
83: if (!((PetscObject)eps->rand)->type_name) {
84: PetscRandomSetFromOptions(eps->rand);
85: }
87: /* Set problem dimensions */
88: STGetNumMatrices(eps->st,&nmat);
89: if (!nmat) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
90: STGetOperators(eps->st,0,&A);
91: MatGetSize(A,&eps->n,NULL);
92: MatGetLocalSize(A,&eps->nloc,NULL);
93: VecDestroy(&eps->t);
94: SlepcMatGetVecsTemplate(A,&eps->t,NULL);
95: PetscLogObjectParent(eps,eps->t);
97: /* Set default problem type */
98: if (!eps->problem_type) {
99: if (nmat==1) {
100: EPSSetProblemType(eps,EPS_NHEP);
101: } else {
102: EPSSetProblemType(eps,EPS_GNHEP);
103: }
104: } else if (nmat==1 && eps->isgeneralized) {
105: PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
106: eps->isgeneralized = PETSC_FALSE;
107: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
108: } else if (nmat>1 && !eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");
109: #if defined(PETSC_USE_COMPLEX)
110: STGetShift(eps->st,&sigma);
111: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
112: #endif
113: if (eps->ishermitian && eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Requesting left eigenvectors not allowed in Hermitian problems");
115: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
116: STGetBilinearForm(eps->st,&B);
117: IPSetMatrix(eps->ip,B);
118: MatDestroy(&B);
119: if (!eps->ispositive) {
120: IPSetType(eps->ip,IPINDEFINITE);
121: }
122: } else {
123: IPSetMatrix(eps->ip,NULL);
124: }
126: if (eps->nev > eps->n) eps->nev = eps->n;
127: if (eps->ncv > eps->n) eps->ncv = eps->n;
129: /* initialization of matrix norms */
130: if (eps->nrma == PETSC_DETERMINE) {
131: MatHasOperation(A,MATOP_NORM,&flg);
132: if (flg) {
133: MatNorm(A,NORM_INFINITY,&eps->nrma);
134: } else eps->nrma = 1.0;
135: }
136: if (eps->nrmb == PETSC_DETERMINE) {
137: if (nmat>1) { STGetOperators(eps->st,1,&B); }
138: MatHasOperation(B,MATOP_NORM,&flg);
139: if (flg) {
140: MatNorm(B,NORM_INFINITY,&eps->nrmb);
141: } else eps->nrmb = 1.0;
142: }
144: if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
146: /* call specific solver setup */
147: (*eps->ops->setup)(eps);
149: /* check extraction */
150: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STPRECOND,STSHIFT,"");
151: if (!flg && eps->extraction && eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
153: /* set tolerance if not yet set */
154: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
156: /* set eigenvalue comparison */
157: switch (eps->which) {
158: case EPS_LARGEST_MAGNITUDE:
159: eps->comparison = SlepcCompareLargestMagnitude;
160: eps->comparisonctx = NULL;
161: break;
162: case EPS_SMALLEST_MAGNITUDE:
163: eps->comparison = SlepcCompareSmallestMagnitude;
164: eps->comparisonctx = NULL;
165: break;
166: case EPS_LARGEST_REAL:
167: eps->comparison = SlepcCompareLargestReal;
168: eps->comparisonctx = NULL;
169: break;
170: case EPS_SMALLEST_REAL:
171: eps->comparison = SlepcCompareSmallestReal;
172: eps->comparisonctx = NULL;
173: break;
174: case EPS_LARGEST_IMAGINARY:
175: eps->comparison = SlepcCompareLargestImaginary;
176: eps->comparisonctx = NULL;
177: break;
178: case EPS_SMALLEST_IMAGINARY:
179: eps->comparison = SlepcCompareSmallestImaginary;
180: eps->comparisonctx = NULL;
181: break;
182: case EPS_TARGET_MAGNITUDE:
183: eps->comparison = SlepcCompareTargetMagnitude;
184: eps->comparisonctx = &eps->target;
185: break;
186: case EPS_TARGET_REAL:
187: eps->comparison = SlepcCompareTargetReal;
188: eps->comparisonctx = &eps->target;
189: break;
190: case EPS_TARGET_IMAGINARY:
191: eps->comparison = SlepcCompareTargetImaginary;
192: eps->comparisonctx = &eps->target;
193: break;
194: case EPS_ALL:
195: eps->comparison = SlepcCompareSmallestReal;
196: eps->comparisonctx = NULL;
197: break;
198: case EPS_WHICH_USER:
199: break;
200: }
202: /* Build balancing matrix if required */
203: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
204: if (!eps->D) {
205: VecDuplicate(eps->V[0],&eps->D);
206: PetscLogObjectParent(eps,eps->D);
207: } else {
208: VecSet(eps->D,1.0);
209: }
210: EPSBuildBalance_Krylov(eps);
211: STSetBalanceMatrix(eps->st,eps->D);
212: }
214: /* Setup ST */
215: STSetUp(eps->st);
217: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
218: if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
220: PetscObjectTypeCompare((PetscObject)eps->st,STFOLD,&flg);
221: if (flg && !eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Fold spectral transformation requires a Hermitian problem");
223: if (eps->nds>0) {
224: if (!eps->ds_ortho) {
225: /* allocate memory and copy deflation basis vectors into defl */
226: VecDuplicateVecs(eps->t,eps->nds,&newDS);
227: for (i=0;i<eps->nds;i++) {
228: VecCopy(eps->defl[i],newDS[i]);
229: VecDestroy(&eps->defl[i]);
230: }
231: PetscFree(eps->defl);
232: eps->defl = newDS;
233: PetscLogObjectParents(eps,eps->nds,eps->defl);
234: /* orthonormalize vectors in defl */
235: k = 0;
236: for (i=0;i<eps->nds;i++) {
237: IPOrthogonalize(eps->ip,0,NULL,k,NULL,eps->defl,eps->defl[k],NULL,&norm,&lindep);
238: if (norm==0.0 || lindep) {
239: PetscInfo(eps,"Linearly dependent deflation vector found, removing...\n");
240: } else {
241: VecScale(eps->defl[k],1.0/norm);
242: k++;
243: }
244: }
245: for (i=k;i<eps->nds;i++) { VecDestroy(&eps->defl[i]); }
246: eps->nds = k;
247: eps->ds_ortho = PETSC_TRUE;
248: }
249: }
250: STCheckNullSpace(eps->st,eps->nds,eps->defl);
252: /* process initial vectors */
253: if (eps->nini<0) {
254: eps->nini = -eps->nini;
255: if (eps->nini>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The number of initial vectors is larger than ncv");
256: IPOrthonormalizeBasis_Private(eps->ip,&eps->nini,&eps->IS,eps->V);
257: }
258: if (eps->ninil<0) {
259: if (!eps->leftvecs) {
260: PetscInfo(eps,"Ignoring initial left vectors\n");
261: } else {
262: eps->ninil = -eps->ninil;
263: if (eps->ninil>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The number of initial left vectors is larger than ncv");
264: IPOrthonormalizeBasis_Private(eps->ip,&eps->ninil,&eps->ISL,eps->W);
265: }
266: }
268: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
269: eps->setupcalled = 1;
270: return(0);
271: }
275: /*@
276: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
278: Collective on EPS and Mat
280: Input Parameters:
281: + eps - the eigenproblem solver context
282: . A - the matrix associated with the eigensystem
283: - B - the second matrix in the case of generalized eigenproblems
285: Notes:
286: To specify a standard eigenproblem, use NULL for parameter B.
288: It must be called after EPSSetUp(). If it is called again after EPSSetUp() then
289: the EPS object is reset.
291: Level: beginner
293: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
294: @*/
295: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
296: {
298: PetscInt m,n,m0,nmat;
299: Mat mat[2];
308: /* Check for square matrices */
309: MatGetSize(A,&m,&n);
310: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
311: if (B) {
312: MatGetSize(B,&m0,&n);
313: if (m0!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
314: if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
315: }
317: if (eps->setupcalled) { EPSReset(eps); }
318: if (!eps->st) { EPSGetST(eps,&eps->st); }
319: mat[0] = A;
320: if (B) {
321: mat[1] = B;
322: nmat = 2;
323: } else nmat = 1;
324: STSetOperators(eps->st,nmat,mat);
325: return(0);
326: }
330: /*@C
331: EPSGetOperators - Gets the matrices associated with the eigensystem.
333: Collective on EPS and Mat
335: Input Parameter:
336: . eps - the EPS context
338: Output Parameters:
339: + A - the matrix associated with the eigensystem
340: - B - the second matrix in the case of generalized eigenproblems
342: Level: intermediate
344: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
345: @*/
346: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
347: {
349: ST st;
350: PetscInt k;
354: EPSGetST(eps,&st);
355: if (A) { STGetOperators(st,0,A); }
356: if (B) {
357: STGetNumMatrices(st,&k);
358: if (k==1) B = NULL;
359: else {
360: STGetOperators(st,1,B);
361: }
362: }
363: return(0);
364: }
368: /*@
369: EPSSetDeflationSpace - Specify a basis of vectors that constitute
370: the deflation space.
372: Collective on EPS and Vec
374: Input Parameter:
375: + eps - the eigenproblem solver context
376: . n - number of vectors
377: - v - set of basis vectors of the deflation space
379: Notes:
380: When a deflation space is given, the eigensolver seeks the eigensolution
381: in the restriction of the problem to the orthogonal complement of this
382: space. This can be used for instance in the case that an invariant
383: subspace is known beforehand (such as the nullspace of the matrix).
385: Basis vectors set by a previous call to EPSSetDeflationSpace() are
386: replaced.
388: The vectors do not need to be mutually orthonormal, since they are explicitly
389: orthonormalized internally.
391: These vectors persist from one EPSSolve() call to the other, use
392: EPSRemoveDeflationSpace() to eliminate them.
394: Level: intermediate
396: .seealso: EPSRemoveDeflationSpace()
397: @*/
398: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)
399: {
401: PetscInt i;
406: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
408: /* free previous vectors */
409: EPSRemoveDeflationSpace(eps);
411: /* get references of passed vectors */
412: if (n>0) {
413: PetscMalloc(n*sizeof(Vec),&eps->defl);
414: PetscLogObjectMemory(eps,n*sizeof(Vec));
415: for (i=0;i<n;i++) {
416: PetscObjectReference((PetscObject)v[i]);
417: eps->defl[i] = v[i];
418: }
419: eps->setupcalled = 0;
420: eps->ds_ortho = PETSC_FALSE;
421: }
423: eps->nds = n;
424: return(0);
425: }
429: /*@
430: EPSRemoveDeflationSpace - Removes the deflation space.
432: Collective on EPS
434: Input Parameter:
435: . eps - the eigenproblem solver context
437: Level: intermediate
439: .seealso: EPSSetDeflationSpace()
440: @*/
441: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
442: {
447: VecDestroyVecs(eps->nds,&eps->defl);
448: eps->nds = 0;
449: eps->setupcalled = 0;
450: eps->ds_ortho = PETSC_FALSE;
451: return(0);
452: }
456: /*@
457: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
458: space, that is, the subspace from which the solver starts to iterate.
460: Collective on EPS and Vec
462: Input Parameter:
463: + eps - the eigenproblem solver context
464: . n - number of vectors
465: - is - set of basis vectors of the initial space
467: Notes:
468: Some solvers start to iterate on a single vector (initial vector). In that case,
469: the other vectors are ignored.
471: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
472: EPSSolve() call to the other, so the initial space should be set every time.
474: The vectors do not need to be mutually orthonormal, since they are explicitly
475: orthonormalized internally.
477: Common usage of this function is when the user can provide a rough approximation
478: of the wanted eigenspace. Then, convergence may be faster.
480: Level: intermediate
482: .seealso: EPSSetInitialSpaceLeft(), EPSSetDeflationSpace()
483: @*/
484: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
485: {
491: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
492: SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS);
493: if (n>0) eps->setupcalled = 0;
494: return(0);
495: }
499: /*@
500: EPSSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
501: left space, that is, the subspace from which the solver starts to iterate for
502: building the left subspace (in methods that work with two subspaces).
504: Collective on EPS and Vec
506: Input Parameter:
507: + eps - the eigenproblem solver context
508: . n - number of vectors
509: - is - set of basis vectors of the initial left space
511: Notes:
512: Some solvers start to iterate on a single vector (initial left vector). In that case,
513: the other vectors are ignored.
515: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
516: EPSSolve() call to the other, so the initial left space should be set every time.
518: The vectors do not need to be mutually orthonormal, since they are explicitly
519: orthonormalized internally.
521: Common usage of this function is when the user can provide a rough approximation
522: of the wanted left eigenspace. Then, convergence may be faster.
524: Level: intermediate
526: .seealso: EPSSetInitialSpace(), EPSSetDeflationSpace()
527: @*/
528: PetscErrorCode EPSSetInitialSpaceLeft(EPS eps,PetscInt n,Vec *is)
529: {
535: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
536: SlepcBasisReference_Private(n,is,&eps->ninil,&eps->ISL);
537: if (n>0) eps->setupcalled = 0;
538: return(0);
539: }