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-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/qepimpl.h> /*I "slepcqep.h" I*/
25: #include <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 (!((PetscObject)qep->rand)->type_name) {
68: PetscRandomSetFromOptions(qep->rand);
69: }
71: /* Check matrices */
72: if (!qep->M || !qep->C || !qep->K)
73: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSetOperators must be called first");
74:
75: /* Set problem dimensions */
76: MatGetSize(qep->M,&qep->n,PETSC_NULL);
77: MatGetLocalSize(qep->M,&qep->nloc,PETSC_NULL);
78: VecDestroy(&qep->t);
79: SlepcMatGetVecsTemplate(qep->M,&qep->t,PETSC_NULL);
81: /* Set default problem type */
82: if (!qep->problem_type) {
83: QEPSetProblemType(qep,QEP_GENERAL);
84: }
86: /* Compute scaling factor if not set by user */
87: if (qep->sfactor==0.0) {
88: MatHasOperation(qep->K,MATOP_NORM,&khas);
89: MatHasOperation(qep->M,MATOP_NORM,&mhas);
90: if (khas && mhas) {
91: MatNorm(qep->K,NORM_INFINITY,&knorm);
92: MatNorm(qep->M,NORM_INFINITY,&mnorm);
93: qep->sfactor = PetscSqrtReal(knorm/mnorm);
94: }
95: else qep->sfactor = 1.0;
96: }
98: /* Call specific solver setup */
99: (*qep->ops->setup)(qep);
101: /* set tolerance if not yet set */
102: if (qep->tol==PETSC_DEFAULT) qep->tol = SLEPC_DEFAULT_TOL;
104: if (qep->ncv > 2*qep->n)
105: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"ncv must be twice the problem size at most");
106: if (qep->nev > qep->ncv)
107: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
109: /* process initial vectors */
110: if (qep->nini<0) {
111: qep->nini = -qep->nini;
112: if (qep->nini>qep->ncv) SETERRQ(((PetscObject)qep)->comm,1,"The number of initial vectors is larger than ncv");
113: k = 0;
114: for (i=0;i<qep->nini;i++) {
115: VecCopy(qep->IS[i],qep->V[k]);
116: VecDestroy(&qep->IS[i]);
117: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->V,qep->V[k],PETSC_NULL,&norm,&lindep);
118: if (norm==0.0 || lindep) PetscInfo(qep,"Linearly dependent initial vector found, removing...\n");
119: else {
120: VecScale(qep->V[k],1.0/norm);
121: k++;
122: }
123: }
124: qep->nini = k;
125: PetscFree(qep->IS);
126: }
127: if (qep->ninil<0) {
128: if (!qep->leftvecs) PetscInfo(qep,"Ignoring initial left vectors\n");
129: else {
130: qep->ninil = -qep->ninil;
131: if (qep->ninil>qep->ncv) SETERRQ(((PetscObject)qep)->comm,1,"The number of initial left vectors is larger than ncv");
132: k = 0;
133: for (i=0;i<qep->ninil;i++) {
134: VecCopy(qep->ISL[i],qep->W[k]);
135: VecDestroy(&qep->ISL[i]);
136: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->W,qep->W[k],PETSC_NULL,&norm,&lindep);
137: if (norm==0.0 || lindep) PetscInfo(qep,"Linearly dependent initial left vector found, removing...\n");
138: else {
139: VecScale(qep->W[k],1.0/norm);
140: k++;
141: }
142: }
143: qep->ninil = k;
144: PetscFree(qep->ISL);
145: }
146: }
147: PetscLogEventEnd(QEP_SetUp,qep,0,0,0);
148: qep->setupcalled = 1;
149: return(0);
150: }
154: /*@
155: QEPSetOperators - Sets the matrices associated with the quadratic eigenvalue problem.
157: Collective on QEP and Mat
159: Input Parameters:
160: + qep - the eigenproblem solver context
161: . M - the first coefficient matrix
162: . C - the second coefficient matrix
163: - K - the third coefficient matrix
165: Notes:
166: The quadratic eigenproblem is defined as (l^2*M + l*C + K)*x = 0, where l is
167: the eigenvalue and x is the eigenvector.
169: Level: beginner
171: .seealso: QEPSolve(), QEPGetOperators()
172: @*/
173: PetscErrorCode QEPSetOperators(QEP qep,Mat M,Mat C,Mat K)
174: {
176: PetscInt m,n,m0;
187: /* Check for square matrices */
188: MatGetSize(M,&m,&n);
189: if (m!=n) { SETERRQ(((PetscObject)qep)->comm,1,"M is a non-square matrix"); }
190: m0=m;
191: MatGetSize(C,&m,&n);
192: if (m!=n) { SETERRQ(((PetscObject)qep)->comm,1,"C is a non-square matrix"); }
193: if (m!=m0) { SETERRQ(((PetscObject)qep)->comm,1,"Dimensions of M and C do not match"); }
194: MatGetSize(K,&m,&n);
195: if (m!=n) { SETERRQ(((PetscObject)qep)->comm,1,"K is a non-square matrix"); }
196: if (m!=m0) { SETERRQ(((PetscObject)qep)->comm,1,"Dimensions of M and K do not match"); }
198: /* Store a copy of the matrices */
199: if (qep->setupcalled) { QEPReset(qep); }
200: PetscObjectReference((PetscObject)M);
201: MatDestroy(&qep->M);
202: qep->M = M;
203: PetscObjectReference((PetscObject)C);
204: MatDestroy(&qep->C);
205: qep->C = C;
206: PetscObjectReference((PetscObject)K);
207: MatDestroy(&qep->K);
208: qep->K = K;
209: return(0);
210: }
214: /*@
215: QEPGetOperators - Gets the matrices associated with the quadratic eigensystem.
217: Collective on QEP and Mat
219: Input Parameter:
220: . qep - the QEP context
222: Output Parameters:
223: + M - the first coefficient matrix
224: . C - the second coefficient matrix
225: - K - the third coefficient matrix
227: Level: intermediate
229: .seealso: QEPSolve(), QEPSetOperators()
230: @*/
231: PetscErrorCode QEPGetOperators(QEP qep,Mat *M,Mat *C,Mat *K)
232: {
238: return(0);
239: }
243: /*@
244: QEPSetInitialSpace - Specify a basis of vectors that constitute the initial
245: space, that is, the subspace from which the solver starts to iterate.
247: Collective on QEP and Vec
249: Input Parameter:
250: + qep - the quadratic eigensolver context
251: . n - number of vectors
252: - is - set of basis vectors of the initial space
254: Notes:
255: Some solvers start to iterate on a single vector (initial vector). In that case,
256: the other vectors are ignored.
258: These vectors do not persist from one QEPSolve() call to the other, so the
259: initial space should be set every time.
261: The vectors do not need to be mutually orthonormal, since they are explicitly
262: orthonormalized internally.
264: Common usage of this function is when the user can provide a rough approximation
265: of the wanted eigenspace. Then, convergence may be faster.
267: Level: intermediate
269: .seealso: QEPSetInitialSpaceLeft()
270: @*/
271: PetscErrorCode QEPSetInitialSpace(QEP qep,PetscInt n,Vec *is)
272: {
274: PetscInt i;
275:
279: if (n<0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
281: /* free previous non-processed vectors */
282: if (qep->nini<0) {
283: for (i=0;i<-qep->nini;i++) {
284: VecDestroy(&qep->IS[i]);
285: }
286: PetscFree(qep->IS);
287: }
289: /* get references of passed vectors */
290: PetscMalloc(n*sizeof(Vec),&qep->IS);
291: for (i=0;i<n;i++) {
292: PetscObjectReference((PetscObject)is[i]);
293: qep->IS[i] = is[i];
294: }
296: qep->nini = -n;
297: qep->setupcalled = 0;
298: return(0);
299: }
303: /*@
304: QEPSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
305: left space, that is, the subspace from which the solver starts to iterate for
306: building the left subspace (in methods that work with two subspaces).
308: Collective on QEP and Vec
310: Input Parameter:
311: + qep - the quadratic eigensolver context
312: . n - number of vectors
313: - is - set of basis vectors of the initial left space
315: Notes:
316: Some solvers start to iterate on a single vector (initial left vector). In that case,
317: the other vectors are ignored.
319: These vectors do not persist from one QEPSolve() call to the other, so the
320: initial left space should be set every time.
322: The vectors do not need to be mutually orthonormal, since they are explicitly
323: orthonormalized internally.
325: Common usage of this function is when the user can provide a rough approximation
326: of the wanted left eigenspace. Then, convergence may be faster.
328: Level: intermediate
330: .seealso: QEPSetInitialSpace()
331: @*/
332: PetscErrorCode QEPSetInitialSpaceLeft(QEP qep,PetscInt n,Vec *is)
333: {
335: PetscInt i;
336:
340: if (n<0) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
342: /* free previous non-processed vectors */
343: if (qep->ninil<0) {
344: for (i=0;i<-qep->ninil;i++) {
345: VecDestroy(&qep->ISL[i]);
346: }
347: PetscFree(qep->ISL);
348: }
350: /* get references of passed vectors */
351: PetscMalloc(n*sizeof(Vec),&qep->ISL);
352: for (i=0;i<n;i++) {
353: PetscObjectReference((PetscObject)is[i]);
354: qep->ISL[i] = is[i];
355: }
357: qep->ninil = -n;
358: qep->setupcalled = 0;
359: return(0);
360: }
364: /*
365: QEPAllocateSolution - Allocate memory storage for common variables such
366: as eigenvalues and eigenvectors. All vectors in V (and W) share a
367: contiguous chunk of memory.
368: */
369: PetscErrorCode QEPAllocateSolution(QEP qep)
370: {
372:
374: if (qep->allocated_ncv != qep->ncv) {
375: QEPFreeSolution(qep);
376: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigr);
377: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigi);
378: PetscMalloc(qep->ncv*sizeof(PetscReal),&qep->errest);
379: PetscMalloc(qep->ncv*sizeof(PetscInt),&qep->perm);
380: VecDuplicateVecs(qep->t,qep->ncv,&qep->V);
381: qep->allocated_ncv = qep->ncv;
382: }
383: return(0);
384: }
388: /*
389: QEPFreeSolution - Free memory storage. This routine is related to
390: QEPAllocateSolution().
391: */
392: PetscErrorCode QEPFreeSolution(QEP qep)
393: {
395:
397: if (qep->allocated_ncv > 0) {
398: PetscFree(qep->eigr);
399: PetscFree(qep->eigi);
400: PetscFree(qep->errest);
401: PetscFree(qep->perm);
402: VecDestroyVecs(qep->allocated_ncv,&qep->V);
403: qep->allocated_ncv = 0;
404: }
405: return(0);
406: }