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 == PETSC_TRUE) {
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,*py;
158: PetscReal rn1,rn2;
159: Vec v0,xr,xi,wr,wi;
160: Mat A;
161: IS isV1,isV2;
162: VecScatter vsV1,vsV2;
163:
165: EPSGetOperators(eps,&A,PETSC_NULL);
166: MatGetVecs(A,&v0,PETSC_NULL);
167: VecDuplicate(v0,&xr);
168: VecDuplicate(v0,&xi);
170: if (explicitmatrix) { /* case 1: x needs to be scattered from the owning processes to the rest */
171: VecDuplicate(qep->V[0],&wr);
172: VecDuplicate(qep->V[0],&wi);
173: VecGetOwnershipRange(qep->V[0],&start,&end);
174: idx = start;
175: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV1);
176: VecScatterCreate(xr,isV1,qep->V[0],PETSC_NULL,&vsV1);
177: idx = start+qep->n;
178: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV2);
179: VecScatterCreate(xr,isV2,qep->V[0],PETSC_NULL,&vsV2);
180: for (i=0;i<qep->nconv;i++) {
181: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
182: qep->eigr[i] *= qep->sfactor;
183: qep->eigi[i] *= qep->sfactor;
184: #if !defined(PETSC_USE_COMPLEX)
185: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
186: VecScatterBegin(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
187: VecScatterEnd(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
188: VecScatterBegin(vsV1,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
189: VecScatterEnd(vsV1,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
190: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
191: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
192: VecCopy(wr,qep->V[i]);
193: VecCopy(wi,qep->V[i+1]);
194: VecScatterBegin(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
195: VecScatterEnd(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
196: VecScatterBegin(vsV2,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
197: VecScatterEnd(vsV2,xi,wi,INSERT_VALUES,SCATTER_FORWARD);
198: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
199: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
200: if (rn1>rn2) {
201: VecCopy(wr,qep->V[i]);
202: VecCopy(wi,qep->V[i+1]);
203: }
204: }
205: else if (qep->eigi[i]==0.0) /* real eigenvalue */
206: #endif
207: {
208: VecScatterBegin(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
209: VecScatterEnd(vsV1,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
210: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
211: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
212: VecCopy(wr,qep->V[i]);
213: VecScatterBegin(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
214: VecScatterEnd(vsV2,xr,wr,INSERT_VALUES,SCATTER_FORWARD);
215: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
216: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
217: if (rn1>rn2) {
218: VecCopy(wr,qep->V[i]);
219: }
220: }
221: }
222: ISDestroy(isV1);
223: ISDestroy(isV2);
224: VecScatterDestroy(vsV1);
225: VecScatterDestroy(vsV2);
226: }
227: else { /* case 2: elements of x are already in the right process */
228: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wr);
229: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wi);
230: for (i=0;i<qep->nconv;i++) {
231: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
232: qep->eigr[i] *= qep->sfactor;
233: qep->eigi[i] *= qep->sfactor;
234: #if !defined(PETSC_USE_COMPLEX)
235: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
236: VecGetArray(xr,&px);
237: VecGetArray(xi,&py);
238: VecPlaceArray(wr,px);
239: VecPlaceArray(wi,py);
240: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
241: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
242: VecCopy(wr,qep->V[i]);
243: VecCopy(wi,qep->V[i+1]);
244: VecResetArray(wr);
245: VecResetArray(wi);
246: VecPlaceArray(wr,px+qep->nloc);
247: VecPlaceArray(wi,py+qep->nloc);
248: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
249: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
250: if (rn1>rn2) {
251: VecCopy(wr,qep->V[i]);
252: VecCopy(wi,qep->V[i+1]);
253: }
254: VecResetArray(wr);
255: VecResetArray(wi);
256: VecRestoreArray(xr,&px);
257: VecRestoreArray(xi,&py);
258: }
259: else if (qep->eigi[i]==0.0) /* real eigenvalue */
260: #endif
261: {
262: VecGetArray(xr,&px);
263: VecPlaceArray(wr,px);
264: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
265: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
266: VecCopy(wr,qep->V[i]);
267: VecResetArray(wr);
268: VecPlaceArray(wr,px+qep->nloc);
269: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
270: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
271: if (rn1>rn2) {
272: VecCopy(wr,qep->V[i]);
273: }
274: VecResetArray(wr);
275: VecRestoreArray(xr,&px);
276: }
277: }
278: }
279: VecDestroy(wr);
280: VecDestroy(wi);
281: VecDestroy(xr);
282: VecDestroy(xi);
284: return(0);
285: }
289: /*
290: QEPLinearSelect_Simple - Auxiliary routine that copies the solution of the
291: linear eigenproblem to the QEP object. The eigenvector of the generalized
292: problem is supposed to be
293: z = [ x ]
294: [ l*x ]
295: If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
296: Finally, x is normalized so that ||x||_2 = 1.
298: If explicitmatrix==PETSC_TRUE then z is partitioned across processors, otherwise x is.
299: */
300: PetscErrorCode QEPLinearSelect_Simple(QEP qep,EPS eps,PetscTruth explicitmatrix)
301: {
303: PetscInt i,start,end,offset,idx;
304: PetscScalar *px;
305: Vec v0,xr,xi,w;
306: Mat A;
307: IS isV1,isV2;
308: VecScatter vsV,vsV1,vsV2;
309:
311: EPSGetOperators(eps,&A,PETSC_NULL);
312: MatGetVecs(A,&v0,PETSC_NULL);
313: VecDuplicate(v0,&xr);
314: VecDuplicate(v0,&xi);
316: if (explicitmatrix) { /* case 1: x needs to be scattered from the owning processes to the rest */
317: VecGetOwnershipRange(qep->V[0],&start,&end);
318: idx = start;
319: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV1);
320: VecScatterCreate(xr,isV1,qep->V[0],PETSC_NULL,&vsV1);
321: idx = start+qep->n;
322: ISCreateBlock(((PetscObject)qep)->comm,end-start,1,&idx,&isV2);
323: VecScatterCreate(xr,isV2,qep->V[0],PETSC_NULL,&vsV2);
324: for (i=0;i<qep->nconv;i++) {
325: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
326: qep->eigr[i] *= qep->sfactor;
327: qep->eigi[i] *= qep->sfactor;
328: if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) vsV = vsV2;
329: else vsV = vsV1;
330: #if !defined(PETSC_USE_COMPLEX)
331: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
332: VecScatterBegin(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
333: VecScatterEnd(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
334: VecScatterBegin(vsV,xi,qep->V[i+1],INSERT_VALUES,SCATTER_FORWARD);
335: VecScatterEnd(vsV,xi,qep->V[i+1],INSERT_VALUES,SCATTER_FORWARD);
336: SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
337: }
338: else if (qep->eigi[i]==0.0) /* real eigenvalue */
339: #endif
340: {
341: VecScatterBegin(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
342: VecScatterEnd(vsV,xr,qep->V[i],INSERT_VALUES,SCATTER_FORWARD);
343: SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
344: }
345: }
346: ISDestroy(isV1);
347: ISDestroy(isV2);
348: VecScatterDestroy(vsV1);
349: VecScatterDestroy(vsV2);
350: }
351: else { /* case 2: elements of x are already in the right process */
352: VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&w);
353: for (i=0;i<qep->nconv;i++) {
354: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
355: qep->eigr[i] *= qep->sfactor;
356: qep->eigi[i] *= qep->sfactor;
357: if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) offset = qep->nloc;
358: else offset = 0;
359: #if !defined(PETSC_USE_COMPLEX)
360: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
361: VecGetArray(xr,&px);
362: VecPlaceArray(w,px+offset);
363: VecCopy(w,qep->V[i]);
364: VecResetArray(w);
365: VecRestoreArray(xr,&px);
366: VecGetArray(xi,&px);
367: VecPlaceArray(w,px+offset);
368: VecCopy(w,qep->V[i+1]);
369: VecResetArray(w);
370: VecRestoreArray(xi,&px);
371: SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
372: }
373: else if (qep->eigi[i]==0.0) /* real eigenvalue */
374: #endif
375: {
376: VecGetArray(xr,&px);
377: VecPlaceArray(w,px+offset);
378: VecCopy(w,qep->V[i]);
379: VecResetArray(w);
380: VecRestoreArray(xr,&px);
381: SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
382: }
383: }
384: VecDestroy(w);
385: }
386: VecDestroy(xr);
387: VecDestroy(xi);
389: return(0);
390: }
394: PetscErrorCode QEPSolve_LINEAR(QEP qep)
395: {
397: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
398: PetscTruth flg=PETSC_FALSE;
401: EPSSolve(ctx->eps);
402: EPSGetConverged(ctx->eps,&qep->nconv);
403: EPSGetIterationNumber(ctx->eps,&qep->its);
404: EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&qep->reason);
405: EPSGetOperationCounters(ctx->eps,&qep->matvecs,PETSC_NULL,&qep->linits);
406: qep->matvecs *= 2; /* convention: count one matvec for each non-trivial block in A */
407: PetscOptionsGetTruth(((PetscObject)qep)->prefix,"-qep_linear_select_simple",&flg,PETSC_NULL);
408: if (flg) {
409: QEPLinearSelect_Simple(qep,ctx->eps,ctx->explicitmatrix);
410: } else {
411: QEPLinearSelect_Norm(qep,ctx->eps,ctx->explicitmatrix);
412: }
413: return(0);
414: }
418: PetscErrorCode EPSMonitor_QEP_LINEAR(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
419: {
420: PetscInt i;
421: QEP qep = (QEP)ctx;
425: nconv = 0;
426: for (i=0;i<nest;i++) {
427: qep->eigr[i] = eigr[i];
428: qep->eigi[i] = eigi[i];
429: qep->errest[i] = errest[i];
430: if (0.0 < errest[i] && errest[i] < qep->tol) nconv++;
431: }
432: STBackTransform(eps->OP,nest,qep->eigr,qep->eigi);
433: QEPMonitor(qep,its,nconv,qep->eigr,qep->eigi,qep->errest,nest);
434: return(0);
435: }
439: PetscErrorCode QEPSetFromOptions_LINEAR(QEP qep)
440: {
442: PetscTruth set,val;
443: PetscInt i;
444: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
445: ST st;
448: PetscOptionsBegin(((PetscObject)qep)->comm,((PetscObject)qep)->prefix,"LINEAR Quadratic Eigenvalue Problem solver Options","QEP");
450: PetscOptionsInt("-qep_linear_cform","Number of the companion form","QEPLinearSetCompanionForm",ctx->cform,&i,&set);
451: if (set) {
452: QEPLinearSetCompanionForm(qep,i);
453: }
455: PetscOptionsTruth("-qep_linear_explicitmatrix","Use explicit matrix in linearization","QEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
456: if (set) {
457: QEPLinearSetExplicitMatrix(qep,val);
458: }
459: if (!ctx->explicitmatrix) {
460: /* use as default an ST with shell matrix and Jacobi */
461: EPSGetST(ctx->eps,&st);
462: STSetMatMode(st,ST_MATMODE_SHELL);
463: }
465: PetscOptionsEnd();
466: ctx->setfromoptionscalled = PETSC_TRUE;
467: return(0);
468: }
473: PetscErrorCode QEPLinearSetCompanionForm_LINEAR(QEP qep,PetscInt cform)
474: {
475: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
478: if (cform==PETSC_IGNORE) return(0);
479: if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
480: else {
481: if (cform!=1 && cform!=2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
482: ctx->cform = cform;
483: }
484: return(0);
485: }
490: /*@
491: QEPLinearSetCompanionForm - Choose between the two companion forms available
492: for the linearization of the quadratic problem.
494: Collective on QEP
496: Input Parameters:
497: + qep - quadratic eigenvalue solver
498: - cform - 1 or 2 (first or second companion form)
500: Options Database Key:
501: . -qep_linear_cform <int> - Choose the companion form
503: Level: advanced
505: .seealso: QEPLinearGetCompanionForm()
506: @*/
507: PetscErrorCode QEPLinearSetCompanionForm(QEP qep,PetscInt cform)
508: {
509: PetscErrorCode ierr, (*f)(QEP,PetscInt);
513: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetCompanionForm_C",(void (**)())&f);
514: if (f) {
515: (*f)(qep,cform);
516: }
517: return(0);
518: }
523: PetscErrorCode QEPLinearGetCompanionForm_LINEAR(QEP qep,PetscInt *cform)
524: {
525: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
529: *cform = ctx->cform;
530: return(0);
531: }
536: /*@
537: QEPLinearGetCompanionForm - Returns the number of the companion form that
538: will be used for the linearization of the quadratic problem.
540: Not collective
542: Input Parameter:
543: . qep - quadratic eigenvalue solver
545: Output Parameter:
546: . cform - the companion form number (1 or 2)
548: Level: advanced
550: .seealso: QEPLinearSetCompanionForm()
551: @*/
552: PetscErrorCode QEPLinearGetCompanionForm(QEP qep,PetscInt *cform)
553: {
554: PetscErrorCode ierr, (*f)(QEP,PetscInt*);
558: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetCompanionForm_C",(void (**)())&f);
559: if (f) {
560: (*f)(qep,cform);
561: }
562: return(0);
563: }
568: PetscErrorCode QEPLinearSetExplicitMatrix_LINEAR(QEP qep,PetscTruth explicitmatrix)
569: {
570: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
573: ctx->explicitmatrix = explicitmatrix;
574: return(0);
575: }
580: /*@
581: QEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the linearization
582: of the quadratic problem must be built explicitly.
584: Collective on QEP
586: Input Parameters:
587: + qep - quadratic eigenvalue solver
588: - explicit - boolean flag indicating if the matrices are built explicitly
590: Options Database Key:
591: . -qep_linear_explicitmatrix <boolean> - Indicates the boolean flag
593: Level: advanced
595: .seealso: QEPLinearGetExplicitMatrix()
596: @*/
597: PetscErrorCode QEPLinearSetExplicitMatrix(QEP qep,PetscTruth explicitmatrix)
598: {
599: PetscErrorCode ierr, (*f)(QEP,PetscTruth);
603: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetExplicitMatrix_C",(void (**)())&f);
604: if (f) {
605: (*f)(qep,explicitmatrix);
606: }
607: return(0);
608: }
613: PetscErrorCode QEPLinearGetExplicitMatrix_LINEAR(QEP qep,PetscTruth *explicitmatrix)
614: {
615: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
619: *explicitmatrix = ctx->explicitmatrix;
620: return(0);
621: }
626: /*@
627: QEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices A and B
628: for the linearization of the quadratic problem are built explicitly.
630: Not collective
632: Input Parameter:
633: . qep - quadratic eigenvalue solver
635: Output Parameter:
636: . explicitmatrix - the mode flag
638: Level: advanced
640: .seealso: QEPLinearSetExplicitMatrix()
641: @*/
642: PetscErrorCode QEPLinearGetExplicitMatrix(QEP qep,PetscTruth *explicitmatrix)
643: {
644: PetscErrorCode ierr, (*f)(QEP,PetscTruth*);
648: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetExplicitMatrix_C",(void (**)())&f);
649: if (f) {
650: (*f)(qep,explicitmatrix);
651: }
652: return(0);
653: }
658: PetscErrorCode QEPLinearSetEPS_LINEAR(QEP qep,EPS eps)
659: {
660: PetscErrorCode ierr;
661: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
666: PetscObjectReference((PetscObject)eps);
667: EPSDestroy(ctx->eps);
668: ctx->eps = eps;
669: qep->setupcalled = 0;
670: return(0);
671: }
676: /*@
677: QEPLinearSetEPS - Associate an eigensolver object (EPS) to the
678: quadratic eigenvalue solver.
680: Collective on QEP
682: Input Parameters:
683: + qep - quadratic eigenvalue solver
684: - eps - the eigensolver object
686: Level: advanced
688: .seealso: QEPLinearGetEPS()
689: @*/
690: PetscErrorCode QEPLinearSetEPS(QEP qep,EPS eps)
691: {
692: PetscErrorCode ierr, (*f)(QEP,EPS);
696: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearSetEPS_C",(void (**)())&f);
697: if (f) {
698: (*f)(qep,eps);
699: }
700: return(0);
701: }
706: PetscErrorCode QEPLinearGetEPS_LINEAR(QEP qep,EPS *eps)
707: {
708: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
712: *eps = ctx->eps;
713: return(0);
714: }
719: /*@
720: QEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
721: to the quadratic eigenvalue solver.
723: Not Collective
725: Input Parameter:
726: . qep - quadratic eigenvalue solver
728: Output Parameter:
729: . eps - the eigensolver object
731: Level: advanced
733: .seealso: QEPLinearSetEPS()
734: @*/
735: PetscErrorCode QEPLinearGetEPS(QEP qep,EPS *eps)
736: {
737: PetscErrorCode ierr, (*f)(QEP,EPS*);
741: PetscObjectQueryFunction((PetscObject)qep,"QEPLinearGetEPS_C",(void (**)())&f);
742: if (f) {
743: (*f)(qep,eps);
744: }
745: return(0);
746: }
750: PetscErrorCode QEPView_LINEAR(QEP qep,PetscViewer viewer)
751: {
752: PetscErrorCode ierr;
753: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
756: if (ctx->explicitmatrix) {
757: PetscViewerASCIIPrintf(viewer,"linearized matrices: explicit\n");
758: } else {
759: PetscViewerASCIIPrintf(viewer,"linearized matrices: implicit\n");
760: }
761: PetscViewerASCIIPrintf(viewer,"companion form: %d\n",ctx->cform);
762: EPSView(ctx->eps,viewer);
763: return(0);
764: }
768: PetscErrorCode QEPDestroy_LINEAR(QEP qep)
769: {
770: PetscErrorCode ierr;
771: QEP_LINEAR *ctx = (QEP_LINEAR *)qep->data;
774: EPSDestroy(ctx->eps);
775: if (ctx->A) {
776: MatDestroy(ctx->A);
777: MatDestroy(ctx->B);
778: }
779: if (ctx->x1) {
780: VecDestroy(ctx->x1);
781: VecDestroy(ctx->x2);
782: VecDestroy(ctx->y1);
783: VecDestroy(ctx->y2);
784: }
786: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","",PETSC_NULL);
787: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","",PETSC_NULL);
788: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","",PETSC_NULL);
789: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","",PETSC_NULL);
790: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","",PETSC_NULL);
791: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","",PETSC_NULL);
793: QEPDestroy_Default(qep);
794: return(0);
795: }
800: PetscErrorCode QEPCreate_LINEAR(QEP qep)
801: {
803: QEP_LINEAR *ctx;
806: PetscNew(QEP_LINEAR,&ctx);
807: PetscLogObjectMemory(qep,sizeof(QEP_LINEAR));
808: qep->data = (void *)ctx;
809: qep->ops->solve = QEPSolve_LINEAR;
810: qep->ops->setup = QEPSetUp_LINEAR;
811: qep->ops->setfromoptions = QEPSetFromOptions_LINEAR;
812: qep->ops->destroy = QEPDestroy_LINEAR;
813: qep->ops->view = QEPView_LINEAR;
814: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","QEPLinearSetCompanionForm_LINEAR",QEPLinearSetCompanionForm_LINEAR);
815: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","QEPLinearGetCompanionForm_LINEAR",QEPLinearGetCompanionForm_LINEAR);
816: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","QEPLinearSetEPS_LINEAR",QEPLinearSetEPS_LINEAR);
817: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","QEPLinearGetEPS_LINEAR",QEPLinearGetEPS_LINEAR);
818: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","QEPLinearSetExplicitMatrix_LINEAR",QEPLinearSetExplicitMatrix_LINEAR);
819: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","QEPLinearGetExplicitMatrix_LINEAR",QEPLinearGetExplicitMatrix_LINEAR);
821: EPSCreate(((PetscObject)qep)->comm,&ctx->eps);
822: EPSSetOptionsPrefix(ctx->eps,((PetscObject)qep)->prefix);
823: EPSAppendOptionsPrefix(ctx->eps,"qep_");
824: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)qep,1);
825: PetscLogObjectParent(qep,ctx->eps);
826: EPSSetIP(ctx->eps,qep->ip);
827: EPSMonitorSet(ctx->eps,EPSMonitor_QEP_LINEAR,qep,PETSC_NULL);
828: ctx->explicitmatrix = PETSC_FALSE;
829: ctx->cform = 1;
830: ctx->A = PETSC_NULL;
831: ctx->B = PETSC_NULL;
832: ctx->x1 = PETSC_NULL;
833: ctx->x2 = PETSC_NULL;
834: ctx->y1 = PETSC_NULL;
835: ctx->y2 = PETSC_NULL;
836: ctx->setfromoptionscalled = PETSC_FALSE;
837: return(0);
838: }