Actual source code: linear.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    Straightforward linearization for quadratic eigenproblems.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc-private/pepimpl.h>         /*I "slepcpep.h" I*/
 25: #include <slepcvec.h>
 26:  #include linearp.h

 30: PetscErrorCode PEPSetUp_Linear(PEP pep)
 31: {
 33:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
 34:   PetscInt       i=0;
 35:   EPSWhich       which;
 36:   PetscBool      trackall,istrivial,flg;
 37:   PetscScalar    sigma;
 38:   /* function tables */
 39:   PetscErrorCode (*fcreate[][2])(MPI_Comm,PEP_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 (!ctx->cform) ctx->cform = 1;
 66:   if (!pep->which) pep->which = PEP_LARGEST_MAGNITUDE;
 67:   if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
 68:   if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver only available for quadratic problems");
 69:   if (pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Diagonal scaling not allowed in PEP linear solver");
 70:   STGetTransform(pep->st,&flg);
 71:   if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"ST transformation flag not allowed for PEP linear solver");

 73:   /* compute scale factor if no set by user */
 74:   PEPComputeScaleFactor(pep);
 75:  
 76:   STGetOperators(pep->st,0,&ctx->K);
 77:   STGetOperators(pep->st,1,&ctx->C);
 78:   STGetOperators(pep->st,2,&ctx->M);
 79:   ctx->sfactor = pep->sfactor;

 81:   MatDestroy(&ctx->A);
 82:   MatDestroy(&ctx->B);
 83:   VecDestroy(&ctx->x1);
 84:   VecDestroy(&ctx->x2);
 85:   VecDestroy(&ctx->y1);
 86:   VecDestroy(&ctx->y2);

 88:   switch (pep->problem_type) {
 89:     case PEP_GENERAL:    i = 0; break;
 90:     case PEP_HERMITIAN:  i = 2; break;
 91:     case PEP_GYROSCOPIC: i = 4; break;
 92:     default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->problem_type");
 93:   }
 94:   i += ctx->cform-1;

 96:   if (ctx->explicitmatrix) {
 97:     ctx->x1 = ctx->x2 = ctx->y1 = ctx->y2 = NULL;
 98:     (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
 99:     (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);
100:   } else {
101:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->x1);
102:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->x2);
103:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->y1);
104:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->y2);
105:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->x1);
106:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->x2);
107:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->y1);
108:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->y2);
109:     MatCreateShell(PetscObjectComm((PetscObject)pep),2*pep->nloc,2*pep->nloc,2*pep->n,2*pep->n,ctx,&ctx->A);
110:     MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))fmult[i][0]);
111:     MatShellSetOperation(ctx->A,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][0]);
112:     MatCreateShell(PetscObjectComm((PetscObject)pep),2*pep->nloc,2*pep->nloc,2*pep->n,2*pep->n,ctx,&ctx->B);
113:     MatShellSetOperation(ctx->B,MATOP_MULT,(void(*)(void))fmult[i][1]);
114:     MatShellSetOperation(ctx->B,MATOP_GET_DIAGONAL,(void(*)(void))fgetdiagonal[i][1]);
115:   }
116:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
117:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);

119:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
120:   EPSSetOperators(ctx->eps,ctx->A,ctx->B);
121:   if (pep->problem_type==PEP_HERMITIAN) {
122:     EPSSetProblemType(ctx->eps,EPS_GHIEP);
123:   } else {
124:     EPSSetProblemType(ctx->eps,EPS_GNHEP);
125:   }
126:   switch (pep->which) {
127:       case PEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
128:       case PEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
129:       case PEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
130:       case PEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
131:       case PEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
132:       case PEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
133:       case PEP_TARGET_MAGNITUDE:   which = EPS_TARGET_MAGNITUDE; break;
134:       case PEP_TARGET_REAL:        which = EPS_TARGET_REAL; break;
135:       case PEP_TARGET_IMAGINARY:   which = EPS_TARGET_IMAGINARY; break;
136:       default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of which");
137:   }
138:   EPSSetWhichEigenpairs(ctx->eps,which);
139:   EPSSetDimensions(ctx->eps,pep->nev,pep->ncv?pep->ncv:PETSC_DEFAULT,pep->mpd?pep->mpd:PETSC_DEFAULT);
140:   EPSSetTolerances(ctx->eps,pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:pep->tol/10.0,pep->max_it?pep->max_it:PETSC_DEFAULT);
141:   RGIsTrivial(pep->rg,&istrivial);
142:   if (!istrivial) { EPSSetRG(ctx->eps,pep->rg); }
143:   /* Transfer the trackall option from pep to eps */
144:   PEPGetTrackAll(pep,&trackall);
145:   EPSSetTrackAll(ctx->eps,trackall);

