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-2012, 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 <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;
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: PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSGD,EPSJD,EPSRQCG,"");
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 (!eps->ds) { EPSGetDS(eps,&eps->ds); }
79: DSReset(eps->ds);
80: if (!((PetscObject)eps->rand)->type_name) {
81: PetscRandomSetFromOptions(eps->rand);
82: }
83:
84: /* Set problem dimensions */
85: STGetOperators(eps->OP,&A,&B);
86: if (!A) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
87: MatGetSize(A,&eps->n,PETSC_NULL);
88: MatGetLocalSize(A,&eps->nloc,PETSC_NULL);
89: VecDestroy(&eps->t);
90: SlepcMatGetVecsTemplate(A,&eps->t,PETSC_NULL);
92: /* Set default problem type */
93: if (!eps->problem_type) {
94: if (B==PETSC_NULL) {
95: EPSSetProblemType(eps,EPS_NHEP);
96: } else {
97: EPSSetProblemType(eps,EPS_GNHEP);
98: }
99: } else if (!B && eps->isgeneralized) {
100: PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
101: eps->isgeneralized = PETSC_FALSE;
102: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
103: } else if (B && !eps->isgeneralized) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");
104: #if defined(PETSC_USE_COMPLEX)
105: STGetShift(eps->OP,&sigma);
106: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
107: #endif
108: if (eps->ishermitian && eps->leftvecs) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Requesting left eigenvectors not allowed in Hermitian problems");
109:
110: if (eps->ispositive || (eps->isgeneralized && eps->ishermitian)) {
111: STGetBilinearForm(eps->OP,&B);
112: IPSetMatrix(eps->ip,B);
113: MatDestroy(&B);
114: if (!eps->ispositive) { IPSetType(eps->ip,IPINDEFINITE); }
115: } else {
116: IPSetMatrix(eps->ip,PETSC_NULL);
117: }
118:
119: if (eps->nev > eps->n) eps->nev = eps->n;
120: if (eps->ncv > eps->n) eps->ncv = eps->n;
122: /* initialization of matrix norms */
123: if (eps->nrma == PETSC_DETERMINE) {
124: MatHasOperation(A,MATOP_NORM,&flg);
125: if (flg) { MatNorm(A,NORM_INFINITY,&eps->nrma); }
126: else eps->nrma = 1.0;
127: }
128: if (eps->nrmb == PETSC_DETERMINE) {
129: MatHasOperation(B,MATOP_NORM,&flg);
130: if (flg) { MatNorm(B,NORM_INFINITY,&eps->nrmb); }
131: else eps->nrmb = 1.0;
132: }
134: if (!eps->balance) eps->balance = EPS_BALANCE_NONE;
136: /* call specific solver setup */
137: (*eps->ops->setup)(eps);
139: /* check extraction */
140: PetscObjectTypeCompareAny((PetscObject)eps->OP,&flg,STPRECOND,STSHIFT,"");
141: if (!flg && eps->extraction && eps->extraction!=EPS_RITZ) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
143: /* set tolerance if not yet set */
144: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
146: /* set eigenvalue comparison */
147: switch (eps->which) {
148: case EPS_LARGEST_MAGNITUDE:
149: eps->which_func = SlepcCompareLargestMagnitude;
150: eps->which_ctx = PETSC_NULL;
151: break;
152: case EPS_SMALLEST_MAGNITUDE:
153: eps->which_func = SlepcCompareSmallestMagnitude;
154: eps->which_ctx = PETSC_NULL;
155: break;
156: case EPS_LARGEST_REAL:
157: eps->which_func = SlepcCompareLargestReal;
158: eps->which_ctx = PETSC_NULL;
159: break;
160: case EPS_SMALLEST_REAL:
161: eps->which_func = SlepcCompareSmallestReal;
162: eps->which_ctx = PETSC_NULL;
163: break;
164: case EPS_LARGEST_IMAGINARY:
165: eps->which_func = SlepcCompareLargestImaginary;
166: eps->which_ctx = PETSC_NULL;
167: break;
168: case EPS_SMALLEST_IMAGINARY:
169: eps->which_func = SlepcCompareSmallestImaginary;
170: eps->which_ctx = PETSC_NULL;
171: break;
172: case EPS_TARGET_MAGNITUDE:
173: eps->which_func = SlepcCompareTargetMagnitude;
174: eps->which_ctx = &eps->target;
175: break;
176: case EPS_TARGET_REAL:
177: eps->which_func = SlepcCompareTargetReal;
178: eps->which_ctx = &eps->target;
179: break;
180: case EPS_TARGET_IMAGINARY:
181: eps->which_func = SlepcCompareTargetImaginary;
182: eps->which_ctx = &eps->target;
183: break;
184: case EPS_ALL:
185: eps->which_func = SlepcCompareSmallestReal;
186: eps->which_ctx = PETSC_NULL;
187: break;
188: case EPS_WHICH_USER:
189: break;
190: }
192: /* Build balancing matrix if required */
193: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
194: if (!eps->D) {
195: VecDuplicate(eps->V[0],&eps->D);
196: } else {
197: VecSet(eps->D,1.0);
198: }
199: EPSBuildBalance_Krylov(eps);
200: STSetBalanceMatrix(eps->OP,eps->D);
201: }
203: /* Setup ST */
204: STSetUp(eps->OP);
205:
206: PetscObjectTypeCompare((PetscObject)eps->OP,STCAYLEY,&flg);
207: if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
209: PetscObjectTypeCompare((PetscObject)eps->OP,STFOLD,&flg);
210: if (flg && !eps->ishermitian) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_SUP,"Fold spectral transformation requires a Hermitian problem");
212: if (eps->nds>0) {
213: if (!eps->ds_ortho) {
214: /* allocate memory and copy deflation basis vectors into defl */
215: VecDuplicateVecs(eps->t,eps->nds,&newDS);
216: for (i=0;i<eps->nds;i++) {
217: VecCopy(eps->defl[i],newDS[i]);
218: VecDestroy(&eps->defl[i]);
219: }
220: PetscFree(eps->defl);
221: eps->defl = newDS;
222: /* orthonormalize vectors in defl */
223: k = 0;
224: for (i=0;i<eps->nds;i++) {
225: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->defl,eps->defl[k],PETSC_NULL,&norm,&lindep);
226: if (norm==0.0 || lindep) {
227: PetscInfo(eps,"Linearly dependent deflation vector found, removing...\n");
228: } else {
229: VecScale(eps->defl[k],1.0/norm);
230: k++;
231: }
232: }
233: for (i=k;i<eps->nds;i++) { VecDestroy(&eps->defl[i]); }
234: eps->nds = k;
235: eps->ds_ortho = PETSC_TRUE;
236: }
237: }
238: STCheckNullSpace(eps->OP,eps->nds,eps->defl);
240: /* process initial vectors */
241: if (eps->nini<0) {
242: eps->nini = -eps->nini;
243: if (eps->nini>eps->ncv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"The number of initial vectors is larger than ncv");
244: k = 0;
245: for (i=0;i<eps->nini;i++) {
246: VecCopy(eps->IS[i],eps->V[k]);
247: VecDestroy(&eps->IS[i]);
248: IPOrthogonalize(eps->ip,eps->nds,eps->defl,k,PETSC_NULL,eps->V,eps->V[k],PETSC_NULL,&norm,&lindep);
249: if (norm==0.0 || lindep) {
250: PetscInfo(eps,"Linearly dependent initial vector found, removing...\n");
251: } else {
252: VecScale(eps->V[k],1.0/norm);
253: k++;
254: }
255: }
256: eps->nini = k;
257: PetscFree(eps->IS);
258: }
259: if (eps->ninil<0) {
260: if (!eps->leftvecs) {
261: PetscInfo(eps,"Ignoring initial left vectors\n");
262: } else {
263: eps->ninil = -eps->ninil;
264: if (eps->ninil>eps->ncv) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"The number of initial left vectors is larger than ncv");
265: k = 0;
266: for (i=0;i<eps->ninil;i++) {
267: VecCopy(eps->ISL[i],eps->W[k]);
268: VecDestroy(&eps->ISL[i]);
269: IPOrthogonalize(eps->ip,0,PETSC_NULL,k,PETSC_NULL,eps->W,eps->W[k],PETSC_NULL,&norm,&lindep);
270: if (norm==0.0 || lindep) {
271: PetscInfo(eps,"Linearly dependent initial left vector found, removing...\n");
272: } else {
273: VecScale(eps->W[k],1.0/norm);
274: k++;
275: }
276: }
277: eps->ninil = k;
278: PetscFree(eps->ISL);
279: }
280: }
282: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
283: eps->setupcalled = 1;
284: return(0);
285: }
289: /*@
290: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
292: Collective on EPS and Mat
294: Input Parameters:
295: + eps - the eigenproblem solver context
296: . A - the matrix associated with the eigensystem
297: - B - the second matrix in the case of generalized eigenproblems
299: Notes:
300: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
302: It must be called after EPSSetUp(). If it is called again after EPSSetUp() then
303: the EPS object is reset.
305: Level: beginner
307: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
308: @*/
309: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
310: {
312: PetscInt m,n,m0;
321: /* Check for square matrices */
322: MatGetSize(A,&m,&n);
323: if (m!=n) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
324: if (B) {
325: MatGetSize(B,&m0,&n);
326: if (m0!=n) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
327: if (m!=m0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
328: }
330: if (eps->setupcalled) { EPSReset(eps); }
331: if (!eps->OP) { EPSGetST(eps,&eps->OP); }
332: STSetOperators(eps->OP,A,B);
333: return(0);
334: }
338: /*@
339: EPSGetOperators - Gets the matrices associated with the eigensystem.
341: Collective on EPS and Mat
343: Input Parameter:
344: . eps - the EPS context
346: Output Parameters:
347: + A - the matrix associated with the eigensystem
348: - B - the second matrix in the case of generalized eigenproblems
350: Level: intermediate
352: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
353: @*/
354: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
355: {
357: ST st;
363: EPSGetST(eps,&st);
364: STGetOperators(st,A,B);
365: return(0);
366: }
370: /*@
371: EPSSetDeflationSpace - Specify a basis of vectors that constitute
372: the deflation space.
374: Collective on EPS and Vec
376: Input Parameter:
377: + eps - the eigenproblem solver context
378: . n - number of vectors
379: - v - set of basis vectors of the deflation space
381: Notes:
382: When a deflation space is given, the eigensolver seeks the eigensolution
383: in the restriction of the problem to the orthogonal complement of this
384: space. This can be used for instance in the case that an invariant
385: subspace is known beforehand (such as the nullspace of the matrix).
387: Basis vectors set by a previous call to EPSSetDeflationSpace() are
388: replaced.
390: The vectors do not need to be mutually orthonormal, since they are explicitly
391: orthonormalized internally.
393: These vectors persist from one EPSSolve() call to the other, use
394: EPSRemoveDeflationSpace() to eliminate them.
396: Level: intermediate
398: .seealso: EPSRemoveDeflationSpace()
399: @*/
400: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)
401: {
403: PetscInt i;
404:
408: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
410: /* free previous vectors */
411: EPSRemoveDeflationSpace(eps);
413: /* get references of passed vectors */
414: if (n>0) {
415: PetscMalloc(n*sizeof(Vec),&eps->defl);
416: for (i=0;i<n;i++) {
417: PetscObjectReference((PetscObject)v[i]);
418: eps->defl[i] = v[i];
419: }
420: eps->setupcalled = 0;
421: eps->ds_ortho = PETSC_FALSE;
422: }
424: eps->nds = n;
425: return(0);
426: }
430: /*@
431: EPSRemoveDeflationSpace - Removes the deflation space.
433: Collective on EPS
435: Input Parameter:
436: . eps - the eigenproblem solver context
438: Level: intermediate
440: .seealso: EPSSetDeflationSpace()
441: @*/
442: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
443: {
445:
448: VecDestroyVecs(eps->nds,&eps->defl);
449: eps->nds = 0;
450: eps->setupcalled = 0;
451: eps->ds_ortho = PETSC_FALSE;
452: return(0);
453: }
457: /*@
458: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
459: space, that is, the subspace from which the solver starts to iterate.
461: Collective on EPS and Vec
463: Input Parameter:
464: + eps - the eigenproblem solver context
465: . n - number of vectors
466: - is - set of basis vectors of the initial space
468: Notes:
469: Some solvers start to iterate on a single vector (initial vector). In that case,
470: the other vectors are ignored.
472: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
473: EPSSolve() call to the other, so the initial space should be set every time.
475: The vectors do not need to be mutually orthonormal, since they are explicitly
476: orthonormalized internally.
478: Common usage of this function is when the user can provide a rough approximation
479: of the wanted eigenspace. Then, convergence may be faster.
481: Level: intermediate
483: .seealso: EPSSetInitialSpaceLeft(), EPSSetDeflationSpace()
484: @*/
485: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
486: {
488: PetscInt i;
489:
493: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
495: /* free previous non-processed vectors */
496: if (eps->nini<0) {
497: for (i=0;i<-eps->nini;i++) {
498: VecDestroy(&eps->IS[i]);
499: }
500: PetscFree(eps->IS);
501: }
503: /* get references of passed vectors */
504: if (n>0) {
505: PetscMalloc(n*sizeof(Vec),&eps->IS);
506: for (i=0;i<n;i++) {
507: PetscObjectReference((PetscObject)is[i]);
508: eps->IS[i] = is[i];
509: }
510: eps->setupcalled = 0;
511: }
513: eps->nini = -n;
514: return(0);
515: }
519: /*@
520: EPSSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
521: left space, that is, the subspace from which the solver starts to iterate for
522: building the left subspace (in methods that work with two subspaces).
524: Collective on EPS and Vec
526: Input Parameter:
527: + eps - the eigenproblem solver context
528: . n - number of vectors
529: - is - set of basis vectors of the initial left space
531: Notes:
532: Some solvers start to iterate on a single vector (initial left vector). In that case,
533: the other vectors are ignored.
535: In contrast to EPSSetDeflationSpace(), these vectors do not persist from one
536: EPSSolve() call to the other, so the initial left space should be set every time.
538: The vectors do not need to be mutually orthonormal, since they are explicitly
539: orthonormalized internally.
541: Common usage of this function is when the user can provide a rough approximation
542: of the wanted left eigenspace. Then, convergence may be faster.
544: Level: intermediate
546: .seealso: EPSSetInitialSpace(), EPSSetDeflationSpace()
547: @*/
548: PetscErrorCode EPSSetInitialSpaceLeft(EPS eps,PetscInt n,Vec *is)
549: {
551: PetscInt i;
552:
556: if (n<0) SETERRQ(((PetscObject)eps)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
558: /* free previous non-processed vectors */
559: if (eps->ninil<0) {
560: for (i=0;i<-eps->ninil;i++) {
561: VecDestroy(&eps->ISL[i]);
562: }
563: PetscFree(eps->ISL);
564: }
566: /* get references of passed vectors */
567: if (n>0) {
568: PetscMalloc(n*sizeof(Vec),&eps->ISL);
569: for (i=0;i<n;i++) {
570: PetscObjectReference((PetscObject)is[i]);
571: eps->ISL[i] = is[i];
572: }
573: eps->setupcalled = 0;
574: }
576: eps->ninil = -n;
577: return(0);
578: }