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-2011, 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 <private/qepimpl.h>         /*I "slepcqep.h" I*/
 26: #include <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,qep->nloc,qep->n,PETSC_NULL,&ctx->x1);
 92:     VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->x2);
 93:     VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&ctx->y1);
 94:     VecCreateMPIWithArray(((PetscObject)qep)->comm,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:   EPSSetProblemType(ctx->eps,EPS_GNHEP);
105:   switch (qep->which) {
106:       case QEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
107:       case QEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
108:       case QEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
109:       case QEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
110:       case QEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
111:       case QEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
112:       default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of which");
113:   }
114:   EPSSetWhichEigenpairs(ctx->eps,which);
115:   EPSSetLeftVectorsWanted(ctx->eps,qep->leftvecs);
116:   EPSSetDimensions(ctx->eps,qep->nev,qep->ncv,qep->mpd);
117:   EPSSetTolerances(ctx->eps,qep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/10.0:qep->tol/10.0,qep->max_it);
118:   /* Transfer the trackall option from qep to eps */
119:   QEPGetTrackAll(qep,&trackall);
120:   EPSSetTrackAll(ctx->eps,trackall);
121:   if (ctx->setfromoptionscalled) {
122:     EPSSetFromOptions(ctx->eps);
123:     ctx->setfromoptionscalled = PETSC_FALSE;
124:   }
125:   EPSSetUp(ctx->eps);
126:   EPSGetDimensions(ctx->eps,PETSC_NULL,&qep->ncv,&qep->mpd);
127:   EPSGetTolerances(ctx->eps,PETSC_NULL,&qep->max_it);
128:   if (qep->nini>0 || qep->ninil>0) PetscInfo(qep,"Ignoring initial vectors\n");
129:   QEPAllocateSolution(qep);
130:   return(0);
131: }

135: /*
136:    QEPLinearSelect_Norm - Auxiliary routine that copies the solution of the
137:    linear eigenproblem to the QEP object. The eigenvector of the generalized 
138:    problem is supposed to be
139:                                z = [  x  ]
140:                                    [ l*x ]
141:    The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
142:    computed residual norm.
143:    Finally, x is normalized so that ||x||_2 = 1.
144: */
145: PetscErrorCode QEPLinearSelect_Norm(QEP qep,EPS eps)
146: {
148:   PetscInt       i;
149:   PetscScalar    *px;
150:   PetscReal      rn1,rn2;
151:   Vec            xr,xi,wr,wi;
152:   Mat            A;
153: #if !defined(PETSC_USE_COMPLEX)
154:   PetscScalar    *py;
155: #endif
156: 
158:   EPSGetOperators(eps,&A,PETSC_NULL);
159:   MatGetVecs(A,&xr,PETSC_NULL);
160:   VecDuplicate(xr,&xi);
161:   VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wr);
162:   VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&wi);
163:   for (i=0;i<qep->nconv;i++) {
164:     EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
165:     qep->eigr[i] *= qep->sfactor;
166:     qep->eigi[i] *= qep->sfactor;
167: #if !defined(PETSC_USE_COMPLEX)
168:     if (qep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
169:       VecGetArray(xr,&px);
170:       VecGetArray(xi,&py);
171:       VecPlaceArray(wr,px);
172:       VecPlaceArray(wi,py);
173:       SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
174:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn1);
175:       VecCopy(wr,qep->V[i]);
176:       VecCopy(wi,qep->V[i+1]);
177:       VecResetArray(wr);
178:       VecResetArray(wi);
179:       VecPlaceArray(wr,px+qep->nloc);
180:       VecPlaceArray(wi,py+qep->nloc);
181:       SlepcVecNormalize(wr,wi,PETSC_TRUE,PETSC_NULL);
182:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,wi,&rn2);
183:       if (rn1>rn2) {
184:         VecCopy(wr,qep->V[i]);
185:         VecCopy(wi,qep->V[i+1]);
186:       }
187:       VecResetArray(wr);
188:       VecResetArray(wi);
189:       VecRestoreArray(xr,&px);
190:       VecRestoreArray(xi,&py);
191:     }
192:     else if (qep->eigi[i]==0.0)   /* real eigenvalue */
193: #endif
194:     {
195:       VecGetArray(xr,&px);
196:       VecPlaceArray(wr,px);
197:       SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
198:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn1);
199:       VecCopy(wr,qep->V[i]);
200:       VecResetArray(wr);
201:       VecPlaceArray(wr,px+qep->nloc);
202:       SlepcVecNormalize(wr,PETSC_NULL,PETSC_FALSE,PETSC_NULL);
203:       QEPComputeResidualNorm_Private(qep,qep->eigr[i],qep->eigi[i],wr,PETSC_NULL,&rn2);
204:       if (rn1>rn2) {
205:         VecCopy(wr,qep->V[i]);
206:       }
207:       VecResetArray(wr);
208:       VecRestoreArray(xr,&px);
209:     }
210:   }
211:   VecDestroy(&wr);
212:   VecDestroy(&wi);
213:   VecDestroy(&xr);
214:   VecDestroy(&xi);
215:   return(0);
216: }