147:   /* temporary change of target */
148:   if (pep->sfactor!=1.0) {
149:     EPSGetTarget(ctx->eps,&sigma);
150:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
151:   }
152:   EPSSetUp(ctx->eps);
153:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
154:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
155:   if (pep->nini>0) { PetscInfo(pep,"Ignoring initial vectors\n"); }
156:   PEPAllocateSolution(pep,0);
157:   return(0);
158: }

162: /*
163:    PEPLinearExtract_Residual - Auxiliary routine that copies the solution of the
164:    linear eigenproblem to the PEP object. The eigenvector of the generalized
165:    problem is supposed to be
166:                                z = [  x  ]
167:                                    [ l*x ]
168:    The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
169:    computed residual norm.
170:    Finally, x is normalized so that ||x||_2 = 1.
171: */
172: static PetscErrorCode PEPLinearExtract_Residual(PEP pep,EPS eps)
173: {
175:   PetscInt       i;
176:   PetscScalar    *px;
177:   PetscReal      rn1,rn2;
178:   Vec            xr,xi,wr,wi;
179:   Mat            A;
180: #if !defined(PETSC_USE_COMPLEX)
181:   PetscScalar    *py;
182: #endif

185:   EPSGetOperators(eps,&A,NULL);
186:   MatGetVecs(A,&xr,NULL);
187:   VecDuplicate(xr,&xi);
188:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&wr);
189:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&wi);
190:   for (i=0;i<pep->nconv;i++) {
191:     EPSGetEigenpair(eps,i,&pep->eigr[i],&pep->eigi[i],xr,xi);
192:     pep->eigr[i] *= pep->sfactor;
193:     pep->eigi[i] *= pep->sfactor;
194: #if !defined(PETSC_USE_COMPLEX)
195:     if (pep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
196:       VecGetArray(xr,&px);
197:       VecGetArray(xi,&py);
198:       VecPlaceArray(wr,px);
199:       VecPlaceArray(wi,py);
200:       SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
201:       PEPComputeResidualNorm_Private(pep,pep->eigr[i],pep->eigi[i],wr,wi,&rn1);
202:       BVInsertVec(pep->V,i,wr);
203:       BVInsertVec(pep->V,i+1,wi);
204:       VecResetArray(wr);
205:       VecResetArray(wi);
206:       VecPlaceArray(wr,px+pep->nloc);
207:       VecPlaceArray(wi,py+pep->nloc);
208:       SlepcVecNormalize(wr,wi,PETSC_TRUE,NULL);
209:       PEPComputeResidualNorm_Private(pep,pep->eigr[i],pep->eigi[i],wr,wi,&rn2);
210:       if (rn1>rn2) {
211:         BVInsertVec(pep->V,i,wr);
212:         BVInsertVec(pep->V,i+1,wi);
213:       }
214:       VecResetArray(wr);
215:       VecResetArray(wi);
216:       VecRestoreArray(xr,&px);
217:       VecRestoreArray(xi,&py);
218:     } else if (pep->eigi[i]==0.0)   /* real eigenvalue */
219: #endif
220:     {
221:       VecGetArray(xr,&px);
222:       VecPlaceArray(wr,px);
223:       SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
224:       PEPComputeResidualNorm_Private(pep,pep->eigr[i],pep->eigi[i],wr,NULL,&rn1);
225:       BVInsertVec(pep->V,i,wr);
226:       VecResetArray(wr);
227:       VecPlaceArray(wr,px+pep->nloc);
228:       SlepcVecNormalize(wr,NULL,PETSC_FALSE,NULL);
229:       PEPComputeResidualNorm_Private(pep,pep->eigr[i],pep->eigi[i],wr,NULL,&rn2);
230:       if (rn1>rn2) {
231:         BVInsertVec(pep->V,i,wr);
232:       }
233:       VecResetArray(wr);
234:       VecRestoreArray(xr,&px);
235:     }
236:   }
237:   VecDestroy(&wr);
238:   VecDestroy(&wi);
239:   VecDestroy(&xr);
240:   VecDestroy(&xi);
241:   return(0);
242: }

