Actual source code: qepsetup.c
1: /*
2: QEP 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/qepimpl.h> /*I "slepcqep.h" I*/
25: #include <slepc-private/ipimpl.h>
29: /*@
30: QEPSetUp - Sets up all the internal data structures necessary for the
31: execution of the QEP solver.
33: Collective on QEP
35: Input Parameter:
36: . qep - solver context
38: Notes:
39: This function need not be called explicitly in most cases, since QEPSolve()
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: QEPCreate(), QEPSolve(), QEPDestroy()
46: @*/
47: PetscErrorCode QEPSetUp(QEP qep)
48: {
50: PetscInt i,k;
51: PetscBool khas,mhas,lindep;
52: PetscReal knorm,mnorm,norm;
53:
56: if (qep->setupcalled) return(0);
57: PetscLogEventBegin(QEP_SetUp,qep,0,0,0);
59: /* Set default solver type (QEPSetFromOptions was not called) */
60: if (!((PetscObject)qep)->type_name) {
61: QEPSetType(qep,QEPLINEAR);
62: }
63: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
64: if (!((PetscObject)qep->ip)->type_name) {
65: IPSetDefaultType_Private(qep->ip);
66: }
67: if (!qep->ds) { QEPGetDS(qep,&qep->ds); }
68: DSReset(qep->ds);
69: if (!((PetscObject)qep->rand)->type_name) {
70: PetscRandomSetFromOptions(qep->rand);
71: }
73: /* Check matrices */
74: if (!qep->M || !qep->C || !qep->K) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSetOperators must be called first");
75:
76: /* Set problem dimensions */
77: MatGetSize(qep->M,&qep->n,PETSC_NULL);
78: MatGetLocalSize(qep->M,&qep->nloc,PETSC_NULL);
79: VecDestroy(&qep->t);
80: SlepcMatGetVecsTemplate(qep->M,&qep->t,PETSC_NULL);
82: /* Set default problem type */
83: if (!qep->problem_type) {
84: QEPSetProblemType(qep,QEP_GENERAL);
85: }
87: /* Compute scaling factor if not set by user */
88: if (qep->sfactor==0.0) {
89: MatHasOperation(qep->K,MATOP_NORM,&khas);
90: MatHasOperation(qep->M,MATOP_NORM,&mhas);
91: if (khas && mhas) {
92: MatNorm(qep->K,NORM_INFINITY,&knorm);
93: MatNorm(qep->M,NORM_INFINITY,&mnorm);
94: qep->sfactor = PetscSqrtReal(knorm/mnorm);
95: }
96: else qep->sfactor = 1.0;
97: }
99: /* Call specific solver setup */
100: (*qep->ops->setup)(qep);
102: /* set tolerance if not yet set */
103: if (qep->tol==PETSC_DEFAULT) qep->tol = SLEPC_DEFAULT_TOL;
105: /* set eigenvalue comparison */
106: switch (qep->which) {
107: case QEP_LARGEST_MAGNITUDE:
108: qep->which_func = SlepcCompareLargestMagnitude;
109: qep->which_ctx = PETSC_NULL;
110: break;
111: case QEP_SMALLEST_MAGNITUDE:
112: qep->which_func = SlepcCompareSmallestMagnitude;
113: qep->which_ctx = PETSC_NULL;
114: break;
115: case QEP_LARGEST_REAL:
116: qep->which_func = SlepcCompareLargestReal;
117: qep->which_ctx = PETSC_NULL;
118: break;
119: case QEP_SMALLEST_REAL:
120: qep->which_func = SlepcCompareSmallestReal;
121: qep->which_ctx = PETSC_NULL;
122: break;
123: case QEP_LARGEST_IMAGINARY:
124: qep->which_func = SlepcCompareLargestImaginary;
125: qep->which_ctx = PETSC_NULL;
126: break;
127: case QEP_SMALLEST_IMAGINARY:
128: qep->which_func = SlepcCompareSmallestImaginary;
129: qep->which_ctx = PETSC_NULL;
130: break;
131: }
133: if (qep->ncv > 2*qep->n) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"ncv must be twice the problem size at most");
134: if (qep->nev > qep->ncv) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
136: /* process initial vectors */
137: if (qep->nini<0) {
138: qep->nini = -qep->nini;
139: if (qep->nini>qep->ncv) SETERRQ(((PetscObject)qep)->comm,1,"The number of initial vectors is larger than ncv");
140: k = 0;
141: for (i=0;i<qep->nini;i++) {
142: VecCopy(qep->IS[i],qep->V[k]);
143: VecDestroy(&qep->IS[i]);
144: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->V,qep->V[k],PETSC_NULL,&norm,&lindep);
145: if (norm==0.0 || lindep) { PetscInfo(qep,"Linearly dependent initial vector found, removing...\n"); }
146: else {
147: VecScale(qep->V[k],1.0/norm);
148: k++;
149: }
150: }
151: qep->nini = k;
152: PetscFree(qep->IS);
153: }
154: if (qep->ninil<0) {
155: if (!qep->leftvecs) { PetscInfo(qep,"Ignoring initial left vectors\n"); }
156: else {
157: qep->ninil = -qep->ninil;
158: if (qep->ninil>qep->ncv) SETERRQ(((PetscObject)qep)->comm,1,"The number of initial left vectors is larger than ncv");
159: k = 0;
160: for (i=0;i<qep->ninil;i++) {
161: VecCopy(qep->ISL[i],qep->W[k]);
162: VecDestroy(&qep->ISL[i]);
163: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->W,qep->W[k],PETSC_NULL,&norm,&lindep);
164: if (norm==0.0 || lindep) { PetscInfo(qep,"Linearly dependent initial left vector found, removing...\n"); }
165: else {
166: VecScale(qep->W[k],1.0/norm);
167: k++;
168: }
169: }
170: qep->ninil = k;
171: PetscFree(qep->ISL);
172: }
173: }
174: PetscLogEventEnd(QEP_SetUp,qep,0,0,0);
175: qep->setupcalled = 1;
176: return(0);
177: }
181: /*@
182: QEPSetOperators - Sets the matrices associated with the quadratic eigenvalue problem.
184: Collective on QEP and Mat
186: Input Parameters:
187: + qep - the eigenproblem solver context
188: . M - the first coefficient matrix
189: . C - the second coefficient matrix
190: - K - the third coefficient matrix
192: Notes:
193: The quadratic eigenproblem is defined as (l^2*M + l*C + K)*x = 0, where l is
194: the eigenvalue and x is the eigenvector.
196: Level: beginner
198: .seealso: QEPSolve(), QEPGetOperators()
199: @*/
200: PetscErrorCode QEPSetOperators(QEP qep,Mat M,Mat C,Mat K)
201: {
203: PetscInt m,n,m0;
214: /* Check for square matrices */
215: MatGetSize(M,&m,&n);
216: if (m!=n) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONG,"M is a non-square matrix");
217: m0=m;
218: MatGetSize(C,&m,&n);
219: if (m!=n) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONG,"C is a non-square matrix");
220: if (m!=m0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of M and C do not match");
221: MatGetSize(K,&m,&n);
222: if (m!=n) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONG,"K is a non-square matrix");
223: if (m!=m0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_INCOMP,"Dimensions of M and K do not match");
225: /* Store a copy of the matrices */
226: if (qep->setupcalled) { QEPReset(qep); }
227: PetscObjectReference((PetscObject)M);
228: MatDestroy(&qep->M);
229: qep->M = M;
230: PetscObjectReference((PetscObject)C);
231: MatDestroy(&qep->C);
232: qep->C = C;
233: PetscObjectReference((PetscObject)K);
234: MatDestroy(&qep->K);
235: qep->K = K;
236: return(0);
237: }
241: /*@
242: QEPGetOperators - Gets the matrices associated with the quadratic eigensystem.
244: Collective on QEP and Mat
246: Input Parameter:
247: . qep - the QEP context
249: Output Parameters:
250: + M - the first coefficient matrix
251: . C - the second coefficient matrix
252: - K - the third coefficient matrix
254: Level: intermediate
256: .seealso: QEPSolve(), QEPSetOperators()
257: @*/
258: PetscErrorCode QEPGetOperators(QEP qep,Mat *M,Mat *C,Mat *K)
259: {
265: return(0);
266: }
270: /*@
271: QEPSetInitialSpace - Specify a basis of vectors that constitute the initial
272: space, that is, the subspace from which the solver starts to iterate.
274: Collective on QEP and Vec
276: Input Parameter:
277: + qep - the quadratic eigensolver context
278: . n - number of vectors
279: - is - set of basis vectors of the initial space
281: Notes:
282: Some solvers start to iterate on a single vector (initial vector). In that case,
283: the other vectors are ignored.
285: These vectors do not persist from one QEPSolve() call to the other, so the
286: initial space should be set every time.
288: The vectors do not need to be mutually orthonormal, since they are explicitly
289: orthonormalized internally.
291: Common usage of this function is when the user can provide a rough approximation
292: of the wanted eigenspace. Then, convergence may be faster.
294: Level: intermediate
296: .seealso: QEPSetInitialSpaceLeft()
297: @*/
298: PetscErrorCode QEPSetInitialSpace(QEP qep,PetscInt n,Vec *is)
299: {
301: PetscInt i;
302:
306: if (n<0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
308: /* free previous non-processed vectors */
309: if (qep->nini<0) {
310: for (i=0;i<-qep->nini;i++) {
311: VecDestroy(&qep->IS[i]);
312: }
313: PetscFree(qep->IS);
314: }
316: /* get references of passed vectors */
317: PetscMalloc(n*sizeof(Vec),&qep->IS);
318: for (i=0;i<n;i++) {
319: PetscObjectReference((PetscObject)is[i]);
320: qep->IS[i] = is[i];
321: }
323: qep->nini = -n;
324: qep->setupcalled = 0;
325: return(0);
326: }
330: /*@
331: QEPSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
332: left space, that is, the subspace from which the solver starts to iterate for
333: building the left subspace (in methods that work with two subspaces).
335: Collective on QEP and Vec
337: Input Parameter:
338: + qep - the quadratic eigensolver context
339: . n - number of vectors
340: - is - set of basis vectors of the initial left space
342: Notes:
343: Some solvers start to iterate on a single vector (initial left vector). In that case,
344: the other vectors are ignored.
346: These vectors do not persist from one QEPSolve() call to the other, so the
347: initial left space should be set every time.
349: The vectors do not need to be mutually orthonormal, since they are explicitly
350: orthonormalized internally.
352: Common usage of this function is when the user can provide a rough approximation
353: of the wanted left eigenspace. Then, convergence may be faster.
355: Level: intermediate
357: .seealso: QEPSetInitialSpace()
358: @*/
359: PetscErrorCode QEPSetInitialSpaceLeft(QEP qep,PetscInt n,Vec *is)
360: {
362: PetscInt i;
363:
367: if (n<0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
369: /* free previous non-processed vectors */
370: if (qep->ninil<0) {
371: for (i=0;i<-qep->ninil;i++) {
372: VecDestroy(&qep->ISL[i]);
373: }
374: PetscFree(qep->ISL);
375: }
377: /* get references of passed vectors */
378: PetscMalloc(n*sizeof(Vec),&qep->ISL);
379: for (i=0;i<n;i++) {
380: PetscObjectReference((PetscObject)is[i]);
381: qep->ISL[i] = is[i];
382: }
384: qep->ninil = -n;
385: qep->setupcalled = 0;
386: return(0);
387: }
391: /*
392: QEPAllocateSolution - Allocate memory storage for common variables such
393: as eigenvalues and eigenvectors. All vectors in V (and W) share a
394: contiguous chunk of memory.
395: */
396: PetscErrorCode QEPAllocateSolution(QEP qep)
397: {
399: PetscInt newc,cnt;
400:
402: if (qep->allocated_ncv != qep->ncv) {
403: newc = PetscMax(0,qep->ncv-qep->allocated_ncv);
404: QEPFreeSolution(qep);
405: cnt = 0;
406: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigr);
407: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigi);
408: cnt += 2*newc*sizeof(PetscScalar);
409: PetscMalloc(qep->ncv*sizeof(PetscReal),&qep->errest);
410: PetscMalloc(qep->ncv*sizeof(PetscInt),&qep->perm);
411: cnt += 2*newc*sizeof(PetscReal);
412: PetscLogObjectMemory(qep,cnt);
413: VecDuplicateVecs(qep->t,qep->ncv,&qep->V);
414: PetscLogObjectParents(qep,qep->ncv,qep->V);
415: qep->allocated_ncv = qep->ncv;
416: }
417: return(0);
418: }
422: /*
423: QEPFreeSolution - Free memory storage. This routine is related to
424: QEPAllocateSolution().
425: */
426: PetscErrorCode QEPFreeSolution(QEP qep)
427: {
429:
431: if (qep->allocated_ncv > 0) {
432: PetscFree(qep->eigr);
433: PetscFree(qep->eigi);
434: PetscFree(qep->errest);
435: PetscFree(qep->perm);
436: VecDestroyVecs(qep->allocated_ncv,&qep->V);
437: qep->allocated_ncv = 0;
438: }
439: return(0);
440: }