220: /*
221:    QEPLinearSelect_Simple - Auxiliary routine that copies the solution of the
222:    linear eigenproblem to the QEP object. The eigenvector of the generalized 
223:    problem is supposed to be
224:                                z = [  x  ]
225:                                    [ l*x ]
226:    If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
227:    Finally, x is normalized so that ||x||_2 = 1.
228: */
229: PetscErrorCode QEPLinearSelect_Simple(QEP qep,EPS eps)
230: {
232:   PetscInt       i,offset;
233:   PetscScalar    *px;
234:   Vec            xr,xi,w;
235:   Mat            A;
236: 
238:   EPSGetOperators(eps,&A,PETSC_NULL);
239:   MatGetVecs(A,&xr,PETSC_NULL);
240:   VecDuplicate(xr,&xi);
241:   VecCreateMPIWithArray(((PetscObject)qep)->comm,qep->nloc,qep->n,PETSC_NULL,&w);
242:   for (i=0;i<qep->nconv;i++) {
243:     EPSGetEigenpair(eps,i,&qep->eigr[i],&qep->eigi[i],xr,xi);
244:     qep->eigr[i] *= qep->sfactor;
245:     qep->eigi[i] *= qep->sfactor;
246:     if (SlepcAbsEigenvalue(qep->eigr[i],qep->eigi[i])>1.0) offset = qep->nloc;
247:     else offset = 0;
248: #if !defined(PETSC_USE_COMPLEX)
249:     if (qep->eigi[i]>0.0) {   /* first eigenvalue of a complex conjugate pair */
250:       VecGetArray(xr,&px);
251:       VecPlaceArray(w,px+offset);
252:       VecCopy(w,qep->V[i]);
253:       VecResetArray(w);
254:       VecRestoreArray(xr,&px);
255:       VecGetArray(xi,&px);
256:       VecPlaceArray(w,px+offset);
257:       VecCopy(w,qep->V[i+1]);
258:       VecResetArray(w);
259:       VecRestoreArray(xi,&px);
260:       SlepcVecNormalize(qep->V[i],qep->V[i+1],PETSC_TRUE,PETSC_NULL);
261:     }
262:     else if (qep->eigi[i]==0.0)   /* real eigenvalue */
263: #endif
264:     {
265:       VecGetArray(xr,&px);
266:       VecPlaceArray(w,px+offset);
267:       VecCopy(w,qep->V[i]);
268:       VecResetArray(w);
269:       VecRestoreArray(xr,&px);
270:       SlepcVecNormalize(qep->V[i],PETSC_NULL,PETSC_FALSE,PETSC_NULL);
271:     }
272:   }
273:   VecDestroy(&w);
274:   VecDestroy(&xr);
275:   VecDestroy(&xi);
276:   return(0);
277: }