246: /*
247:    PEPLinearExtract_Norm - Auxiliary routine that copies the solution of the
248:    linear eigenproblem to the PEP object. The eigenvector of the generalized
249:    problem is supposed to be
250:                                z = [  x  ]
251:                                    [ l*x ]
252:    If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
253:    Finally, x is normalized so that ||x||_2 = 1.
254: */
255: static PetscErrorCode PEPLinearExtract_Norm(PEP pep,EPS eps)
256: {
258:   PetscInt       i,offset;
259:   PetscScalar    *px;
260:   Vec            xr,xi,w,vi;
261: #if !defined(PETSC_USE_COMPLEX)
262:   Vec            vi1;
263: #endif
264:   Mat            A;

267:   EPSGetOperators(eps,&A,NULL);
268:   MatGetVecs(A,&xr,NULL);
269:   VecDuplicate(xr,&xi);
270:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&w);
271:   for (i=0;i<pep->nconv;i++) {
272:     EPSGetEigenpair(eps,i,&pep->eigr[i],&pep->eigi[i],xr,xi);
273:     pep->eigr[i] *= pep->sfactor;
274:     pep->eigi[i] *= pep->sfactor;
275:     if (SlepcAbsEigenvalue(pep->eigr[i],pep->eigi[i])>1.0) offset = pep->nloc;
276:     else offset = 0;
277: #if !defined(PETSC_USE_COMPLEX)
278:     if (pep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
279:       VecGetArray(xr,&px);
280:       VecPlaceArray(w,px+offset);
281:       BVInsertVec(pep->V,i,w);
282:       VecResetArray(w);
283:       VecRestoreArray(xr,&px);
284:       VecGetArray(xi,&px);
285:       VecPlaceArray(w,px+offset);
286:       BVInsertVec(pep->V,i+1,w);
287:       VecResetArray(w);
288:       VecRestoreArray(xi,&px);
289:       BVGetColumn(pep->V,i,&vi);
290:       BVGetColumn(pep->V,i+1,&vi1);
291:       SlepcVecNormalize(vi,vi1,PETSC_TRUE,NULL);
292:       BVRestoreColumn(pep->V,i,&vi);
293:       BVRestoreColumn(pep->V,i+1,&vi1);
294:     } else if (pep->eigi[i]==0.0)   /* real eigenvalue */
295: #endif
296:     {
297:       VecGetArray(xr,&px);
298:       VecPlaceArray(w,px+offset);
299:       BVInsertVec(pep->V,i,w);
300:       VecResetArray(w);
301:       VecRestoreArray(xr,&px);
302:       BVGetColumn(pep->V,i,&vi);
303:       SlepcVecNormalize(vi,NULL,PETSC_FALSE,NULL);
304:       BVRestoreColumn(pep->V,i,&vi);
305:     }
306:   }
307:   VecDestroy(&w);
308:   VecDestroy(&xr);
309:   VecDestroy(&xi);
310:   return(0);
311: }

315: PetscErrorCode PEPSolve_Linear(PEP pep)
316: {
318:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
319:   PetscScalar    sigma;

322:   EPSSolve(ctx->eps);
323:   EPSGetConverged(ctx->eps,&pep->nconv);
324:   EPSGetIterationNumber(ctx->eps,&pep->its);
325:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);
326:   /* restore target */
327:   EPSGetTarget(ctx->eps,&sigma);
328:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

330:   switch (pep->extract) {
331:   case PEP_EXTRACT_NORM:
332:     PEPLinearExtract_Norm(pep,ctx->eps);
333:     break;
334:   case PEP_EXTRACT_RESIDUAL:
335:     PEPLinearExtract_Residual(pep,ctx->eps);
336:     break;
337:   default:
338:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
339:   }
340:   return(0);
341: }

