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-2013, Universitat Politecnica de Valencia, Spain
  9:    This file is part of SLEPc.
 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:   PetscScalar    sigma;
 39:   /* function tables */
 40:   PetscErrorCode (*fcreate[][2])(MPI_Comm,QEP_LINEAR*,Mat*) = {
 41:     { MatCreateExplicit_Linear_N1A, MatCreateExplicit_Linear_N1B },   /* N1 */
 42:     { MatCreateExplicit_Linear_N2A, MatCreateExplicit_Linear_N2B },   /* N2 */
 43:     { MatCreateExplicit_Linear_S1A, MatCreateExplicit_Linear_S1B },   /* S1 */
 44:     { MatCreateExplicit_Linear_S2A, MatCreateExplicit_Linear_S2B },   /* S2 */
 45:     { MatCreateExplicit_Linear_H1A, MatCreateExplicit_Linear_H1B },   /* H1 */
 46:     { MatCreateExplicit_Linear_H2A, MatCreateExplicit_Linear_H2B }    /* H2 */
 47:   };
 48:   PetscErrorCode (*fmult[][2])(Mat,Vec,Vec) = {
 49:     { MatMult_Linear_N1A, MatMult_Linear_N1B },
 50:     { MatMult_Linear_N2A, MatMult_Linear_N2B },
 51:     { MatMult_Linear_S1A, MatMult_Linear_S1B },
 52:     { MatMult_Linear_S2A, MatMult_Linear_S2B },
 53:     { MatMult_Linear_H1A, MatMult_Linear_H1B },
 54:     { MatMult_Linear_H2A, MatMult_Linear_H2B }
 55:   };
 56:   PetscErrorCode (*fgetdiagonal[][2])(Mat,Vec) = {
 57:     { MatGetDiagonal_Linear_N1A, MatGetDiagonal_Linear_N1B },
 58:     { MatGetDiagonal_Linear_N2A, MatGetDiagonal_Linear_N2B },
 59:     { MatGetDiagonal_Linear_S1A, MatGetDiagonal_Linear_S1B },
 60:     { MatGetDiagonal_Linear_S2A, MatGetDiagonal_Linear_S2B },
 61:     { MatGetDiagonal_Linear_H1A, MatGetDiagonal_Linear_H1B },
 62:     { MatGetDiagonal_Linear_H2A, MatGetDiagonal_Linear_H2B }
 63:   };
 66:   if (!ctx->cform) ctx->cform = 1;
 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:   MatDestroy(&ctx->A);
 74:   MatDestroy(&ctx->B);
 75:   VecDestroy(&ctx->x1);
 76:   VecDestroy(&ctx->x2);
 77:   VecDestroy(&ctx->y1);
 78:   VecDestroy(&ctx->y2);
 80:   switch (qep->problem_type) {
 81:     case QEP_GENERAL:    i = 0; break;
 82:     case QEP_HERMITIAN:  i = 2; break;
 83:     case QEP_GYROSCOPIC: i = 4; break;
 84:     default: SETERRQ(PetscObjectComm((PetscObject)qep),1,"Wrong value of qep->problem_type");
 85:   }
 86:   i += ctx->cform-1;
 88:   if (ctx->explicitmatrix) {
 89:     ctx->x1 = ctx->x2 = ctx->y1 = ctx->y2 = NULL;
 90:     (*fcreate[i][0])(PetscObjectComm((PetscObject)qep),ctx,&ctx->A);
 91:     (*fcreate[i][1])(PetscObjectComm((PetscObject)qep),ctx,&ctx->B);
 92:   } else {
 93:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&ctx->x1);
 94:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&ctx->x2);
 95:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&ctx->y1);
 96:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&ctx->y2);
 97:     PetscLogObjectParent(qep,ctx->x1);
 98:     PetscLogObjectParent(qep,ctx->x2);
 99:     PetscLogObjectParent(qep,ctx->y1);