281: PetscErrorCode QEPSolve_Linear(QEP qep)
282: {
284:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
285:   PetscBool      flg=PETSC_FALSE;

288:   EPSSolve(ctx->eps);
289:   EPSGetConverged(ctx->eps,&qep->nconv);
290:   EPSGetIterationNumber(ctx->eps,&qep->its);
291:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&qep->reason);
292:   EPSGetOperationCounters(ctx->eps,&qep->matvecs,PETSC_NULL,&qep->linits);
293:   qep->matvecs *= 2;  /* convention: count one matvec for each non-trivial block in A */
294:   PetscOptionsGetBool(((PetscObject)qep)->prefix,"-qep_linear_select_simple",&flg,PETSC_NULL);
295:   if (flg) {
296:     QEPLinearSelect_Simple(qep,ctx->eps);
297:   } else {
298:     QEPLinearSelect_Norm(qep,ctx->eps);
299:   }
300:   return(0);
301: }

305: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
306: {
307:   PetscInt       i;
308:   QEP            qep = (QEP)ctx;

312:   nconv = 0;
313:   for (i=0;i<PetscMin(nest,qep->ncv);i++) {
314:     qep->eigr[i] = eigr[i];
315:     qep->eigi[i] = eigi[i];
316:     qep->errest[i] = errest[i];
317:     if (0.0 < errest[i] && errest[i] < qep->tol) nconv++;
318:   }
319:   STBackTransform(eps->OP,nest,qep->eigr,qep->eigi);
320:   QEPMonitor(qep,its,nconv,qep->eigr,qep->eigi,qep->errest,nest);
321:   return(0);
322: }