345: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
346: {
347:   PetscInt       i;
348:   PEP            pep = (PEP)ctx;
349:   ST             st;

353:   for (i=0;i<PetscMin(nest,pep->ncv);i++) {
354:     pep->eigr[i] = eigr[i];
355:     pep->eigi[i] = eigi[i];
356:     pep->errest[i] = errest[i];
357:   }
358:   EPSGetST(eps,&st);
359:   STBackTransform(st,nest,pep->eigr,pep->eigi);
360:   PEPMonitor(pep,its,nconv,pep->eigr,pep->eigi,pep->errest,nest);
361:   return(0);
362: }

366: PetscErrorCode PEPSetFromOptions_Linear(PEP pep)
367: {
369:   PetscBool      set,val;
370:   PetscInt       i;
371:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
372:   ST             st;

375:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
376:   EPSSetFromOptions(ctx->eps);
377:   PetscOptionsHead("PEP Linear Options");
378:   PetscOptionsInt("-pep_linear_cform","Number of the companion form","PEPLinearSetCompanionForm",ctx->cform,&i,&set);
379:   if (set) {
380:     PEPLinearSetCompanionForm(pep,i);
381:   }
382:   PetscOptionsBool("-pep_linear_explicitmatrix","Use explicit matrix in linearization","PEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
383:   if (set) {
384:     PEPLinearSetExplicitMatrix(pep,val);
385:   }
386:   if (!ctx->explicitmatrix) {
387:     /* use as default an ST with shell matrix and Jacobi */
388:     if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
389:     EPSGetST(ctx->eps,&st);
390:     STSetMatMode(st,ST_MATMODE_SHELL);
391:   }
392:   PetscOptionsTail();
393:   return(0);
394: }

398: static PetscErrorCode PEPLinearSetCompanionForm_Linear(PEP pep,PetscInt cform)
399: {
400:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

403:   if (!cform) return(0);
404:   if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
405:   else {
406:     if (cform!=1 && cform!=2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
407:     ctx->cform = cform;
408:   }
409:   return(0);
410: }

414: /*@
415:    PEPLinearSetCompanionForm - Choose between the two companion forms available
416:    for the linearization of a quadratic eigenproblem.

418:    Logically Collective on PEP

420:    Input Parameters:
421: +  pep   - polynomial eigenvalue solver
422: -  cform - 1 or 2 (first or second companion form)

424:    Options Database Key:
425: .  -pep_linear_cform <int> - Choose the companion form

427:    Level: advanced

429: .seealso: PEPLinearGetCompanionForm()
430: @*/
431: PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform)
432: {

438:   PetscTryMethod(pep,"PEPLinearSetCompanionForm_C",(PEP,PetscInt),(pep,cform));
439:   return(0);
440: }

444: static PetscErrorCode PEPLinearGetCompanionForm_Linear(PEP pep,PetscInt *cform)
445: {
446:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

449:   *cform = ctx->cform;
450:   return(0);
451: }

455: /*@
456:    PEPLinearGetCompanionForm - Returns the number of the companion form that
457:    will be used for the linearization of a quadratic eigenproblem.

459:    Not Collective

461:    Input Parameter:
462: .  pep  - polynomial eigenvalue solver

464:    Output Parameter:
465: .  cform - the companion form number (1 or 2)

467:    Level: advanced

469: .seealso: PEPLinearSetCompanionForm()
470: @*/
471: PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform)
472: {

478:   PetscTryMethod(pep,"PEPLinearGetCompanionForm_C",(PEP,PetscInt*),(pep,cform));
479:   return(0);
480: }

484: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
485: {
486:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

489:   ctx->explicitmatrix = explicitmatrix;
490:   return(0);
491: }

495: /*@
496:    PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
497:    linearization of the problem must be built explicitly.

499:    Logically Collective on PEP

501:    Input Parameters:
502: +  pep      - polynomial eigenvalue solver
503: -  explicit - boolean flag indicating if the matrices are built explicitly

505:    Options Database Key:
506: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

508:    Level: advanced

510: .seealso: PEPLinearGetExplicitMatrix()
511: @*/
512: PetscErrorCode PEPLinearSetExplicitMatrix(PEP pep,PetscBool explicitmatrix)
513: {

519:   PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmatrix));
520:   return(0);
521: }

