Actual source code: linear.c
1: /*
3: Straightforward linearization for quadratic eigenproblems.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include private/qepimpl.h
26: #include private/epsimpl.h
27: #include slepceps.h
28: #include linearp.h
32: PetscErrorCode QEPSetUp_LINEAR(QEP qep)
33: {
34: PetscErrorCode ierr;
35: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
36: PetscInt i=0;
37: EPSWhich which;
38: PetscTruth trackall;
40: /* function tables */
41: PetscErrorCode (*fcreate[][2])(MPI_Comm,QEP_LINEAR*,Mat*) = {
42: { MatCreateExplicit_QEPLINEAR_N1A, MatCreateExplicit_QEPLINEAR_N1B }, /* N1 */
43: { MatCreateExplicit_QEPLINEAR_N2A, MatCreateExplicit_QEPLINEAR_N2B }, /* N2 */
44: { MatCreateExplicit_QEPLINEAR_S1A, MatCreateExplicit_QEPLINEAR_S1B }, /* S1 */
45: { MatCreateExplicit_QEPLINEAR_S2A, MatCreateExplicit_QEPLINEAR_S2B }, /* S2 */
46: { MatCreateExplicit_QEPLINEAR_H1A, MatCreateExplicit_QEPLINEAR_H1B }, /* H1 */
47: { MatCreateExplicit_QEPLINEAR_H2A, MatCreateExplicit_QEPLINEAR_H2B } /* H2 */
48: };
49: PetscErrorCode (*fmult[][2])(Mat,Vec,Vec) = {
50: { MatMult_QEPLINEAR_N1A, MatMult_QEPLINEAR_N1B },
51: { MatMult_QEPLINEAR_N2A, MatMult_QEPLINEAR_N2B },
52: { MatMult_QEPLINEAR_S1A, MatMult_QEPLINEAR_S1B },
53: { MatMult_QEPLINEAR_S2A, MatMult_QEPLINEAR_S2B },
54: { MatMult_QEPLINEAR_H1A, MatMult_QEPLINEAR_H1B },
55: { MatMult_QEPLINEAR_H2A, MatMult_QEPLINEAR_H2B }
56: };
57: PetscErrorCode (*fgetdiagonal[][2])(Mat,Vec) = {
58: { MatGetDiagonal_QEPLINEAR_N1A, MatGetDiagonal_QEPLINEAR_N1B },
59: { MatGetDiagonal_QEPLINEAR_N2A, MatGetDiagonal_QEPLINEAR_N2B },
60: { MatGetDiagonal_QEPLINEAR_S1A, MatGetDiagonal_QEPLINEAR_S1B },
61: { MatGetDiagonal_QEPLINEAR_S2A, MatGetDiagonal_QEPLINEAR_S2B },
62: { MatGetDiagonal_QEPLINEAR_H1A, MatGetDiagonal_QEPLINEAR_H1B },
63: { MatGetDiagonal_QEPLINEAR_H2A, MatGetDiagonal_QEPLINEAR_H2B }
64: };
67: if (!qep->which) qep->which = QEP_LARGEST_MAGNITUDE;
68: ctx->M = qep->M;
69: ctx->C = qep->C;
70: ctx->K = qep->K;
71: ctx->sfactor = qep->sfactor;
73: if (ctx->A) {
74: MatDestroy(ctx->A);
75: MatDestroy(ctx->B);
76: }
77: if (ctx->x1) {
78: VecDestroy(ctx->x1);
79: VecDestroy(ctx->x2);
80: VecDestroy(ctx->y1);
81: VecDestroy(ctx->y2);
82: }
84: switch (qep->problem_type) {
85: case QEP_GENERAL: i = 0; break;
86: case QEP_HERMITIAN: i = 2; break;
87: case QEP_GYROSCOPIC: i = 4; break;
88: default: SETERRQ(1,"Wrong value of qep->problem_type");
89: }
90: i += ctx->cform-1;
92: if (ctx->explicitmatrix) {
93: ctx->x1 = ctx->x2 = ctx->y1 = ctx->y2 = PETSC_NULL;
94: (*fcreate[i][0])(((PetscObject)qep)->comm,ctx,&ctx->A);
95: (*fcreate[i][1])(((PetscObject)qep)->comm,ctx,&ctx->B);
96: } else {
97: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->x1);
98: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->x2);
99: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->y1);
100: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->y2);
101: MatCreateShell(((PetscObject)qep)->comm,2*qep->nloc,2*qep->nloc,2*qep->n,2*qep->n,ctx,&ctx->A);
102: MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))fmult[i][0]);
103: MatShellSetOperation(ctx->A,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][0]);
104: MatCreateShell(((PetscObject)qep)->comm,2*qep->nloc,2*qep->nloc,2*qep->n,2*qep->n,ctx,&ctx->B);
105: MatShellSetOperation(ctx->B,MATOP_MULT,(void(*)(void))fmult[i][1]);
106: MatShellSetOperation(ctx->B,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][1]);
107: }
109: EPSSetOperators(ctx->eps,ctx->A,ctx->B);
110: EPSSetProblemType(ctx->eps,EPS_GNHEP);
111: switch (qep->which) {
112: case QEP_LARGEST_MAGNITUDE: which = EPS_LARGEST_MAGNITUDE; break;
113: case QEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
114: case QEP_LARGEST_REAL: which = EPS_LARGEST_REAL; break;
115: case QEP_SMALLEST_REAL: which = EPS_SMALLEST_REAL; break;
116: case QEP_LARGEST_IMAGINARY: which = EPS_LARGEST_IMAGINARY; break;
117: case QEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
118: default: SETERRQ(1,"Wrong value of which");
119: }
120: EPSSetWhichEigenpairs(ctx->eps,which);
121: EPSSetLeftVectorsWanted(ctx->eps,qep->leftvecs);
122: EPSSetDimensions(ctx->eps,qep->nev,qep->ncv,qep->mpd);
123: EPSSetTolerances(ctx->eps,qep->tol,qep->max_it);
124: /* Transfer the trackall option from qep to eps */
125: QEPGetTrackAll(qep,&trackall);
126: EPSSetTrackAll(ctx->eps,trackall);
127: if (ctx->setfromoptionscalled) {
128: EPSSetFromOptions(ctx->eps);
129: ctx->setfromoptionscalled = PETSC_FALSE;
130: }
131: EPSSetUp(ctx->eps);
132: EPSGetDimensions(ctx->eps,PETSC_NULL,&qep->ncv,&qep->mpd);
133: EPSGetTolerances(ctx->eps,&qep->tol,&qep->max_it);
134: if (qep->nini>0 || qep->ninil>0) PetscInfo(qep,"Ignoring initial vectors\n");
136: return(0);
137: }
141: /*
142: QEPLinearSelect_Norm - Auxiliary routine that copies the solution of the
143: linear eigenproblem to the QEP object. The eigenvector of the generalized
144: problem is supposed to be
145: z = [ x ]
146: [ l*x ]
147: The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
148: computed residual norm.
149: Finally, x is normalized so that ||x||_2 = 1.
151: If explicitmatrix==PETSC_TRUE then z is partitioned across processors, otherwise x is.
152: */
153: PetscErrorCode QEPLinearSelect_Norm(QEP qep,EPS eps,PetscTruth explicitmatrix)
154: {
156: PetscInt i,start,end,idx;
157: PetscScalar *px;
158: PetscReal rn1,rn2;
159: Vec v0,xr,xi,wr,wi;
160: Mat A;
161: IS isV1,isV2;
162: VecScatter vsV1,vsV2;
163: #if !defined(PETSC_USE_COMPLEX)
164: PetscScalar *py;
165: #endif
166:
168: EPSGetOperators(eps,&A,PETSC_NULL);
169: MatGetVecs(A,&v0,PETSC_NULL);
170: VecDuplicate(v0,&xr);
171: VecDuplicate(v0,&xi);
173: if (explicitmatrix) { /* case 1: x needs to be scattered from the owning processes to the rest */
174: VecDuplicate(qep->V[0],&wr);
175: VecDuplicate(qep->V[0],&wi);
176: VecGetOwnershipRange(qep->V[0],&start,&end);
177: idx = start;
178: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV1);
179: VecScatterCreate(xr,isV1,qep->V[0],PETSC_NULL,&vsV1);
180: idx = start+qep->n;
181: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV2);
182: VecScatterCreate(xr,isV2,qep->V[0],PETSC_NULL,&vsV2);
183: for (i=0;i<qep->nconv;i++) {
184: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
185: qep->eigr[i] *= qep->sfactor;
186: qep->eigi[i] *= qep->sfactor;
187: #if !defined(PETSC_USE_COMPLEX)
188: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
189: VecScatterBegin(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
190: VecScatterEnd(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
191: VecScatterBegin(vsV1,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
192: VecScatterEnd(vsV1,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
193: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
194: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
195: VecCopy(wr,qep->V[i]);
196: VecCopy(wi,qep->V[i+1]);
197: VecScatterBegin(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
198: VecScatterEnd(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
199: VecScatterBegin(vsV2,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
200: VecScatterEnd(vsV2,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
201: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
202: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
203: if (rn1>rn2) {
204: VecCopy(wr,qep->V[i]);
205: VecCopy(wi,qep->V[i+1]);
206: }
207: }
208: else if (qep->eigi[i]==0.0) /* real eigenvalue */
209: #endif
210: {
211: VecScatterBegin(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
212: VecScatterEnd(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
213: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
214: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
215: VecCopy(wr,qep->V[i]);
216: VecScatterBegin(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
217: VecScatterEnd(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
218: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
219: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
220: if (rn1>rn2) {
221: VecCopy(wr,qep->V[i]);
222: }
223: }
224: }
225: ISDestroy(isV1);
226: ISDestroy(isV2);
227: VecScatterDestroy(vsV1);
228: VecScatterDestroy(vsV2);
229: }
230: else { /* case 2: elements of x are already in the right process */
231: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wr);
232: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wi);
233: for (i=0;i<qep->nconv;i++) {
234: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
235: qep->eigr[i] *= qep->sfactor;
236: qep->eigi[i] *= qep->sfactor;
237: #if !defined(PETSC_USE_COMPLEX)
238: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
239: VecGetArray(xr,&px);
240: VecGetArray(xi,&py);
241: VecPlaceArray(wr,px);
242: VecPlaceArray(wi,py);
243: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
244: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
245: VecCopy(wr,qep->V[i]);
246: VecCopy(wi,qep->V[i+1]);
247: VecResetArray(wr);
248: VecResetArray(wi);
249: VecPlaceArray(wr,px+qep->nloc);
250: VecPlaceArray(wi,py+qep->nloc);
251: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
252: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
253: if (rn1>rn2) {
254: VecCopy(wr,qep->V[i]);
255: VecCopy(wi,qep->V[i+1]);
256: }
257: VecResetArray(wr);
258: VecResetArray(wi);
259: VecRestoreArray(xr,&px);
260: VecRestoreArray(xi,&py);
261: }
262: else if (qep->eigi[i]==0.0) /* real eigenvalue */
263: #endif
264: {
265: VecGetArray(xr,&px);
266: VecPlaceArray(wr,px);
267: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
268: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
269: VecCopy(wr,qep->V[i]);
270: VecResetArray(wr);
271: VecPlaceArray(wr,px+qep->nloc);
272: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
273: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
274: if (rn1>rn2) {
275: VecCopy(wr,qep->V[i]);
276: }
277: VecResetArray(wr);
278: VecRestoreArray(xr,&px);
279: }
280: }
281: }
282: VecDestroy(wr);
283: VecDestroy(wi);
284: VecDestroy(xr);
285: VecDestroy(xi);
287: return(0);
288: }
292: /*
293: QEPLinearSelect_Simple - Auxiliary routine that copies the solution of the
294: linear eigenproblem to the QEP object. The eigenvector of the generalized
295: problem is supposed to be
296: z = [ x ]
297: [ l*x ]
298: If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
299: Finally, x is normalized so that ||x||_2 = 1.
301: If explicitmatrix==PETSC_TRUE then z is partitioned across processors, otherwise x is.
302: */
303: PetscErrorCode QEPLinearSelect_Simple(QEP qep,EPS eps,PetscTruth explicitmatrix)
304: {
306: PetscInt i,start,end,offset,idx;
307: PetscScalar *px;
308: Vec v0,xr,xi,w;
309: Mat A;
310: IS isV1,isV2;
311: VecScatter vsV,vsV1,vsV2;
312:
314: EPSGetOperators(eps,&A,PETSC_NULL);
315: MatGetVecs(A,&v0,PETSC_NULL);
316: VecDuplicate(v0,&xr);
317: VecDuplicate(v0,&xi);
319: if (explicitmatrix) { /* case 1: x needs to be scattered from the owning processes to the rest */
320: VecGetOwnershipRange(qep->V[0],&start,&end);
321: idx = start;
322: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV1);
323: VecScatterCreate(xr,isV1,qep->V[0],PETSC_NULL,&vsV1);
324: idx = start+qep->n;
325: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV2);
326: VecScatterCreate(xr,isV2,qep->V[0],PETSC_NULL,&vsV2);
327: for (i=0;i<qep->nconv;i++) {
328: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
329: qep->eigr[i] *= qep->sfactor;
330: qep->eigi[i] *= qep->sfactor;
331: if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) vsV = vsV2;
332: else vsV = vsV1;
333: #if !defined(PETSC_USE_COMPLEX)
334: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
335: VecScatterBegin(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
336: VecScatterEnd(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
337: VecScatterBegin(vsV,xi,qep->V[i+1],INSERT_VALUES,SCATTER_FORWARD);
338: VecScatterEnd(vsV,xi,qep->V[i+1],INSERT_VALUES,SCATTER_FORWARD);
339: SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
340: }
341: else if (qep->eigi[i]==0.0) /* real eigenvalue */
342: #endif
343: {
344: VecScatterBegin(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
345: VecScatterEnd(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
346: SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
347: }
348: }
349: ISDestroy(isV1);
350: ISDestroy(isV2);
351: VecScatterDestroy(vsV1);
352: VecScatterDestroy(vsV2);
353: }
354: else { /* case 2: elements of x are already in the right process */
355: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&w);
356: for (i=0;i<qep->nconv;i++) {
357: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
358: qep->eigr[i] *= qep->sfactor;
359: qep->eigi[i] *= qep->sfactor;
360: if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) offset = qep->nloc;
361: else offset = 0;
362: #if !defined(PETSC_USE_COMPLEX)
363: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
364: VecGetArray(xr,&px);
365: VecPlaceArray(w,px+offset);
366: VecCopy(w,qep->V[i]);
367: VecResetArray(w);
368: VecRestoreArray(xr,&px);
369: VecGetArray(xi,&px);
370: VecPlaceArray(w,px+offset);
371: VecCopy(w,qep->V[i+1]);
372: VecResetArray(w);
373: VecRestoreArray(xi,&px);
374: SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
375: }
376: else if (qep->eigi[i]==0.0) /* real eigenvalue */
377: #endif
378: {
379: VecGetArray(xr,&px);
380: VecPlaceArray(w,px+offset);
381: VecCopy(w,qep->V[i]);
382: VecResetArray(w);
383: VecRestoreArray(xr,&px);
384: SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
385: }
386: }
387: VecDestroy(w);
388: }
389: VecDestroy(xr);
390: VecDestroy(xi);
392: return(0);
393: }
397: PetscErrorCode QEPSolve_LINEAR(QEP qep)
398: {
400: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
401: PetscTruth flg=PETSC_FALSE;
404: EPSSolve(ctx->eps);
405: EPSGetConverged(ctx->eps,&qep->nconv);
406: EPSGetIterationNumber(ctx->eps,&qep->its);
407: EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&qep->reason);
408: EPSGetOperationCounters(ctx->eps,&qep->matvecs,PETSC_NULL,&qep->linits);
409: qep->matvecs *= 2; /* convention: count one matvec for each non-trivial block in A */
410: PetscOptionsGetTruth(((PetscObject)qep)->prefix,"-qep_linear_select_simple",&flg,PETSC_NULL);
411: if (flg) {
412: QEPLinearSelect_Simple(qep,ctx->eps,ctx->explicitmatrix);
413: } else {
414: QEPLinearSelect_Norm(qep,ctx->eps,ctx->explicitmatrix);
415: }
416: return(0);
417: }
421: PetscErrorCode EPSMonitor_QEP_LINEAR(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
422: {
423: PetscInt i;
424: QEP qep = (QEP)ctx;
428: nconv = 0;
429: for (i=0;i<nest;i++) {
430: qep->eigr[i] = eigr[i];
431: qep->eigi[i] = eigi[i];
432: qep->errest[i] = errest[i];
433: if (0.0 < errest[i] && errest[i] < qep->tol) nconv++;
434: }
435: STBackTransform(eps->OP,nest,qep->eigr,qep->eigi);
436: QEPMonitor(qep,its,nconv,qep->eigr,qep->eigi,qep->errest,nest);
437: return(0);
438: }
442: PetscErrorCode QEPSetFromOptions_LINEAR(QEP qep)
443: {
445: PetscTruth set,val;
446: PetscInt i;
447: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
448: ST st;
451: PetscOptionsBegin(((PetscObject)qep)->comm,((PetscObject)qep)->prefix,"LINEAR Quadratic Eigenvalue Problem solver Options","QEP");
453: PetscOptionsInt("-qep_linear_cform","Number of the companion form","QEPLinearSetCompanionForm",ctx->cform,&i,&set);
454: if (set) {
455: QEPLinearSetCompanionForm(qep,i);
456: }
458: PetscOptionsTruth("-qep_linear_explicitmatrix","Use explicit matrix in linearization","QEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
459: if (set) {
460: QEPLinearSetExplicitMatrix(qep,val);
461: }
462: if (!ctx->explicitmatrix) {
463: /* use as default an ST with shell matrix and Jacobi */
464: EPSGetST(ctx->eps,&st);
465: STSetMatMode(st,ST_MATMODE_SHELL);
466: }
468: PetscOptionsEnd();
469: ctx->setfromoptionscalled = PETSC_TRUE;
470: return(0);
471: }
476: PetscErrorCode QEPLinearSetCompanionForm_LINEAR(QEP qep,PetscInt cform)
477: {
478: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
481: if (cform==PETSC_IGNORE) return(0);
482: if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
483: else {
484: if (cform!=1 && cform!=2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
485: ctx->cform = cform;
486: }
487: return(0);
488: }
493: /*@
494: QEPLinearSetCompanionForm - Choose between the two companion forms available
495: for the linearization of the quadratic problem.
497: Collective on QEP
499: Input Parameters:
500: + qep - quadratic eigenvalue solver
501: - cform - 1 or 2 (first or second companion form)
503: Options Database Key:
504: . -qep_linear_cform <int> - Choose the companion form
506: Level: advanced
508: .seealso: QEPLinearGetCompanionForm()
509: @*/
510: PetscErrorCode QEPLinearSetCompanionForm(QEP qep,PetscInt cform)
511: {
512: PetscErrorCode ierr, (*f)(QEP,PetscInt);
516: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetCompanionForm_C",(void (**)())&f);
517: if (f) {
518: (*f)(qep,cform);
519: }
520: return(0);
521: }
526: PetscErrorCode QEPLinearGetCompanionForm_LINEAR(QEP qep,PetscInt *cform)
527: {
528: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
532: *cform = ctx->cform;
533: return(0);
534: }
539: /*@
540: QEPLinearGetCompanionForm - Returns the number of the companion form that
541: will be used for the linearization of the quadratic problem.
543: Not collective
545: Input Parameter:
546: . qep - quadratic eigenvalue solver
548: Output Parameter:
549: . cform - the companion form number (1 or 2)
551: Level: advanced
553: .seealso: QEPLinearSetCompanionForm()
554: @*/
555: PetscErrorCode QEPLinearGetCompanionForm(QEP qep,PetscInt *cform)
556: {
557: PetscErrorCode ierr, (*f)(QEP,PetscInt*);
561: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetCompanionForm_C",(void (**)())&f);
562: if (f) {
563: (*f)(qep,cform);
564: }
565: return(0);
566: }
571: PetscErrorCode QEPLinearSetExplicitMatrix_LINEAR(QEP qep,PetscTruth explicitmatrix)
572: {
573: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
576: ctx->explicitmatrix = explicitmatrix;
577: return(0);
578: }
583: /*@
584: QEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the linearization
585: of the quadratic problem must be built explicitly.
587: Collective on QEP
589: Input Parameters:
590: + qep - quadratic eigenvalue solver
591: - explicit - boolean flag indicating if the matrices are built explicitly
593: Options Database Key:
594: . -qep_linear_explicitmatrix <boolean> - Indicates the boolean flag
596: Level: advanced
598: .seealso: QEPLinearGetExplicitMatrix()
599: @*/
600: PetscErrorCode QEPLinearSetExplicitMatrix(QEP qep,PetscTruth explicitmatrix)
601: {
602: PetscErrorCode ierr, (*f)(QEP,PetscTruth);
606: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetExplicitMatrix_C",(void (**)())&f);
607: if (f) {
608: (*f)(qep,explicitmatrix);
609: }
610: return(0);
611: }
616: PetscErrorCode QEPLinearGetExplicitMatrix_LINEAR(QEP qep,PetscTruth *explicitmatrix)
617: {
618: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
622: *explicitmatrix = ctx->explicitmatrix;
623: return(0);
624: }
629: /*@
630: QEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices A and B
631: for the linearization of the quadratic problem are built explicitly.
633: Not collective
635: Input Parameter:
636: . qep - quadratic eigenvalue solver
638: Output Parameter:
639: . explicitmatrix - the mode flag
641: Level: advanced
643: .seealso: QEPLinearSetExplicitMatrix()
644: @*/
645: PetscErrorCode QEPLinearGetExplicitMatrix(QEP qep,PetscTruth *explicitmatrix)
646: {
647: PetscErrorCode ierr, (*f)(QEP,PetscTruth*);
651: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetExplicitMatrix_C",(void (**)())&f);
652: if (f) {
653: (*f)(qep,explicitmatrix);
654: }
655: return(0);
656: }
661: PetscErrorCode QEPLinearSetEPS_LINEAR(QEP qep,EPS eps)
662: {
663: PetscErrorCode ierr;
664: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
669: PetscObjectReference((PetscObject)eps);
670: EPSDestroy(ctx->eps);
671: ctx->eps = eps;
672: qep->setupcalled = 0;
673: return(0);
674: }
679: /*@
680: QEPLinearSetEPS - Associate an eigensolver object (EPS) to the
681: quadratic eigenvalue solver.
683: Collective on QEP
685: Input Parameters:
686: + qep - quadratic eigenvalue solver
687: - eps - the eigensolver object
689: Level: advanced
691: .seealso: QEPLinearGetEPS()
692: @*/
693: PetscErrorCode QEPLinearSetEPS(QEP qep,EPS eps)
694: {
695: PetscErrorCode ierr, (*f)(QEP,EPS);
699: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetEPS_C",(void (**)())&f);
700: if (f) {
701: (*f)(qep,eps);
702: }
703: return(0);
704: }
709: PetscErrorCode QEPLinearGetEPS_LINEAR(QEP qep,EPS *eps)
710: {
711: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
715: *eps = ctx->eps;
716: return(0);
717: }
722: /*@
723: QEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
724: to the quadratic eigenvalue solver.
726: Not Collective
728: Input Parameter:
729: . qep - quadratic eigenvalue solver
731: Output Parameter:
732: . eps - the eigensolver object
734: Level: advanced
736: .seealso: QEPLinearSetEPS()
737: @*/
738: PetscErrorCode QEPLinearGetEPS(QEP qep,EPS *eps)
739: {
740: PetscErrorCode ierr, (*f)(QEP,EPS*);
744: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetEPS_C",(void (**)())&f);
745: if (f) {
746: (*f)(qep,eps);
747: }
748: return(0);
749: }
753: PetscErrorCode QEPView_LINEAR(QEP qep,PetscViewer viewer)
754: {
755: PetscErrorCode ierr;
756: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
759: if (ctx->explicitmatrix) {
760: PetscViewerASCIIPrintf(viewer,"linearized matrices: explicit\n");
761: } else {
762: PetscViewerASCIIPrintf(viewer,"linearized matrices: implicit\n");
763: }
764: PetscViewerASCIIPrintf(viewer,"companion form: %d\n",ctx->cform);
765: EPSView(ctx->eps,viewer);
766: return(0);
767: }
771: PetscErrorCode QEPDestroy_LINEAR(QEP qep)
772: {
773: PetscErrorCode ierr;
774: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
777: EPSDestroy(ctx->eps);
778: if (ctx->A) {
779: MatDestroy(ctx->A);
780: MatDestroy(ctx->B);
781: }
782: if (ctx->x1) {
783: VecDestroy(ctx->x1);
784: VecDestroy(ctx->x2);
785: VecDestroy(ctx->y1);
786: VecDestroy(ctx->y2);
787: }
789: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","",PETSC_NULL);
790: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","",PETSC_NULL);
791: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","",PETSC_NULL);
792: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","",PETSC_NULL);
793: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","",PETSC_NULL);
794: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","",PETSC_NULL);
796: QEPDestroy_Default(qep);
797: return(0);
798: }
803: PetscErrorCode QEPCreate_LINEAR(QEP qep)
804: {
806: QEP_LINEAR *ctx;
809: PetscNew(QEP_LINEAR,&ctx);
810: PetscLogObjectMemory(qep,sizeof(QEP_LINEAR));
811: qep->data = (void *)ctx;
812: qep->ops->solve = QEPSolve_LINEAR;
813: qep->ops->setup = QEPSetUp_LINEAR;
814: qep->ops->setfromoptions = QEPSetFromOptions_LINEAR;
815: qep->ops->destroy = QEPDestroy_LINEAR;
816: qep->ops->view = QEPView_LINEAR;
817: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","QEPLinearSetCompanionForm_LINEAR",QEPLinearSetCompanionForm_LINEAR);
818: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","QEPLinearGetCompanionForm_LINEAR",QEPLinearGetCompanionForm_LINEAR);
819: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","QEPLinearSetEPS_LINEAR",QEPLinearSetEPS_LINEAR);
820: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","QEPLinearGetEPS_LINEAR",QEPLinearGetEPS_LINEAR);
821: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","QEPLinearSetExplicitMatrix_LINEAR",QEPLinearSetExplicitMatrix_LINEAR);
822: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","QEPLinearGetExplicitMatrix_LINEAR",QEPLinearGetExplicitMatrix_LINEAR);
824: EPSCreate(((PetscObject)qep)->comm,&ctx->eps);
825: EPSSetOptionsPrefix(ctx->eps,((PetscObject)qep)->prefix);
826: EPSAppendOptionsPrefix(ctx->eps,"qep_");
827: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)qep,1);
828: PetscLogObjectParent(qep,ctx->eps);
829: EPSSetIP(ctx->eps,qep->ip);
830: EPSMonitorSet(ctx->eps,EPSMonitor_QEP_LINEAR,qep,PETSC_NULL);
831: ctx->explicitmatrix = PETSC_FALSE;
832: ctx->cform = 1;
833: ctx->A = PETSC_NULL;
834: ctx->B = PETSC_NULL;
835: ctx->x1 = PETSC_NULL;
836: ctx->x2 = PETSC_NULL;
837: ctx->y1 = PETSC_NULL;
838: ctx->y2 = PETSC_NULL;
839: ctx->setfromoptionscalled = PETSC_FALSE;
840: return(0);
841: }