100:     PetscLogObjectParent(qep,ctx->y2);
101:     MatCreateShell(PetscObjectComm((PetscObject)qep),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(PetscObjectComm((PetscObject)qep),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:   }
108:   PetscLogObjectParent(qep,ctx->A);
109:   PetscLogObjectParent(qep,ctx->B);
111:   if (!ctx->eps) { QEPLinearGetEPS(qep,&ctx->eps); }
112:   EPSSetOperators(ctx->eps,ctx->A,ctx->B);
113:   if (qep->problem_type==QEP_HERMITIAN) {
114:     EPSSetProblemType(ctx->eps,EPS_GHIEP);
115:   } else {
116:     EPSSetProblemType(ctx->eps,EPS_GNHEP);
117:   }
118:   switch (qep->which) {
119:       case QEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
120:       case QEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
121:       case QEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
122:       case QEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
123:       case QEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
124:       case QEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
125:       default: SETERRQ(PetscObjectComm((PetscObject)qep),1,"Wrong value of which");
126:   }
127:   EPSSetWhichEigenpairs(ctx->eps,which);
128:   EPSSetLeftVectorsWanted(ctx->eps,qep->leftvecs);
129:   EPSSetDimensions(ctx->eps,qep->nev,qep->ncv,qep->mpd);
130:   EPSSetTolerances(ctx->eps,qep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:qep->tol/10.0,qep->max_it);
131:   /* Transfer the trackall option from qep to eps */
132:   QEPGetTrackAll(qep,&trackall);
133:   EPSSetTrackAll(ctx->eps,trackall);
134:   if (ctx->setfromoptionscalled) {
135:     EPSSetFromOptions(ctx->eps);
136:     ctx->setfromoptionscalled = PETSC_FALSE;
137:   }
138:   /* temporary change of target */
139:   if (qep->sfactor!=1.0) {
140:     EPSGetTarget(ctx->eps,&sigma);
141:     EPSSetTarget(ctx->eps,sigma/qep->sfactor);
142:   }
143:   EPSSetUp(ctx->eps);
144:   EPSGetDimensions(ctx->eps,NULL,&qep->ncv,&qep->mpd);
145:   EPSGetTolerances(ctx->eps,NULL,&qep->max_it);
146:   if (qep->nini>0 || qep->ninil>0) { PetscInfo(qep,"Ignoring initial vectors\n"); }
147:   QEPAllocateSolution(qep);
148:   return(0);
149: }
153: /*
154:    QEPLinearSelect_Norm - Auxiliary routine that copies the solution of the
155:    linear eigenproblem to the QEP object. The eigenvector of the generalized
156:    problem is supposed to be
157:                                z = [  x  ]
158:                                    [ l*x ]
159:    The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
160:    computed residual norm.
161:    Finally, x is normalized so that ||x||_2 = 1.
162: */
163: static PetscErrorCode QEPLinearSelect_Norm(QEP qep,EPS eps)
164: {
166:   PetscInt       i;
167:   PetscScalar    *px;
168:   PetscReal      rn1,rn2;
169:   Vec            xr,xi,wr,wi;
170:   Mat            A;
171: #if !defined(PETSC_USE_COMPLEX)
172:   PetscScalar    *py;
173: #endif
176:   EPSGetOperators(eps,&A,NULL);
177:   MatGetVecs(A,&xr,NULL);
178:   VecDuplicate(xr,&xi);
179:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&wr);
180:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&wi);
181:   for (i=0;i<qep->nconv;i++) {
182:     EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
183:     qep->eigr[i] *= qep->sfactor;
184:     qep->eigi[i] *= qep->sfactor;
185: #if !defined(PETSC_USE_COMPLEX)
186:     if (qep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
187:       VecGetArray(xr,&px);
188:       VecGetArray(xi,&py);
189:       VecPlaceArray(wr,px);
190:       VecPlaceArray(wi,py);
191:       SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
192:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
193:       VecCopy(wr,qep->V[i]);
194:       VecCopy(wi,qep->V[i+1]);
195:       VecResetArray(wr);
196:       VecResetArray(wi);
197:       VecPlaceArray(wr,px+qep->nloc);
198:       VecPlaceArray(wi,py+qep->nloc);
199:       SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
200:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
201:       if (rn1>rn2) {
202:         VecCopy(wr,qep->V[i]);
203:         VecCopy(wi,qep->V[i+1]);
204:       }
205:       VecResetArray(wr);
206:       VecResetArray(wi);
207:       VecRestoreArray(xr,&px);
208:       VecRestoreArray(xi,&py);
209:     } else if (qep->eigi[i]==0.0)   /* real eigenvalue */
210: #endif
211:     {
212:       VecGetArray(xr,&px);
213:       VecPlaceArray(wr,px);
214:       SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
215:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,NULL,&rn1);
216:       VecCopy(wr,qep->V[i]);
217:       VecResetArray(wr);
218:       VecPlaceArray(wr,px+qep->nloc);
219:       SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
220:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,NULL,&rn2);
221:       if (rn1>rn2) {
222:         VecCopy(wr,qep->V[i]);
223:       }
224:       VecResetArray(wr);
225:       VecRestoreArray(xr,&px);
226:     }
227:   }
228:   VecDestroy(&wr);
229:   VecDestroy(&wi);
230:   VecDestroy(&xr);
231:   VecDestroy(&xi);
232:   return(0);
233: }
237: /*
238:    QEPLinearSelect_Simple - Auxiliary routine that copies the solution of the
239:    linear eigenproblem to the QEP object. The eigenvector of the generalized
240:    problem is supposed to be
241:                                z = [  x  ]
242:                                    [ l*x ]
243:    If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
244:    Finally, x is normalized so that ||x||_2 = 1.
245: */
246: static PetscErrorCode QEPLinearSelect_Simple(QEP qep,EPS eps)
247: {
249:   PetscInt       i,offset;
250:   PetscScalar    *px;
251:   Vec            xr,xi,w;
252:   Mat            A;
255:   EPSGetOperators(eps,&A,NULL);
256:   MatGetVecs(A,&xr,NULL);
257:   VecDuplicate(xr,&xi);
258:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)qep),1,qep->nloc,qep->n,NULL,&w);
259:   for (i=0;i<qep->nconv;i++) {
260:     EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
261:     qep->eigr[i] *= qep->sfactor;
262:     qep->eigi[i] *= qep->sfactor;
263:     if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) offset = qep->nloc;
264:     else offset = 0;
265: #if !defined(PETSC_USE_COMPLEX)
266:     if (qep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
267:       VecGetArray(xr,&px);
268:       VecPlaceArray(w,px+offset);
269:       VecCopy(w,qep->V[i]);
270:       VecResetArray(w);
271:       VecRestoreArray(xr,&px);
272:       VecGetArray(xi,&px);
273:       VecPlaceArray(w,px+offset);
274:       VecCopy(w,qep->V[i+1]);
275:       VecResetArray(w);
276:       VecRestoreArray(xi,&px);
277:       SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,NULL);
278:     } else if (qep->eigi[i]==0.0)   /* real eigenvalue */
279: #endif
280:     {
281:       VecGetArray(xr,&px);
282:       VecPlaceArray(w,px+offset);
283:       VecCopy(w,qep->V[i]);
284:       VecResetArray(w);
285:       VecRestoreArray(xr,&px);
286:       SlepcVecNormalize(qep->V[i],NULL,PETSC_FALSE,NULL);
287:     }
288:   }
289:   VecDestroy(&w);
290:   VecDestroy(&xr);
291:   VecDestroy(&xi);
292:   return(0);
293: }
297: PetscErrorCode QEPSolve_Linear(QEP qep)
298: {
300:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
301:   PetscBool      flg=PETSC_FALSE;
302:   PetscScalar    sigma;
305:   EPSSolve(ctx->eps);
306:   EPSGetConverged(ctx->eps,&qep->nconv);
307:   EPSGetIterationNumber(ctx->eps,&qep->its);
308:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&qep->reason);
309:   EPSGetOperationCounters(ctx->eps,&qep->matvecs,NULL,&qep->linits);
310:   /* restore target */
311:   EPSGetTarget(ctx->eps,&sigma);
312:   EPSSetTarget(ctx->eps,sigma*qep->sfactor);
314:   qep->matvecs *= 2;  /* convention: count one matvec for each non-trivial block in A */
315:   PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_linear_select_simple",&flg,NULL);
316:   if (flg) {
317:     QEPLinearSelect_Simple(qep,ctx->eps);
318:   } else {
319:     QEPLinearSelect_Norm(qep,ctx->eps);
320:   }
321:   return(0);
322: }
326: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
327: {
328:   PetscInt       i;
329:   QEP            qep = (QEP)ctx;
333:   nconv = 0;
334:   for (i=0;i<PetscMin(nest,qep->ncv);i++) {
335:     qep->eigr[i] = eigr[i];
336:     qep->eigi[i] = eigi[i];
337:     qep->errest[i] = errest[i];
338:     if (0.0 < errest[i] && errest[i] < qep->tol) nconv++;
339:   }
340:   STBackTransform(eps->st,nest,qep->eigr,qep->eigi);
341:   QEPMonitor(qep,its,nconv,qep->eigr,qep->eigi,qep->errest,nest);
342:   return(0);
343: }
347: PetscErrorCode QEPSetFromOptions_Linear(QEP qep)
348: {
350:   PetscBool      set,val;
351:   PetscInt       i;
352:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
353:   ST             st;
356:   ctx->setfromoptionscalled = PETSC_TRUE;
357:   PetscOptionsHead("QEP Linear Options");
358:   PetscOptionsInt("-qep_linear_cform","Number of the companion form","QEPLinearSetCompanionForm",ctx->cform,&i,&set);
359:   if (set) {
360:     QEPLinearSetCompanionForm(qep,i);
361:   }
362:   PetscOptionsBool("-qep_linear_explicitmatrix","Use explicit matrix in linearization","QEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
363:   if (set) {
364:     QEPLinearSetExplicitMatrix(qep,val);
365:   }
366:   if (!ctx->explicitmatrix) {
367:     /* use as default an ST with shell matrix and Jacobi */
368:     if (!ctx->eps) { QEPLinearGetEPS(qep,&ctx->eps); }
369:     EPSGetST(ctx->eps,&st);
370:     STSetMatMode(st,ST_MATMODE_SHELL);
371:   }
372:   PetscOptionsTail();
373:   return(0);
374: }
378: static PetscErrorCode QEPLinearSetCompanionForm_Linear(QEP qep,PetscInt cform)
379: {
380:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
383:   if (!cform) return(0);
384:   if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
385:   else {
386:     if (cform!=1 && cform!=2) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
387:     ctx->cform = cform;
388:   }
389:   return(0);
390: }
394: /*@
395:    QEPLinearSetCompanionForm - Choose between the two companion forms available
396:    for the linearization of the quadratic problem.
398:    Logically Collective on QEP
400:    Input Parameters:
401: +  qep   - quadratic eigenvalue solver
402: -  cform - 1 or 2 (first or second companion form)
404:    Options Database Key:
405: .  -qep_linear_cform <int> - Choose the companion form
407:    Level: advanced
409: .seealso: QEPLinearGetCompanionForm()
410: @*/
411: PetscErrorCode QEPLinearSetCompanionForm(QEP qep,PetscInt cform)
412: {
418:   PetscTryMethod(qep,"QEPLinearSetCompanionForm_C",(QEP,PetscInt),(qep,cform));
419:   return(0);
420: }
424: static PetscErrorCode QEPLinearGetCompanionForm_Linear(QEP qep,PetscInt *cform)
425: {
426:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
429:   *cform = ctx->cform;
430:   return(0);
431: }
435: /*@
436:    QEPLinearGetCompanionForm - Returns the number of the companion form that
437:    will be used for the linearization of the quadratic problem.
439:    Not Collective
441:    Input Parameter:
442: .  qep  - quadratic eigenvalue solver
444:    Output Parameter:
445: .  cform - the companion form number (1 or 2)
447:    Level: advanced
449: .seealso: QEPLinearSetCompanionForm()
450: @*/
451: PetscErrorCode QEPLinearGetCompanionForm(QEP qep,PetscInt *cform)
452: {
458:   PetscTryMethod(qep,"QEPLinearGetCompanionForm_C",(QEP,PetscInt*),(qep,cform));
459:   return(0);
460: }
464: static PetscErrorCode QEPLinearSetExplicitMatrix_Linear(QEP qep,PetscBool explicitmatrix)
465: {
466:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
469:   ctx->explicitmatrix = explicitmatrix;
470:   return(0);
471: }
475: /*@
476:    QEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the linearization
477:    of the quadratic problem must be built explicitly.
479:    Logically Collective on QEP
481:    Input Parameters:
482: +  qep      - quadratic eigenvalue solver
483: -  explicit - boolean flag indicating if the matrices are built explicitly
485:    Options Database Key:
486: .  -qep_linear_explicitmatrix <boolean> - Indicates the boolean flag
488:    Level: advanced
490: .seealso: QEPLinearGetExplicitMatrix()
491: @*/
492: PetscErrorCode QEPLinearSetExplicitMatrix(QEP qep,PetscBool explicitmatrix)
493: {
499:   PetscTryMethod(qep,"QEPLinearSetExplicitMatrix_C",(QEP,PetscBool),(qep,explicitmatrix));
500:   return(0);
501: }
505: static PetscErrorCode QEPLinearGetExplicitMatrix_Linear(QEP qep,PetscBool *explicitmatrix)
506: {
507:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;
510:   *explicitmatrix = ctx->explicitmatrix;
511:   return(0);
512: }
516: /*@
517:    QEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices A and B
518:    for the linearization of the quadratic problem are built explicitly.
520:    Not Collective
522:    Input Parameter:
523: .  qep  - quadratic eigenvalue solver
525:    Output Parameter:
526: .  explicitmatrix - the mode flag
528:    Level: advanced
530: .seealso: QEPLinearSetExplicitMatrix()
531: @*/
532: PetscErrorCode QEPLinearGetExplicitMatrix(QEP qep,PetscBool *explicitmatrix)
533: {
539:   PetscTryMethod(qep,"QEPLinearGetExplicitMatrix_C",(QEP,PetscBool*),(qep,explicitmatrix));
540:   return(0);
541: }
545: static PetscErrorCode QEPLinearSetEPS_Linear(QEP qep,EPS eps)
546: {
548:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
551:   PetscObjectReference((PetscObject)eps);
552:   EPSDestroy(&ctx->eps);
553:   ctx->eps = eps;
554:   PetscLogObjectParent(qep,ctx->eps);
555:   qep->setupcalled = 0;
556:   return(0);
557: }
561: /*@
562:    QEPLinearSetEPS - Associate an eigensolver object (EPS) to the
563:    quadratic eigenvalue solver.
565:    Collective on QEP
567:    Input Parameters:
568: +  qep - quadratic eigenvalue solver
569: -  eps - the eigensolver object
571:    Level: advanced
573: .seealso: QEPLinearGetEPS()
574: @*/
575: PetscErrorCode QEPLinearSetEPS(QEP qep,EPS eps)
576: {
583:   PetscTryMethod(qep,"QEPLinearSetEPS_C",(QEP,EPS),(qep,eps));
584:   return(0);
585: }
589: static PetscErrorCode QEPLinearGetEPS_Linear(QEP qep,EPS *eps)
590: {
592:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
595:   if (!ctx->eps) {
596:     EPSCreate(PetscObjectComm((PetscObject)qep),&ctx->eps);
597:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)qep)->prefix);
598:     EPSAppendOptionsPrefix(ctx->eps,"qep_");
599:     STSetOptionsPrefix(ctx->eps->st,((PetscObject)ctx->eps)->prefix);
600:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)qep,1);
601:     PetscLogObjectParent(qep,ctx->eps);
602:     if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
603:     EPSSetIP(ctx->eps,qep->ip);
604:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,qep,NULL);
605:   }
606:   *eps = ctx->eps;
607:   return(0);
608: }
612: /*@
613:    QEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
614:    to the quadratic eigenvalue solver.
616:    Not Collective
618:    Input Parameter:
619: .  qep - quadratic eigenvalue solver
621:    Output Parameter:
622: .  eps - the eigensolver object
624:    Level: advanced
626: .seealso: QEPLinearSetEPS()
627: @*/
628: PetscErrorCode QEPLinearGetEPS(QEP qep,EPS *eps)
629: {
635:   PetscTryMethod(qep,"QEPLinearGetEPS_C",(QEP,EPS*),(qep,eps));
636:   return(0);
637: }
641: PetscErrorCode QEPView_Linear(QEP qep,PetscViewer viewer)
642: {
644:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
647:   if (!ctx->eps) { QEPLinearGetEPS(qep,&ctx->eps); }
648:   PetscViewerASCIIPrintf(viewer,"  Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
649:   PetscViewerASCIIPrintf(viewer,"  Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
650:   PetscViewerASCIIPushTab(viewer);
651:   EPSView(ctx->eps,viewer);
652:   PetscViewerASCIIPopTab(viewer);
653:   return(0);
654: }
658: PetscErrorCode QEPReset_Linear(QEP qep)
659: {
661:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
664:   if (!ctx->eps) { EPSReset(ctx->eps); }
665:   MatDestroy(&ctx->A);
666:   MatDestroy(&ctx->B);
667:   VecDestroy(&ctx->x1);
668:   VecDestroy(&ctx->x2);
669:   VecDestroy(&ctx->y1);
670:   VecDestroy(&ctx->y2);
671:   QEPReset_Default(qep);
672:   return(0);
673: }
677: PetscErrorCode QEPDestroy_Linear(QEP qep)
678: {
680:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
683:   EPSDestroy(&ctx->eps);
684:   PetscFree(qep->data);
685:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetCompanionForm_C",NULL);
686:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetCompanionForm_C",NULL);
687:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetEPS_C",NULL);
688:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetEPS_C",NULL);
689:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetExplicitMatrix_C",NULL);
690:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetExplicitMatrix_C",NULL);
691:   return(0);
692: }
696: PETSC_EXTERN PetscErrorCode QEPCreate_Linear(QEP qep)
697: {
699:   QEP_LINEAR     *ctx;
702:   PetscNewLog(qep,QEP_LINEAR,&ctx);
703:   qep->data                      = (void*)ctx;
704:   qep->ops->solve                = QEPSolve_Linear;
705:   qep->ops->setup                = QEPSetUp_Linear;
706:   qep->ops->setfromoptions       = QEPSetFromOptions_Linear;
707:   qep->ops->destroy              = QEPDestroy_Linear;
708:   qep->ops->reset                = QEPReset_Linear;
709:   qep->ops->view                 = QEPView_Linear;
710:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetCompanionForm_C",QEPLinearSetCompanionForm_Linear);
711:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetCompanionForm_C",QEPLinearGetCompanionForm_Linear);
712:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetEPS_C",QEPLinearSetEPS_Linear);
713:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetEPS_C",QEPLinearGetEPS_Linear);
714:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearSetExplicitMatrix_C",QEPLinearSetExplicitMatrix_Linear);
715:   PetscObjectComposeFunction((PetscObject)qep,"QEPLinearGetExplicitMatrix_C",QEPLinearGetExplicitMatrix_Linear);
716:   return(0);
717: }