525: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmatrix)
526: {
527:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

530:   *explicitmatrix = ctx->explicitmatrix;
531:   return(0);
532: }

536: /*@
537:    PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
538:    A and B for the linearization are built explicitly.

540:    Not Collective

542:    Input Parameter:
543: .  pep  - polynomial eigenvalue solver

545:    Output Parameter:
546: .  explicitmatrix - the mode flag

548:    Level: advanced

550: .seealso: PEPLinearSetExplicitMatrix()
551: @*/
552: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmatrix)
553: {

559:   PetscTryMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmatrix));
560:   return(0);
561: }

565: static PetscErrorCode PEPLinearSetEPS_Linear(PEP pep,EPS eps)
566: {
568:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

571:   PetscObjectReference((PetscObject)eps);
572:   EPSDestroy(&ctx->eps);
573:   ctx->eps = eps;
574:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
575:   pep->state = PEP_STATE_INITIAL;
576:   return(0);
577: }

581: /*@
582:    PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
583:    polynomial eigenvalue solver.

585:    Collective on PEP

587:    Input Parameters:
588: +  pep - polynomial eigenvalue solver
589: -  eps - the eigensolver object

591:    Level: advanced

593: .seealso: PEPLinearGetEPS()
594: @*/
595: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
596: {

603:   PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
604:   return(0);
605: }

609: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
610: {
612:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
613:   ST             st;

616:   if (!ctx->eps) {
617:     EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
618:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
619:     EPSAppendOptionsPrefix(ctx->eps,"pep_");
620:     EPSGetST(ctx->eps,&st);
621:     STSetOptionsPrefix(st,((PetscObject)ctx->eps)->prefix);
622:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
623:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
624:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
625:   }
626:   *eps = ctx->eps;
627:   return(0);
628: }

632: /*@
633:    PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
634:    to the polynomial eigenvalue solver.

636:    Not Collective

638:    Input Parameter:
639: .  pep - polynomial eigenvalue solver

641:    Output Parameter:
642: .  eps - the eigensolver object

644:    Level: advanced

646: .seealso: PEPLinearSetEPS()
647: @*/
648: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
649: {

655:   PetscTryMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
656:   return(0);
657: }

661: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
662: {
664:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

667:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
668:   PetscViewerASCIIPrintf(viewer,"  Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
669:   PetscViewerASCIIPrintf(viewer,"  Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
670:   PetscViewerASCIIPushTab(viewer);
671:   EPSView(ctx->eps,viewer);
672:   PetscViewerASCIIPopTab(viewer);
673:   return(0);
674: }

678: PetscErrorCode PEPReset_Linear(PEP pep)
679: {
681:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

684:   if (!ctx->eps) { EPSReset(ctx->eps); }
685:   MatDestroy(&ctx->A);
686:   MatDestroy(&ctx->B);
687:   VecDestroy(&ctx->x1);
688:   VecDestroy(&ctx->x2);
689:   VecDestroy(&ctx->y1);
690:   VecDestroy(&ctx->y2);
691:   return(0);
692: }

696: PetscErrorCode PEPDestroy_Linear(PEP pep)
697: {
699:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

702:   EPSDestroy(&ctx->eps);
703:   PetscFree(pep->data);
704:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",NULL);
705:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",NULL);
706:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
707:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
708:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
709:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
710:   return(0);
711: }

715: PETSC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
716: {
718:   PEP_LINEAR     *ctx;

721:   PetscNewLog(pep,&ctx);
722:   ctx->explicitmatrix = PETSC_TRUE;
723:   pep->data = (void*)ctx;

725:   pep->ops->solve                = PEPSolve_Linear;
726:   pep->ops->setup                = PEPSetUp_Linear;
727:   pep->ops->setfromoptions       = PEPSetFromOptions_Linear;
728:   pep->ops->destroy              = PEPDestroy_Linear;
729:   pep->ops->reset                = PEPReset_Linear;
730:   pep->ops->view                 = PEPView_Linear;
731:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",PEPLinearSetCompanionForm_Linear);
732:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",PEPLinearGetCompanionForm_Linear);
733:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
734:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
735:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
736:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
737:   return(0);
738: }