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-2012, Universitat 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 <slepc-private/qepimpl.h> /*I "slepcqep.h" I*/
26: #include <slepc-private/epsimpl.h> /*I "slepceps.h" I*/
27: #include linearp.h
31: PetscErrorCode QEPSetUp_Linear(QEP qep)
32: {
34: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
35: PetscInt i=0;
36: EPSWhich which;
37: PetscBool trackall;
38: /* function tables */
39: PetscErrorCode (*fcreate[][2])(MPI_Comm,QEP_LINEAR*,Mat*) = {
40: { MatCreateExplicit_Linear_N1A, MatCreateExplicit_Linear_N1B }, /* N1 */
41: { MatCreateExplicit_Linear_N2A, MatCreateExplicit_Linear_N2B }, /* N2 */
42: { MatCreateExplicit_Linear_S1A, MatCreateExplicit_Linear_S1B }, /* S1 */
43: { MatCreateExplicit_Linear_S2A, MatCreateExplicit_Linear_S2B }, /* S2 */
44: { MatCreateExplicit_Linear_H1A, MatCreateExplicit_Linear_H1B }, /* H1 */
45: { MatCreateExplicit_Linear_H2A, MatCreateExplicit_Linear_H2B } /* H2 */
46: };
47: PetscErrorCode (*fmult[][2])(Mat,Vec,Vec) = {
48: { MatMult_Linear_N1A, MatMult_Linear_N1B },
49: { MatMult_Linear_N2A, MatMult_Linear_N2B },
50: { MatMult_Linear_S1A, MatMult_Linear_S1B },
51: { MatMult_Linear_S2A, MatMult_Linear_S2B },
52: { MatMult_Linear_H1A, MatMult_Linear_H1B },
53: { MatMult_Linear_H2A, MatMult_Linear_H2B }
54: };
55: PetscErrorCode (*fgetdiagonal[][2])(Mat,Vec) = {
56: { MatGetDiagonal_Linear_N1A, MatGetDiagonal_Linear_N1B },
57: { MatGetDiagonal_Linear_N2A, MatGetDiagonal_Linear_N2B },
58: { MatGetDiagonal_Linear_S1A, MatGetDiagonal_Linear_S1B },
59: { MatGetDiagonal_Linear_S2A, MatGetDiagonal_Linear_S2B },
60: { MatGetDiagonal_Linear_H1A, MatGetDiagonal_Linear_H1B },
61: { MatGetDiagonal_Linear_H2A, MatGetDiagonal_Linear_H2B }
62: };
65: if (!qep->which) qep->which = QEP_LARGEST_MAGNITUDE;
66: ctx->M = qep->M;
67: ctx->C = qep->C;
68: ctx->K = qep->K;
69: ctx->sfactor = qep->sfactor;
71: MatDestroy(&ctx->A);
72: MatDestroy(&ctx->B);
73: VecDestroy(&ctx->x1);
74: VecDestroy(&ctx->x2);
75: VecDestroy(&ctx->y1);
76: VecDestroy(&ctx->y2);
78: switch (qep->problem_type) {
79: case QEP_GENERAL: i = 0; break;
80: case QEP_HERMITIAN: i = 2; break;
81: case QEP_GYROSCOPIC: i = 4; break;
82: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
83: }
84: i += ctx->cform-1;
86: if (ctx->explicitmatrix) {
87: ctx->x1 = ctx->x2 = ctx->y1 = ctx->y2 = PETSC_NULL;
88: (*fcreate[i][0])(((PetscObject)qep)->comm,ctx,&ctx->A);
89: (*fcreate[i][1])(((PetscObject)qep)->comm,ctx,&ctx->B);
90: } else {
91: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&ctx->x1);
92: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&ctx->x2);
93: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&ctx->y1);
94: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&ctx->y2);
95: MatCreateShell(((PetscObject)qep)->comm,2*qep->nloc,2*qep->nloc,2*qep->n,2*qep->n,ctx,&ctx->A);
96: MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))fmult[i][0]);
97: MatShellSetOperation(ctx->A,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][0]);
98: MatCreateShell(((PetscObject)qep)->comm,2*qep->nloc,2*qep->nloc,2*qep->n,2*qep->n,ctx,&ctx->B);
99: MatShellSetOperation(ctx->B,MATOP_MULT,(void(*)(void))fmult[i][1]);
100: MatShellSetOperation(ctx->B,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][1]);
101: }
103: EPSSetOperators(ctx->eps,ctx->A,ctx->B);
104: if (qep->problem_type==QEP_HERMITIAN) {
105: EPSSetProblemType(ctx->eps,EPS_GHIEP);
106: } else {
107: EPSSetProblemType(ctx->eps,EPS_GNHEP);
108: }
109: switch (qep->which) {
110: case QEP_LARGEST_MAGNITUDE: which = EPS_LARGEST_MAGNITUDE; break;
111: case QEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
112: case QEP_LARGEST_REAL: which = EPS_LARGEST_REAL; break;
113: case QEP_SMALLEST_REAL: which = EPS_SMALLEST_REAL; break;
114: case QEP_LARGEST_IMAGINARY: which = EPS_LARGEST_IMAGINARY; break;
115: case QEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
116: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of which");
117: }
118: EPSSetWhichEigenpairs(ctx->eps,which);
119: EPSSetLeftVectorsWanted(ctx->eps,qep->leftvecs);
120: EPSSetDimensions(ctx->eps,qep->nev,qep->ncv,qep->mpd);
121: EPSSetTolerances(ctx->eps,qep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:qep->tol/10.0,qep->max_it);
122: /* Transfer the trackall option from qep to eps */
123: QEPGetTrackAll(qep,&trackall);
124: EPSSetTrackAll(ctx->eps,trackall);
125: if (ctx->setfromoptionscalled) {
126: EPSSetFromOptions(ctx->eps);
127: ctx->setfromoptionscalled = PETSC_FALSE;
128: }
129: EPSSetUp(ctx->eps);
130: EPSGetDimensions(ctx->eps,PETSC_NULL,&qep->ncv,&qep->mpd);
131: EPSGetTolerances(ctx->eps,PETSC_NULL,&qep->max_it);
132: if (qep->nini>0 || qep->ninil>0) { PetscInfo(qep,"Ignoring initial vectors\n"); }
133: QEPAllocateSolution(qep);
134: return(0);
135: }
139: /*
140: QEPLinearSelect_Norm - Auxiliary routine that copies the solution of the
141: linear eigenproblem to the QEP object. The eigenvector of the generalized
142: problem is supposed to be
143: z = [ x ]
144: [ l*x ]
145: The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
146: computed residual norm.
147: Finally, x is normalized so that ||x||_2 = 1.
148: */
149: PetscErrorCode QEPLinearSelect_Norm(QEP qep,EPS eps)
150: {
152: PetscInt i;
153: PetscScalar *px;
154: PetscReal rn1,rn2;
155: Vec xr,xi,wr,wi;
156: Mat A;
157: #if !defined(PETSC_USE_COMPLEX)
158: PetscScalar *py;
159: #endif
160:
162: EPSGetOperators(eps,&A,PETSC_NULL);
163: MatGetVecs(A,&xr,PETSC_NULL);
164: VecDuplicate(xr,&xi);
165: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&wr);
166: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&wi);
167: for (i=0;i<qep->nconv;i++) {
168: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
169: qep->eigr[i] *= qep->sfactor;
170: qep->eigi[i] *= qep->sfactor;
171: #if !defined(PETSC_USE_COMPLEX)
172: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
173: VecGetArray(xr,&px);
174: VecGetArray(xi,&py);
175: VecPlaceArray(wr,px);
176: VecPlaceArray(wi,py);
177: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
178: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
179: VecCopy(wr,qep->V[i]);
180: VecCopy(wi,qep->V[i+1]);
181: VecResetArray(wr);
182: VecResetArray(wi);
183: VecPlaceArray(wr,px+qep->nloc);
184: VecPlaceArray(wi,py+qep->nloc);
185: SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
186: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
187: if (rn1>rn2) {
188: VecCopy(wr,qep->V[i]);
189: VecCopy(wi,qep->V[i+1]);
190: }
191: VecResetArray(wr);
192: VecResetArray(wi);
193: VecRestoreArray(xr,&px);
194: VecRestoreArray(xi,&py);
195: }
196: else if (qep->eigi[i]==0.0) /* real eigenvalue */
197: #endif
198: {
199: VecGetArray(xr,&px);
200: VecPlaceArray(wr,px);
201: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
202: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
203: VecCopy(wr,qep->V[i]);
204: VecResetArray(wr);
205: VecPlaceArray(wr,px+qep->nloc);
206: SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
207: QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
208: if (rn1>rn2) {
209: VecCopy(wr,qep->V[i]);
210: }
211: VecResetArray(wr);
212: VecRestoreArray(xr,&px);
213: }
214: }
215: VecDestroy(&wr);
216: VecDestroy(&wi);
217: VecDestroy(&xr);
218: VecDestroy(&xi);
219: return(0);
220: }
224: /*
225: QEPLinearSelect_Simple - Auxiliary routine that copies the solution of the
226: linear eigenproblem to the QEP object. The eigenvector of the generalized
227: problem is supposed to be
228: z = [ x ]
229: [ l*x ]
230: If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
231: Finally, x is normalized so that ||x||_2 = 1.
232: */
233: PetscErrorCode QEPLinearSelect_Simple(QEP qep,EPS eps)
234: {
236: PetscInt i,offset;
237: PetscScalar *px;
238: Vec xr,xi,w;
239: Mat A;
240:
242: EPSGetOperators(eps,&A,PETSC_NULL);
243: MatGetVecs(A,&xr,PETSC_NULL);
244: VecDuplicate(xr,&xi);
245: VecCreateMPIWithArray(((PetscObject)qep)->comm,1,qep->nloc,qep->n,PETSC_NULL,&w);
246: for (i=0;i<qep->nconv;i++) {
247: EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
248: qep->eigr[i] *= qep->sfactor;
249: qep->eigi[i] *= qep->sfactor;
250: if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) offset = qep->nloc;
251: else offset = 0;
252: #if !defined(PETSC_USE_COMPLEX)
253: if (qep->eigi[i]>0.0) { /* first eigenvalue of a complex conjugate pair */
254: VecGetArray(xr,&px);
255: VecPlaceArray(w,px+offset);
256: VecCopy(w,qep->V[i]);
257: VecResetArray(w);
258: VecRestoreArray(xr,&px);
259: VecGetArray(xi,&px);
260: VecPlaceArray(w,px+offset);
261: VecCopy(w,qep->V[i+1]);
262: VecResetArray(w);
263: VecRestoreArray(xi,&px);
264: SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
265: }
266: else if (qep->eigi[i]==0.0) /* real eigenvalue */
267: #endif
268: {
269: VecGetArray(xr,&px);
270: VecPlaceArray(w,px+offset);
271: VecCopy(w,qep->V[i]);
272: VecResetArray(w);
273: VecRestoreArray(xr,&px);
274: SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
275: }
276: }
277: VecDestroy(&w);
278: VecDestroy(&xr);
279: VecDestroy(&xi);
280: return(0);
281: }
285: PetscErrorCode QEPSolve_Linear(QEP qep)
286: {
288: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
289: PetscBool flg=PETSC_FALSE;
292: EPSSolve(ctx->eps);
293: EPSGetConverged(ctx->eps,&qep->nconv);
294: EPSGetIterationNumber(ctx->eps,&qep->its);
295: EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&qep->reason);
296: EPSGetOperationCounters(ctx->eps,&qep->matvecs,PETSC_NULL,&qep->linits);
297: qep->matvecs *= 2; /* convention: count one matvec for each non-trivial block in A */
298: PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_linear_select_simple",&flg,PETSC_NULL);
299: if (flg) {
300: QEPLinearSelect_Simple(qep,ctx->eps);
301: } else {
302: QEPLinearSelect_Norm(qep,ctx->eps);
303: }
304: return(0);
305: }
309: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
310: {
311: PetscInt i;
312: QEP qep = (QEP)ctx;
316: nconv = 0;
317: for (i=0;i<PetscMin(nest,qep->ncv);i++) {
318: qep->eigr[i] = eigr[i];
319: qep->eigi[i] = eigi[i];
320: qep->errest[i] = errest[i];
321: if (0.0 < errest[i] && errest[i] < qep->tol) nconv++;
322: }
323: STBackTransform(eps->OP,nest,qep->eigr,qep->eigi);
324: QEPMonitor(qep,its,nconv,qep->eigr,qep->eigi,qep->errest,nest);
325: return(0);
326: }
330: PetscErrorCode QEPSetFromOptions_Linear(QEP qep)
331: {
333: PetscBool set,val;
334: PetscInt i;
335: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
336: ST st;
339: ctx->setfromoptionscalled = PETSC_TRUE;
340: PetscOptionsHead("QEP Linear Options");
341: PetscOptionsInt("-qep_linear_cform","Number of the companion form","QEPLinearSetCompanionForm",ctx->cform,&i,&set);
342: if (set) {
343: QEPLinearSetCompanionForm(qep,i);
344: }
345: PetscOptionsBool("-qep_linear_explicitmatrix","Use explicit matrix in linearization","QEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
346: if (set) {
347: QEPLinearSetExplicitMatrix(qep,val);
348: }
349: if (!ctx->explicitmatrix) {
350: /* use as default an ST with shell matrix and Jacobi */
351: EPSGetST(ctx->eps,&st);
352: STSetMatMode(st,ST_MATMODE_SHELL);
353: }
354: PetscOptionsTail();
355: return(0);
356: }
358: EXTERN_C_BEGIN
361: PetscErrorCode QEPLinearSetCompanionForm_Linear(QEP qep,PetscInt cform)
362: {
363: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
366: if (cform==PETSC_IGNORE) return(0);
367: if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
368: else {
369: if (cform!=1 && cform!=2) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
370: ctx->cform = cform;
371: }
372: return(0);
373: }
374: EXTERN_C_END
378: /*@
379: QEPLinearSetCompanionForm - Choose between the two companion forms available
380: for the linearization of the quadratic problem.
382: Logically Collective on QEP
384: Input Parameters:
385: + qep - quadratic eigenvalue solver
386: - cform - 1 or 2 (first or second companion form)
388: Options Database Key:
389: . -qep_linear_cform <int> - Choose the companion form
391: Level: advanced
393: .seealso: QEPLinearGetCompanionForm()
394: @*/
395: PetscErrorCode QEPLinearSetCompanionForm(QEP qep,PetscInt cform)
396: {
402: PetscTryMethod(qep,"QEPLinearSetCompanionForm_C",(QEP,PetscInt),(qep,cform));
403: return(0);
404: }
406: EXTERN_C_BEGIN
409: PetscErrorCode QEPLinearGetCompanionForm_Linear(QEP qep,PetscInt *cform)
410: {
411: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
414: *cform = ctx->cform;
415: return(0);
416: }
417: EXTERN_C_END
421: /*@
422: QEPLinearGetCompanionForm - Returns the number of the companion form that
423: will be used for the linearization of the quadratic problem.
425: Not Collective
427: Input Parameter:
428: . qep - quadratic eigenvalue solver
430: Output Parameter:
431: . cform - the companion form number (1 or 2)
433: Level: advanced
435: .seealso: QEPLinearSetCompanionForm()
436: @*/
437: PetscErrorCode QEPLinearGetCompanionForm(QEP qep,PetscInt *cform)
438: {
444: PetscTryMethod(qep,"QEPLinearGetCompanionForm_C",(QEP,PetscInt*),(qep,cform));
445: return(0);
446: }
448: EXTERN_C_BEGIN
451: PetscErrorCode QEPLinearSetExplicitMatrix_Linear(QEP qep,PetscBool explicitmatrix)
452: {
453: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
456: ctx->explicitmatrix = explicitmatrix;
457: return(0);
458: }
459: EXTERN_C_END
463: /*@
464: QEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the linearization
465: of the quadratic problem must be built explicitly.
467: Logically Collective on QEP
469: Input Parameters:
470: + qep - quadratic eigenvalue solver
471: - explicit - boolean flag indicating if the matrices are built explicitly
473: Options Database Key:
474: . -qep_linear_explicitmatrix <boolean> - Indicates the boolean flag
476: Level: advanced
478: .seealso: QEPLinearGetExplicitMatrix()
479: @*/
480: PetscErrorCode QEPLinearSetExplicitMatrix(QEP qep,PetscBool explicitmatrix)
481: {
487: PetscTryMethod(qep,"QEPLinearSetExplicitMatrix_C",(QEP,PetscBool),(qep,explicitmatrix));
488: return(0);
489: }
491: EXTERN_C_BEGIN
494: PetscErrorCode QEPLinearGetExplicitMatrix_Linear(QEP qep,PetscBool *explicitmatrix)
495: {
496: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
499: *explicitmatrix = ctx->explicitmatrix;
500: return(0);
501: }
502: EXTERN_C_END
506: /*@
507: QEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices A and B
508: for the linearization of the quadratic problem are built explicitly.
510: Not Collective
512: Input Parameter:
513: . qep - quadratic eigenvalue solver
515: Output Parameter:
516: . explicitmatrix - the mode flag
518: Level: advanced
520: .seealso: QEPLinearSetExplicitMatrix()
521: @*/
522: PetscErrorCode QEPLinearGetExplicitMatrix(QEP qep,PetscBool *explicitmatrix)
523: {
529: PetscTryMethod(qep,"QEPLinearGetExplicitMatrix_C",(QEP,PetscBool*),(qep,explicitmatrix));
530: return(0);
531: }
533: EXTERN_C_BEGIN
536: PetscErrorCode QEPLinearSetEPS_Linear(QEP qep,EPS eps)
537: {
539: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
542: PetscObjectReference((PetscObject)eps);
543: EPSDestroy(&ctx->eps);
544: ctx->eps = eps;
545: PetscLogObjectParent(qep,ctx->eps);
546: qep->setupcalled = 0;
547: return(0);
548: }
549: EXTERN_C_END
553: /*@
554: QEPLinearSetEPS - Associate an eigensolver object (EPS) to the
555: quadratic eigenvalue solver.
557: Collective on QEP
559: Input Parameters:
560: + qep - quadratic eigenvalue solver
561: - eps - the eigensolver object
563: Level: advanced
565: .seealso: QEPLinearGetEPS()
566: @*/
567: PetscErrorCode QEPLinearSetEPS(QEP qep,EPS eps)
568: {
575: PetscTryMethod(qep,"QEPLinearSetEPS_C",(QEP,EPS),(qep,eps));
576: return(0);
577: }
579: EXTERN_C_BEGIN
582: PetscErrorCode QEPLinearGetEPS_Linear(QEP qep,EPS *eps)
583: {
584: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
587: *eps = ctx->eps;
588: return(0);
589: }
590: EXTERN_C_END
594: /*@
595: QEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
596: to the quadratic eigenvalue solver.
598: Not Collective
600: Input Parameter:
601: . qep - quadratic eigenvalue solver
603: Output Parameter:
604: . eps - the eigensolver object
606: Level: advanced
608: .seealso: QEPLinearSetEPS()
609: @*/
610: PetscErrorCode QEPLinearGetEPS(QEP qep,EPS *eps)
611: {
617: PetscTryMethod(qep,"QEPLinearGetEPS_C",(QEP,EPS*),(qep,eps));
618: return(0);
619: }
623: PetscErrorCode QEPView_Linear(QEP qep,PetscViewer viewer)
624: {
626: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
629: PetscViewerASCIIPrintf(viewer," Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
630: PetscViewerASCIIPrintf(viewer," Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
631: PetscViewerASCIIPushTab(viewer);
632: EPSView(ctx->eps,viewer);
633: PetscViewerASCIIPopTab(viewer);
634: return(0);
635: }
639: PetscErrorCode QEPReset_Linear(QEP qep)
640: {
642: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
645: EPSReset(ctx->eps);
646: MatDestroy(&ctx->A);
647: MatDestroy(&ctx->B);
648: VecDestroy(&ctx->x1);
649: VecDestroy(&ctx->x2);
650: VecDestroy(&ctx->y1);
651: VecDestroy(&ctx->y2);
652: QEPDefaultFreeWork(qep);
653: QEPFreeSolution(qep);
654: return(0);
655: }
659: PetscErrorCode QEPDestroy_Linear(QEP qep)
660: {
662: QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
665: EPSDestroy(&ctx->eps);
666: PetscFree(qep->data);
667: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","",PETSC_NULL);
668: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","",PETSC_NULL);
669: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","",PETSC_NULL);
670: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","",PETSC_NULL);
671: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","",PETSC_NULL);
672: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","",PETSC_NULL);
673: return(0);
674: }
676: EXTERN_C_BEGIN
679: PetscErrorCode QEPCreate_Linear(QEP qep)
680: {
682: QEP_LINEAR *ctx;
685: PetscNewLog(qep,QEP_LINEAR,&ctx);
686: qep->data = (void *)ctx;
687: qep->ops->solve = QEPSolve_Linear;
688: qep->ops->setup = QEPSetUp_Linear;
689: qep->ops->setfromoptions = QEPSetFromOptions_Linear;
690: qep->ops->destroy = QEPDestroy_Linear;
691: qep->ops->reset = QEPReset_Linear;
692: qep->ops->view = QEPView_Linear;
693: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","QEPLinearSetCompanionForm_Linear",QEPLinearSetCompanionForm_Linear);
694: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","QEPLinearGetCompanionForm_Linear",QEPLinearGetCompanionForm_Linear);
695: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","QEPLinearSetEPS_Linear",QEPLinearSetEPS_Linear);
696: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","QEPLinearGetEPS_Linear",QEPLinearGetEPS_Linear);
697: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","QEPLinearSetExplicitMatrix_Linear",QEPLinearSetExplicitMatrix_Linear);
698: PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","QEPLinearGetExplicitMatrix_Linear",QEPLinearGetExplicitMatrix_Linear);
700: EPSCreate(((PetscObject)qep)->comm,&ctx->eps);
701: EPSSetOptionsPrefix(ctx->eps,((PetscObject)qep)->prefix);
702: EPSAppendOptionsPrefix(ctx->eps,"qep_");
703: STSetOptionsPrefix(ctx->eps->OP,((PetscObject)ctx->eps)->prefix);
704: PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)qep,1);
705: PetscLogObjectParent(qep,ctx->eps);
706: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
707: EPSSetIP(ctx->eps,qep->ip);
708: EPSMonitorSet(ctx->eps,EPSMonitor_Linear,qep,PETSC_NULL);
709: ctx->explicitmatrix = PETSC_FALSE;
710: ctx->cform = 1;
711: ctx->A = PETSC_NULL;
712: ctx->B = PETSC_NULL;
713: ctx->x1 = PETSC_NULL;
714: ctx->x2 = PETSC_NULL;
715: ctx->y1 = PETSC_NULL;
716: ctx->y2 = PETSC_NULL;
717: ctx->setfromoptionscalled = PETSC_FALSE;
718: return(0);
719: }
720: EXTERN_C_END