326: PetscErrorCode QEPSetFromOptions_Linear(QEP qep)
327: {
329:   PetscBool      set,val;
330:   PetscInt       i;
331:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;
332:   ST             st;

335:   ctx->setfromoptionscalled = PETSC_TRUE;
336:   PetscOptionsHead("QEP Linear Options");
337:   PetscOptionsInt("-qep_linear_cform","Number of the companion form","QEPLinearSetCompanionForm",ctx->cform,&i,&set);
338:   if (set) {
339:     QEPLinearSetCompanionForm(qep,i);
340:   }
341:   PetscOptionsBool("-qep_linear_explicitmatrix","Use explicit matrix in linearization","QEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
342:   if (set) {
343:     QEPLinearSetExplicitMatrix(qep,val);
344:   }
345:   if (!ctx->explicitmatrix) {
346:     /* use as default an ST with shell matrix and Jacobi */
347:     EPSGetST(ctx->eps,&st);
348:     STSetMatMode(st,ST_MATMODE_SHELL);
349:   }
350:   PetscOptionsTail();
351:   return(0);
352: }

357: PetscErrorCode QEPLinearSetCompanionForm_Linear(QEP qep,PetscInt cform)
358: {
359:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;

362:   if (cform==PETSC_IGNORE) return(0);
363:   if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
364:   else {
365:     if (cform!=1 && cform!=2) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
366:     ctx->cform = cform;
367:   }
368:   return(0);
369: }

374: /*@
375:    QEPLinearSetCompanionForm - Choose between the two companion forms available
376:    for the linearization of the quadratic problem.

378:    Logically Collective on QEP

380:    Input Parameters:
381: +  qep   - quadratic eigenvalue solver
382: -  cform - 1 or 2 (first or second companion form)

384:    Options Database Key:
385: .  -qep_linear_cform <int> - Choose the companion form

387:    Level: advanced

389: .seealso: QEPLinearGetCompanionForm()
390: @*/
391: PetscErrorCode QEPLinearSetCompanionForm(QEP qep,PetscInt cform)
392: {

398:   PetscTryMethod(qep,"QEPLinearSetCompanionForm_C",(QEP,PetscInt),(qep,cform));
399:   return(0);
400: }

405: PetscErrorCode QEPLinearGetCompanionForm_Linear(QEP qep,PetscInt *cform)
406: {
407:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;

410:   *cform = ctx->cform;
411:   return(0);
412: }

417: /*@
418:    QEPLinearGetCompanionForm - Returns the number of the companion form that
419:    will be used for the linearization of the quadratic problem.

421:    Not Collective

423:    Input Parameter:
424: .  qep  - quadratic eigenvalue solver

426:    Output Parameter:
427: .  cform - the companion form number (1 or 2)

429:    Level: advanced

431: .seealso: QEPLinearSetCompanionForm()
432: @*/
433: PetscErrorCode QEPLinearGetCompanionForm(QEP qep,PetscInt *cform)
434: {

440:   PetscTryMethod(qep,"QEPLinearGetCompanionForm_C",(QEP,PetscInt*),(qep,cform));
441:   return(0);
442: }

447: PetscErrorCode QEPLinearSetExplicitMatrix_Linear(QEP qep,PetscBool explicitmatrix)
448: {
449:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;

452:   ctx->explicitmatrix = explicitmatrix;
453:   return(0);
454: }

459: /*@
460:    QEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the linearization
461:    of the quadratic problem must be built explicitly.

463:    Logically Collective on QEP

465:    Input Parameters:
466: +  qep      - quadratic eigenvalue solver
467: -  explicit - boolean flag indicating if the matrices are built explicitly

469:    Options Database Key:
470: .  -qep_linear_explicitmatrix <boolean> - Indicates the boolean flag

472:    Level: advanced

474: .seealso: QEPLinearGetExplicitMatrix()
475: @*/
476: PetscErrorCode QEPLinearSetExplicitMatrix(QEP qep,PetscBool explicitmatrix)
477: {

483:   PetscTryMethod(qep,"QEPLinearSetExplicitMatrix_C",(QEP,PetscBool),(qep,explicitmatrix));
484:   return(0);
485: }

490: PetscErrorCode QEPLinearGetExplicitMatrix_Linear(QEP qep,PetscBool *explicitmatrix)
491: {
492:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;

495:   *explicitmatrix = ctx->explicitmatrix;
496:   return(0);
497: }

502: /*@
503:    QEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices A and B
504:    for the linearization of the quadratic problem are built explicitly.

506:    Not Collective

508:    Input Parameter:
509: .  qep  - quadratic eigenvalue solver

511:    Output Parameter:
512: .  explicitmatrix - the mode flag

514:    Level: advanced

516: .seealso: QEPLinearSetExplicitMatrix()
517: @*/
518: PetscErrorCode QEPLinearGetExplicitMatrix(QEP qep,PetscBool *explicitmatrix)
519: {

525:   PetscTryMethod(qep,"QEPLinearGetExplicitMatrix_C",(QEP,PetscBool*),(qep,explicitmatrix));
526:   return(0);
527: }

532: PetscErrorCode QEPLinearSetEPS_Linear(QEP qep,EPS eps)
533: {
535:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;

538:   PetscObjectReference((PetscObject)eps);
539:   EPSDestroy(&ctx->eps);
540:   ctx->eps = eps;
541:   PetscLogObjectParent(qep,ctx->eps);
542:   qep->setupcalled = 0;
543:   return(0);
544: }

549: /*@
550:    QEPLinearSetEPS - Associate an eigensolver object (EPS) to the
551:    quadratic eigenvalue solver. 

553:    Collective on QEP

555:    Input Parameters:
556: +  qep - quadratic eigenvalue solver
557: -  eps - the eigensolver object

559:    Level: advanced

561: .seealso: QEPLinearGetEPS()
562: @*/
563: PetscErrorCode QEPLinearSetEPS(QEP qep,EPS eps)
564: {

571:   PetscTryMethod(qep,"QEPLinearSetEPS_C",(QEP,EPS),(qep,eps));
572:   return(0);
573: }

578: PetscErrorCode QEPLinearGetEPS_Linear(QEP qep,EPS *eps)
579: {
580:   QEP_LINEAR *ctx = (QEP_LINEAR*)qep->data;

583:   *eps = ctx->eps;
584:   return(0);
585: }

590: /*@
591:    QEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
592:    to the quadratic eigenvalue solver.

594:    Not Collective

596:    Input Parameter:
597: .  qep - quadratic eigenvalue solver

599:    Output Parameter:
600: .  eps - the eigensolver object

602:    Level: advanced

604: .seealso: QEPLinearSetEPS()
605: @*/
606: PetscErrorCode QEPLinearGetEPS(QEP qep,EPS *eps)
607: {

613:   PetscTryMethod(qep,"QEPLinearGetEPS_C",(QEP,EPS*),(qep,eps));
614:   return(0);
615: }

619: PetscErrorCode QEPView_Linear(QEP qep,PetscViewer viewer)
620: {
622:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;

625:   PetscViewerASCIIPrintf(viewer,"  Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
626:   PetscViewerASCIIPrintf(viewer,"  Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
627:   PetscViewerASCIIPushTab(viewer);
628:   EPSView(ctx->eps,viewer);
629:   PetscViewerASCIIPopTab(viewer);
630:   return(0);
631: }

635: PetscErrorCode QEPReset_Linear(QEP qep)
636: {
638:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;

641:   EPSReset(ctx->eps);
642:   MatDestroy(&ctx->A);
643:   MatDestroy(&ctx->B);
644:   VecDestroy(&ctx->x1);
645:   VecDestroy(&ctx->x2);
646:   VecDestroy(&ctx->y1);
647:   VecDestroy(&ctx->y2);
648:   QEPDefaultFreeWork(qep);
649:   QEPFreeSolution(qep);
650:   return(0);
651: }

655: PetscErrorCode QEPDestroy_Linear(QEP qep)
656: {
658:   QEP_LINEAR     *ctx = (QEP_LINEAR*)qep->data;

661:   EPSDestroy(&ctx->eps);
662:   PetscFree(qep->data);
663:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","",PETSC_NULL);
664:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","",PETSC_NULL);
665:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","",PETSC_NULL);
666:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","",PETSC_NULL);
667:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","",PETSC_NULL);
668:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","",PETSC_NULL);
669:   return(0);
670: }

675: PetscErrorCode QEPCreate_Linear(QEP qep)
676: {
678:   QEP_LINEAR     *ctx;

681:   PetscNewLog(qep,QEP_LINEAR,&ctx);
682:   qep->data                      = (void *)ctx;
683:   qep->ops->solve                = QEPSolve_Linear;
684:   qep->ops->setup                = QEPSetUp_Linear;
685:   qep->ops->setfromoptions       = QEPSetFromOptions_Linear;
686:   qep->ops->destroy              = QEPDestroy_Linear;
687:   qep->ops->reset                = QEPReset_Linear;
688:   qep->ops->view                 = QEPView_Linear;
689:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetCompanionForm_C","QEPLinearSetCompanionForm_Linear",QEPLinearSetCompanionForm_Linear);
690:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetCompanionForm_C","QEPLinearGetCompanionForm_Linear",QEPLinearGetCompanionForm_Linear);
691:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetEPS_C","QEPLinearSetEPS_Linear",QEPLinearSetEPS_Linear);
692:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetEPS_C","QEPLinearGetEPS_Linear",QEPLinearGetEPS_Linear);
693:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearSetExplicitMatrix_C","QEPLinearSetExplicitMatrix_Linear",QEPLinearSetExplicitMatrix_Linear);
694:   PetscObjectComposeFunctionDynamic((PetscObject)qep,"QEPLinearGetExplicitMatrix_C","QEPLinearGetExplicitMatrix_Linear",QEPLinearGetExplicitMatrix_Linear);

696:   EPSCreate(((PetscObject)qep)->comm,&ctx->eps);
697:   EPSSetOptionsPrefix(ctx->eps,((PetscObject)qep)->prefix);
698:   EPSAppendOptionsPrefix(ctx->eps,"qep_");
699:   PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)qep,1);
700:   PetscLogObjectParent(qep,ctx->eps);
701:   if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
702:   EPSSetIP(ctx->eps,qep->ip);
703:   EPSMonitorSet(ctx->eps,EPSMonitor_Linear,qep,PETSC_NULL);
704:   ctx->explicitmatrix = PETSC_FALSE;
705:   ctx->cform = 1;
706:   ctx->A = PETSC_NULL;
707:   ctx->B = PETSC_NULL;
708:   ctx->x1 = PETSC_NULL;
709:   ctx->x2 = PETSC_NULL;
710:   ctx->y1 = PETSC_NULL;
711:   ctx->y2 = PETSC_NULL;
712:   ctx->setfromoptionscalled = PETSC_FALSE;
713:   return(0);
714: }