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-2010, Universidad 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
28: /*@
29: QEPSetUp - Sets up all the internal data structures necessary for the
30: execution of the QEP solver.
32: Collective on QEP
34: Input Parameter:
35: . qep - solver context
37: Notes:
38: This function need not be called explicitly in most cases, since QEPSolve()
39: calls it. It can be useful when one wants to measure the set-up time
40: separately from the solve time.
42: Level: advanced
44: .seealso: QEPCreate(), QEPSolve(), QEPDestroy()
45: @*/
46: PetscErrorCode QEPSetUp(QEP qep)
47: {
49: PetscInt i,k;
50: PetscScalar *pV;
51: PetscTruth khas,mhas,lindep;
52: PetscReal knorm,mnorm,norm;
53:
57: if (qep->setupcalled) return(0);
59: PetscLogEventBegin(QEP_SetUp,qep,0,0,0);
61: /* Set default solver type */
62: if (!((PetscObject)qep)->type_name) {
63: QEPSetType(qep,QEPLINEAR);
64: }
66: /* Check matrices */
67: if (!qep->M || !qep->C || !qep->K)
68: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "QEPSetOperators must be called first");
69:
70: /* Set problem dimensions */
71: MatGetSize(qep->M,&qep->n,PETSC_NULL);
72: MatGetLocalSize(qep->M,&qep->nloc,PETSC_NULL);
74: /* Set default problem type */
75: if (!qep->problem_type) {
76: QEPSetProblemType(qep,QEP_GENERAL);
77: }
79: /* Compute scaling factor if not set by user */
80: if (qep->sfactor==0.0) {
81: MatHasOperation(qep->K,MATOP_NORM,&khas);
82: MatHasOperation(qep->M,MATOP_NORM,&mhas);
83: if (khas && mhas) {
84: MatNorm(qep->K,NORM_INFINITY,&knorm);
85: MatNorm(qep->M,NORM_INFINITY,&mnorm);
86: qep->sfactor = sqrt(knorm/mnorm);
87: }
88: else qep->sfactor = 1.0;
89: }
91: /* initialize the random number generator */
92: PetscRandomCreate(((PetscObject)qep)->comm,&qep->rand);
93: PetscRandomSetFromOptions(qep->rand);
95: /* Call specific solver setup */
96: (*qep->ops->setup)(qep);
98: if (qep->ncv > 2*qep->n)
99: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"ncv must be twice the problem size at most");
100: if (qep->nev > qep->ncv)
101: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
103: /* Free memory for previous solution */
104: if (qep->eigr) {
105: PetscFree(qep->eigr);
106: PetscFree(qep->eigi);
107: PetscFree(qep->perm);
108: PetscFree(qep->errest);
109: VecGetArray(qep->V[0],&pV);
110: for (i=0;i<qep->ncv;i++) {
111: VecDestroy(qep->V[i]);
112: }
113: PetscFree(pV);
114: PetscFree(qep->V);
115: }
117: /* Allocate memory for next solution */
118: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigr);
119: PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigi);
120: PetscMalloc(qep->ncv*sizeof(PetscInt),&qep->perm);
121: PetscMalloc(qep->ncv*sizeof(PetscReal),&qep->errest);
122: PetscMalloc(qep->ncv*sizeof(Vec),&qep->V);
123: PetscMalloc(qep->ncv*qep->nloc*sizeof(PetscScalar),&pV);
124: for (i=0;i<qep->ncv;i++) {
125: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,PETSC_DECIDE,pV+i*qep->nloc,&qep->V[i]);
126: }
128: /* process initial vectors */
129: if (qep->nini<0) {
130: qep->nini = -qep->nini;
131: if (qep->nini>qep->ncv) SETERRQ(1,"The number of initial vectors is larger than ncv")
132: k = 0;
133: for (i=0;i<qep->nini;i++) {
134: VecCopy(qep->IS[i],qep->V[k]);
135: VecDestroy(qep->IS[i]);
136: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->V,qep->V[k],PETSC_NULL,&norm,&lindep);
137: if (norm==0.0 || lindep) PetscInfo(qep,"Linearly dependent initial vector found, removing...\n");
138: else {
139: VecScale(qep->V[k],1.0/norm);
140: k++;
141: }
142: }
143: qep->nini = k;
144: PetscFree(qep->IS);
145: }
146: if (qep->ninil<0) {
147: if (!qep->leftvecs) PetscInfo(qep,"Ignoring initial left vectors\n");
148: else {
149: qep->ninil = -qep->ninil;
150: if (qep->ninil>qep->ncv) SETERRQ(1,"The number of initial left vectors is larger than ncv")
151: k = 0;
152: for (i=0;i<qep->ninil;i++) {
153: VecCopy(qep->ISL[i],qep->W[k]);
154: VecDestroy(qep->ISL[i]);
155: IPOrthogonalize(qep->ip,0,PETSC_NULL,k,PETSC_NULL,qep->W,qep->W[k],PETSC_NULL,&norm,&lindep);
156: if (norm==0.0 || lindep) PetscInfo(qep,"Linearly dependent initial left vector found, removing...\n");
157: else {
158: VecScale(qep->W[k],1.0/norm);
159: k++;
160: }
161: }
162: qep->ninil = k;
163: PetscFree(qep->ISL);
164: }
165: }
167: PetscLogEventEnd(QEP_SetUp,qep,0,0,0);
168: qep->setupcalled = 1;
169: return(0);
170: }
174: /*@
175: QEPSetOperators - Sets the matrices associated with the quadratic eigenvalue problem.
177: Collective on QEP and Mat
179: Input Parameters:
180: + qep - the eigenproblem solver context
181: . M - the first coefficient matrix
182: . C - the second coefficient matrix
183: - K - the third coefficient matrix
185: Notes:
186: The quadratic eigenproblem is defined as (l^2*M + l*C + K)*x = 0, where l is
187: the eigenvalue and x is the eigenvector.
189: Level: beginner
191: .seealso: QEPSolve(), QEPGetOperators()
192: @*/
193: PetscErrorCode QEPSetOperators(QEP qep,Mat M,Mat C,Mat K)
194: {
196: PetscInt m,n,m0;
207: /* Check for square matrices */
208: MatGetSize(M,&m,&n);
209: if (m!=n) { SETERRQ(1,"M is a non-square matrix"); }
210: m0=m;
211: MatGetSize(C,&m,&n);
212: if (m!=n) { SETERRQ(1,"C is a non-square matrix"); }
213: if (m!=m0) { SETERRQ(1,"Dimensions of M and C do not match"); }
214: MatGetSize(K,&m,&n);
215: if (m!=n) { SETERRQ(1,"K is a non-square matrix"); }
216: if (m!=m0) { SETERRQ(1,"Dimensions of M and K do not match"); }
218: /* Store a copy of the matrices */
219: PetscObjectReference((PetscObject)M);
220: if (qep->M) {
221: MatDestroy(qep->M);
222: }
223: qep->M = M;
224: PetscObjectReference((PetscObject)C);
225: if (qep->C) {
226: MatDestroy(qep->C);
227: }
228: qep->C = C;
229: PetscObjectReference((PetscObject)K);
230: if (qep->K) {
231: MatDestroy(qep->K);
232: }
233: qep->K = K;
235: qep->setupcalled = 0;
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:
305: if (n<0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
307: /* free previous non-processed vectors */
308: if (qep->nini<0) {
309: for (i=0;i<-qep->nini;i++) {
310: VecDestroy(qep->IS[i]);
311: }
312: PetscFree(qep->IS);
313: }
315: /* get references of passed vectors */
316: PetscMalloc(n*sizeof(Vec),&qep->IS);
317: for (i=0;i<n;i++) {
318: PetscObjectReference((PetscObject)is[i]);
319: qep->IS[i] = is[i];
320: }
322: qep->nini = -n;
323: qep->setupcalled = 0;
324: return(0);
325: }
329: /*@
330: QEPSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
331: left space, that is, the subspace from which the solver starts to iterate for
332: building the left subspace (in methods that work with two subspaces).
334: Collective on QEP and Vec
336: Input Parameter:
337: + qep - the quadratic eigensolver context
338: . n - number of vectors
339: - is - set of basis vectors of the initial left space
341: Notes:
342: Some solvers start to iterate on a single vector (initial left vector). In that case,
343: the other vectors are ignored.
345: These vectors do not persist from one QEPSolve() call to the other, so the
346: initial left space should be set every time.
348: The vectors do not need to be mutually orthonormal, since they are explicitly
349: orthonormalized internally.
351: Common usage of this function is when the user can provide a rough approximation
352: of the wanted left eigenspace. Then, convergence may be faster.
354: Level: intermediate
356: .seealso: QEPSetInitialSpace()
357: @*/
358: PetscErrorCode QEPSetInitialSpaceLeft(QEP qep,PetscInt n,Vec *is)
359: {
361: PetscInt i;
362:
365: if (n<0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
367: /* free previous non-processed vectors */
368: if (qep->ninil<0) {
369: for (i=0;i<-qep->ninil;i++) {
370: VecDestroy(qep->ISL[i]);
371: }
372: PetscFree(qep->ISL);
373: }
375: /* get references of passed vectors */
376: PetscMalloc(n*sizeof(Vec),&qep->ISL);
377: for (i=0;i<n;i++) {
378: PetscObjectReference((PetscObject)is[i]);
379: qep->ISL[i] = is[i];
380: }
382: qep->ninil = -n;
383: qep->setupcalled = 0;
384: return(0);
385: }