Actual source code: linear.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    Explicit linearization for polynomial eigenproblems
 12: */

 14: #include <slepc/private/pepimpl.h>         /*I "slepcpep.h" I*/
 15: #include "linearp.h"

 17: static PetscErrorCode MatMult_Linear_Shift(Mat M,Vec x,Vec y)
 18: {
 19:   PetscErrorCode    ierr;
 20:   PEP_LINEAR        *ctx;
 21:   PEP               pep;
 22:   const PetscScalar *px;
 23:   PetscScalar       *py,a,sigma=0.0;
 24:   PetscInt          nmat,deg,i,m;
 25:   Vec               x1,x2,x3,y1,aux;
 26:   PetscReal         *ca,*cb,*cg;
 27:   PetscBool         flg;

 30:   MatShellGetContext(M,(void**)&ctx);
 31:   pep = ctx->pep;
 32:   STGetTransform(pep->st,&flg);
 33:   if (!flg) {
 34:     STGetShift(pep->st,&sigma);
 35:   }
 36:   nmat = pep->nmat;
 37:   deg = nmat-1;
 38:   m = pep->nloc;
 39:   ca = pep->pbc;
 40:   cb = pep->pbc+nmat;
 41:   cg = pep->pbc+2*nmat;
 42:   x1=ctx->w[0];x2=ctx->w[1];x3=ctx->w[2];y1=ctx->w[3];aux=ctx->w[4];

 44:   VecSet(y,0.0);
 45:   VecGetArrayRead(x,&px);
 46:   VecGetArray(y,&py);
 47:   a = 1.0;

 49:   /* first block */
 50:   VecPlaceArray(x2,px);
 51:   VecPlaceArray(x3,px+m);
 52:   VecPlaceArray(y1,py);
 53:   VecAXPY(y1,cb[0]-sigma,x2);
 54:   VecAXPY(y1,ca[0],x3);
 55:   VecResetArray(x2);
 56:   VecResetArray(x3);
 57:   VecResetArray(y1);

 59:   /* inner blocks */
 60:   for (i=1;i<deg-1;i++) {
 61:     VecPlaceArray(x1,px+(i-1)*m);
 62:     VecPlaceArray(x2,px+i*m);
 63:     VecPlaceArray(x3,px+(i+1)*m);
 64:     VecPlaceArray(y1,py+i*m);
 65:     VecAXPY(y1,cg[i],x1);
 66:     VecAXPY(y1,cb[i]-sigma,x2);
 67:     VecAXPY(y1,ca[i],x3);
 68:     VecResetArray(x1);
 69:     VecResetArray(x2);
 70:     VecResetArray(x3);
 71:     VecResetArray(y1);
 72:   }

 74:   /* last block */
 75:   VecPlaceArray(y1,py+(deg-1)*m);
 76:   for (i=0;i<deg;i++) {
 77:     VecPlaceArray(x1,px+i*m);
 78:     STMatMult(pep->st,i,x1,aux);
 79:     VecAXPY(y1,a,aux);
 80:     VecResetArray(x1);
 81:     a *= pep->sfactor;
 82:   }
 83:   VecCopy(y1,aux);
 84:   STMatSolve(pep->st,aux,y1);
 85:   VecScale(y1,-ca[deg-1]/a);
 86:   VecPlaceArray(x1,px+(deg-2)*m);
 87:   VecPlaceArray(x2,px+(deg-1)*m);
 88:   VecAXPY(y1,cg[deg-1],x1);
 89:   VecAXPY(y1,cb[deg-1]-sigma,x2);
 90:   VecResetArray(x1);
 91:   VecResetArray(x2);
 92:   VecResetArray(y1);

 94:   VecRestoreArrayRead(x,&px);
 95:   VecRestoreArray(y,&py);
 96:   return(0);
 97: }

 99: static PetscErrorCode MatMult_Linear_Sinvert(Mat M,Vec x,Vec y)
100: {
101:   PetscErrorCode    ierr;
102:   PEP_LINEAR        *ctx;
103:   PEP               pep;
104:   const PetscScalar *px;
105:   PetscScalar       *py,a,sigma,t=1.0,tp=0.0,tt;
106:   PetscInt          nmat,deg,i,m;
107:   Vec               x1,y1,y2,y3,aux,aux2;
108:   PetscReal         *ca,*cb,*cg;

111:   MatShellGetContext(M,(void**)&ctx);
112:   pep = ctx->pep;
113:   nmat = pep->nmat;
114:   deg = nmat-1;
115:   m = pep->nloc;
116:   ca = pep->pbc;
117:   cb = pep->pbc+nmat;
118:   cg = pep->pbc+2*nmat;
119:   x1=ctx->w[0];y1=ctx->w[1];y2=ctx->w[2];y3=ctx->w[3];aux=ctx->w[4];aux2=ctx->w[5];
120:   EPSGetTarget(ctx->eps,&sigma);
121:   VecSet(y,0.0);
122:   VecGetArrayRead(x,&px);
123:   VecGetArray(y,&py);
124:   a = pep->sfactor;

126:   /* first block */
127:   VecPlaceArray(x1,px);
128:   VecPlaceArray(y1,py+m);
129:   VecCopy(x1,y1);
130:   VecScale(y1,1.0/ca[0]);
131:   VecResetArray(x1);
132:   VecResetArray(y1);

134:   /* second block */
135:   if (deg>2) {
136:     VecPlaceArray(x1,px+m);
137:     VecPlaceArray(y1,py+m);
138:     VecPlaceArray(y2,py+2*m);
139:     VecCopy(x1,y2);
140:     VecAXPY(y2,sigma-cb[1],y1);
141:     VecScale(y2,1.0/ca[1]);
142:     VecResetArray(x1);
143:     VecResetArray(y1);
144:     VecResetArray(y2);
145:   }

147:   /* inner blocks */
148:   for (i=2;i<deg-1;i++) {
149:     VecPlaceArray(x1,px+i*m);
150:     VecPlaceArray(y1,py+(i-1)*m);
151:     VecPlaceArray(y2,py+i*m);
152:     VecPlaceArray(y3,py+(i+1)*m);
153:     VecCopy(x1,y3);
154:     VecAXPY(y3,sigma-cb[i],y2);
155:     VecAXPY(y3,-cg[i],y1);
156:     VecScale(y3,1.0/ca[i]);
157:     VecResetArray(x1);
158:     VecResetArray(y1);
159:     VecResetArray(y2);
160:     VecResetArray(y3);
161:   }

163:   /* last block */
164:   VecPlaceArray(y1,py);
165:   for (i=0;i<deg-2;i++) {
166:     VecPlaceArray(y2,py+(i+1)*m);
167:     STMatMult(pep->st,i+1,y2,aux);
168:     VecAXPY(y1,a,aux);
169:     VecResetArray(y2);
170:     a *= pep->sfactor;
171:   }
172:   i = deg-2;
173:   VecPlaceArray(y2,py+(i+1)*m);
174:   VecPlaceArray(y3,py+i*m);
175:   VecCopy(y2,aux2);
176:   VecAXPY(aux2,cg[i+1]/ca[i+1],y3);
177:   STMatMult(pep->st,i+1,aux2,aux);
178:   VecAXPY(y1,a,aux);
179:   VecResetArray(y2);
180:   VecResetArray(y3);
181:   a *= pep->sfactor;
182:   i = deg-1;
183:   VecPlaceArray(x1,px+i*m);
184:   VecPlaceArray(y3,py+i*m);
185:   VecCopy(x1,aux2);
186:   VecAXPY(aux2,sigma-cb[i],y3);
187:   VecScale(aux2,1.0/ca[i]);
188:   STMatMult(pep->st,i+1,aux2,aux);
189:   VecAXPY(y1,a,aux);
190:   VecResetArray(x1);
191:   VecResetArray(y3);

193:   VecCopy(y1,aux);
194:   STMatSolve(pep->st,aux,y1);
195:   VecScale(y1,-1.0);

197:   /* final update */
198:   for (i=1;i<deg;i++) {
199:     VecPlaceArray(y2,py+i*m);
200:     tt = t;
201:     t = ((sigma-cb[i-1])*t-cg[i-1]*tp)/ca[i-1]; /* i-th basis polynomial */
202:     tp = tt;
203:     VecAXPY(y2,t,y1);
204:     VecResetArray(y2);
205:   }
206:   VecResetArray(y1);

208:   VecRestoreArrayRead(x,&px);
209:   VecRestoreArray(y,&py);
210:   return(0);
211: }

213: static PetscErrorCode BackTransform_Linear(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
214: {
216:   PEP_LINEAR     *ctx;
217:   ST             stctx;

220:   STShellGetContext(st,(void**)&ctx);
221:   PEPGetST(ctx->pep,&stctx);
222:   STBackTransform(stctx,n,eigr,eigi);
223:   return(0);
224: }

226: /*
227:    Dummy backtransform operation
228:  */
229: static PetscErrorCode BackTransform_Skip(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
230: {
232:   return(0);
233: }

235: static PetscErrorCode Apply_Linear(ST st,Vec x,Vec y)
236: {
238:   PEP_LINEAR     *ctx;

241:   STShellGetContext(st,(void**)&ctx);
242:   MatMult(ctx->A,x,y);
243:   return(0);
244: }

246: PetscErrorCode PEPSetUp_Linear(PEP pep)
247: {
249:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
250:   ST             st;
251:   PetscInt       i=0,deg=pep->nmat-1;
252:   EPSWhich       which = EPS_LARGEST_MAGNITUDE;
253:   EPSProblemType ptype;
254:   PetscBool      trackall,istrivial,transf,shift,sinv,ks;
255:   PetscScalar    sigma,*epsarray,*peparray;
256:   Vec            veps,w=NULL;
257:   /* function tables */
258:   PetscErrorCode (*fcreate[][2])(MPI_Comm,PEP_LINEAR*,Mat*) = {
259:     { MatCreateExplicit_Linear_N1A, MatCreateExplicit_Linear_N1B },   /* N1 */
260:     { MatCreateExplicit_Linear_N2A, MatCreateExplicit_Linear_N2B },   /* N2 */
261:     { MatCreateExplicit_Linear_S1A, MatCreateExplicit_Linear_S1B },   /* S1 */
262:     { MatCreateExplicit_Linear_S2A, MatCreateExplicit_Linear_S2B },   /* S2 */
263:     { MatCreateExplicit_Linear_H1A, MatCreateExplicit_Linear_H1B },   /* H1 */
264:     { MatCreateExplicit_Linear_H2A, MatCreateExplicit_Linear_H2B }    /* H2 */
265:   };

268:   if (pep->stopping!=PEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"User-defined stopping test not supported");
269:   pep->lineariz = PETSC_TRUE;
270:   if (!ctx->cform) ctx->cform = 1;
271:   STGetTransform(pep->st,&transf);
272:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
273:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
274:   if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
275:   if (!pep->which) {
276:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
277:     else pep->which = PEP_LARGEST_MAGNITUDE;
278:   }
279:   STSetUp(pep->st);
280:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
281:   EPSGetST(ctx->eps,&st);
282:   if (!transf && !ctx->usereps) { EPSSetTarget(ctx->eps,pep->target); }
283:   if (sinv && !transf && !ctx->usereps) { STSetDefaultShift(st,pep->target); }
284:   /* compute scale factor if not set by user */
285:   PEPComputeScaleFactor(pep);

287:   if (ctx->explicitmatrix) {
288:     if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option is not implemented with st-transform flag active");
289:     if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option only available for quadratic problems");
290:     if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option not implemented for non-monomial bases");
291:     if (pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Diagonal scaling not allowed in PEPLINEAR with explicit matrices");
292:     if (sinv && !transf) { STSetType(st,STSINVERT); }
293:     RGPushScale(pep->rg,1.0/pep->sfactor);
294:     STGetMatrixTransformed(pep->st,0,&ctx->K);
295:     STGetMatrixTransformed(pep->st,1,&ctx->C);
296:     STGetMatrixTransformed(pep->st,2,&ctx->M);
297:     ctx->sfactor = pep->sfactor;
298:     ctx->dsfactor = pep->dsfactor;

300:     MatDestroy(&ctx->A);
301:     MatDestroy(&ctx->B);
302:     VecDestroy(&ctx->w[0]);
303:     VecDestroy(&ctx->w[1]);
304:     VecDestroy(&ctx->w[2]);
305:     VecDestroy(&ctx->w[3]);

307:     switch (pep->problem_type) {
308:       case PEP_GENERAL:    i = 0; break;
309:       case PEP_HERMITIAN:  i = 2; break;
310:       case PEP_GYROSCOPIC: i = 4; break;
311:     }
312:     i += ctx->cform-1;

314:     (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
315:     (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);
316:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
317:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);

319:   } else {   /* implicit matrix */
320:     if (pep->problem_type!=PEP_GENERAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Must use the explicit matrix option if problem type is not general");
321:     if (!((PetscObject)(ctx->eps))->type_name) {
322:       EPSSetType(ctx->eps,EPSKRYLOVSCHUR);
323:     } else {
324:       PetscObjectTypeCompare((PetscObject)ctx->eps,EPSKRYLOVSCHUR,&ks);
325:       if (!ks) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option only implemented for Krylov-Schur");
326:     }
327:     if (ctx->cform!=1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option not available for 2nd companion form");
328:     STSetType(st,STSHELL);
329:     STShellSetContext(st,(PetscObject)ctx);
330:     if (!transf) { STShellSetBackTransform(st,BackTransform_Linear); }
331:     else { STShellSetBackTransform(st,BackTransform_Skip); }
332:     MatCreateVecsEmpty(pep->A[0],&ctx->w[0],&ctx->w[1]);
333:     MatCreateVecsEmpty(pep->A[0],&ctx->w[2],&ctx->w[3]);
334:     MatCreateVecs(pep->A[0],&ctx->w[4],&ctx->w[5]);
335:     PetscLogObjectParents(pep,6,ctx->w);
336:     MatCreateShell(PetscObjectComm((PetscObject)pep),deg*pep->nloc,deg*pep->nloc,deg*pep->n,deg*pep->n,ctx,&ctx->A);
337:     if (sinv && !transf) {
338:       MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Sinvert);
339:     } else {
340:       MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Shift);
341:     }
342:     STShellSetApply(st,Apply_Linear);
343:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
344:     ctx->pep = pep;

346:     PEPBasisCoefficients(pep,pep->pbc);
347:     if (!transf) {
348:       PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
349:       PetscLogObjectMemory((PetscObject)pep,pep->nmat*sizeof(PetscScalar));
350:       if (sinv) {
351:         PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
352:       } else {
353:         for (i=0;i<deg;i++) pep->solvematcoeffs[i] = 0.0;
354:         pep->solvematcoeffs[deg] = 1.0;
355:       }
356:       STScaleShift(pep->st,1.0/pep->sfactor);
357:       RGPushScale(pep->rg,1.0/pep->sfactor);
358:     }
359:     if (pep->sfactor!=1.0) {
360:       for (i=0;i<pep->nmat;i++) {
361:         pep->pbc[pep->nmat+i] /= pep->sfactor;
362:         pep->pbc[2*pep->nmat+i] /= pep->sfactor*pep->sfactor;
363:       }
364:     }
365:   }

367:   EPSSetOperators(ctx->eps,ctx->A,ctx->B);
368:   EPSGetProblemType(ctx->eps,&ptype);
369:   if (!ptype) {
370:     if (ctx->explicitmatrix) {
371:       EPSSetProblemType(ctx->eps,EPS_GNHEP);
372:     } else {
373:       EPSSetProblemType(ctx->eps,EPS_NHEP);
374:     }
375:   }
376:   if (!ctx->usereps) {
377:     if (transf) which = EPS_LARGEST_MAGNITUDE;
378:     else {
379:       switch (pep->which) {
380:           case PEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
381:           case PEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
382:           case PEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
383:           case PEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
384:           case PEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
385:           case PEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
386:           case PEP_TARGET_MAGNITUDE:   which = EPS_TARGET_MAGNITUDE; break;
387:           case PEP_TARGET_REAL:        which = EPS_TARGET_REAL; break;
388:           case PEP_TARGET_IMAGINARY:   which = EPS_TARGET_IMAGINARY; break;
389:           case PEP_WHICH_USER:         which = EPS_WHICH_USER;
390:             EPSSetEigenvalueComparison(ctx->eps,pep->sc->comparison,pep->sc->comparisonctx);
391:             break;
392:       }
393:     }
394:     EPSSetWhichEigenpairs(ctx->eps,which);

396:     EPSSetDimensions(ctx->eps,pep->nev,pep->ncv?pep->ncv:PETSC_DEFAULT,pep->mpd?pep->mpd:PETSC_DEFAULT);
397:     EPSSetTolerances(ctx->eps,pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,pep->max_it?pep->max_it:PETSC_DEFAULT);
398:   }
399:   RGIsTrivial(pep->rg,&istrivial);
400:   if (!istrivial) {
401:     if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEPLINEAR does not support a nontrivial region with st-transform");
402:     EPSSetRG(ctx->eps,pep->rg);
403:   }
404:   /* Transfer the trackall option from pep to eps */
405:   PEPGetTrackAll(pep,&trackall);
406:   EPSSetTrackAll(ctx->eps,trackall);

408:   /* temporary change of target */
409:   if (pep->sfactor!=1.0) {
410:     EPSGetTarget(ctx->eps,&sigma);
411:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
412:   }

414:   /* process initial vector */
415:   if (pep->nini<0) {
416:     VecCreateMPI(PetscObjectComm((PetscObject)ctx->eps),deg*pep->nloc,deg*pep->n,&veps);
417:     VecGetArray(veps,&epsarray);
418:     for (i=0;i<deg;i++) {
419:       if (i<-pep->nini) {
420:         VecGetArray(pep->IS[i],&peparray);
421:         PetscMemcpy(epsarray+i*pep->nloc,peparray,pep->nloc*sizeof(PetscScalar));
422:         VecRestoreArray(pep->IS[i],&peparray);
423:       } else {
424:         if (!w) { VecDuplicate(pep->IS[0],&w); }
425:         VecSetRandom(w,NULL);
426:         VecGetArray(w,&peparray);
427:         PetscMemcpy(epsarray+i*pep->nloc,peparray,pep->nloc*sizeof(PetscScalar));
428:         VecRestoreArray(w,&peparray);
429:       }
430:     }
431:     VecRestoreArray(veps,&epsarray);
432:     EPSSetInitialSpace(ctx->eps,1,&veps);
433:     VecDestroy(&veps);
434:     VecDestroy(&w);
435:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
436:   }

438:   EPSSetUp(ctx->eps);
439:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
440:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
441:   PEPAllocateSolution(pep,0);
442:   return(0);
443: }

445: /*
446:    PEPLinearExtract_Residual - Auxiliary routine that copies the solution of the
447:    linear eigenproblem to the PEP object. The eigenvector of the generalized
448:    problem is supposed to be
449:                                z = [  x  ]
450:                                    [ l*x ]
451:    The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
452:    computed residual norm.
453:    Finally, x is normalized so that ||x||_2 = 1.
454: */
455: static PetscErrorCode PEPLinearExtract_Residual(PEP pep,EPS eps)
456: {
457:   PetscErrorCode    ierr;
458:   PetscInt          i,k;
459:   const PetscScalar *px;
460:   PetscScalar       *er=pep->eigr,*ei=pep->eigi;
461:   PetscReal         rn1,rn2;
462:   Vec               xr,xi=NULL,wr;
463:   Mat               A;
464: #if !defined(PETSC_USE_COMPLEX)
465:   Vec               wi;
466:   const PetscScalar *py;
467: #endif

470: #if defined(PETSC_USE_COMPLEX)
471:   PEPSetWorkVecs(pep,2);
472: #else
473:   PEPSetWorkVecs(pep,4);
474: #endif
475:   EPSGetOperators(eps,&A,NULL);
476:   MatCreateVecs(A,&xr,NULL);
477:   MatCreateVecsEmpty(pep->A[0],&wr,NULL);
478: #if !defined(PETSC_USE_COMPLEX)
479:   VecDuplicate(xr,&xi);
480:   VecDuplicateEmpty(wr,&wi);
481: #endif
482:   for (i=0;i<pep->nconv;i++) {
483:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
484: #if !defined(PETSC_USE_COMPLEX)
485:     if (ei[i]!=0.0) {   /* complex conjugate pair */
486:       VecGetArrayRead(xr,&px);
487:       VecGetArrayRead(xi,&py);
488:       VecPlaceArray(wr,px);
489:       VecPlaceArray(wi,py);
490:       VecNormalizeComplex(wr,wi,PETSC_TRUE,NULL);
491:       PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn1);
492:       BVInsertVec(pep->V,i,wr);
493:       BVInsertVec(pep->V,i+1,wi);
494:       for (k=1;k<pep->nmat-1;k++) {
495:         VecResetArray(wr);
496:         VecResetArray(wi);
497:         VecPlaceArray(wr,px+k*pep->nloc);
498:         VecPlaceArray(wi,py+k*pep->nloc);
499:         VecNormalizeComplex(wr,wi,PETSC_TRUE,NULL);
500:         PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn2);
501:         if (rn1>rn2) {
502:           BVInsertVec(pep->V,i,wr);
503:           BVInsertVec(pep->V,i+1,wi);
504:           rn1 = rn2;
505:         }
506:       }
507:       VecResetArray(wr);
508:       VecResetArray(wi);
509:       VecRestoreArrayRead(xr,&px);
510:       VecRestoreArrayRead(xi,&py);
511:       i++;
512:     } else   /* real eigenvalue */
513: #endif
514:     {
515:       VecGetArrayRead(xr,&px);
516:       VecPlaceArray(wr,px);
517:       VecNormalizeComplex(wr,NULL,PETSC_FALSE,NULL);
518:       PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn1);
519:       BVInsertVec(pep->V,i,wr);
520:       for (k=1;k<pep->nmat-1;k++) {
521:         VecResetArray(wr);
522:         VecPlaceArray(wr,px+k*pep->nloc);
523:         VecNormalizeComplex(wr,NULL,PETSC_FALSE,NULL);
524:         PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn2);
525:         if (rn1>rn2) {
526:           BVInsertVec(pep->V,i,wr);
527:           rn1 = rn2;
528:         }
529:       }
530:       VecResetArray(wr);
531:       VecRestoreArrayRead(xr,&px);
532:     }
533:   }
534:   VecDestroy(&wr);
535:   VecDestroy(&xr);
536: #if !defined(PETSC_USE_COMPLEX)
537:   VecDestroy(&wi);
538:   VecDestroy(&xi);
539: #endif
540:   return(0);
541: }

543: /*
544:    PEPLinearExtract_None - Same as PEPLinearExtract_Norm but always takes
545:    the first block.
546: */
547: static PetscErrorCode PEPLinearExtract_None(PEP pep,EPS eps)
548: {
549:   PetscErrorCode    ierr;
550:   PetscInt          i;
551:   const PetscScalar *px;
552:   Mat               A;
553:   Vec               xr,xi,w;
554: #if !defined(PETSC_USE_COMPLEX)
555:   PetscScalar       *ei=pep->eigi;
556: #endif

559:   EPSGetOperators(eps,&A,NULL);
560:   MatCreateVecs(A,&xr,NULL);
561:   VecDuplicate(xr,&xi);
562:   MatCreateVecsEmpty(pep->A[0],&w,NULL);
563:   for (i=0;i<pep->nconv;i++) {
564:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
565: #if !defined(PETSC_USE_COMPLEX)
566:     if (ei[i]!=0.0) {   /* complex conjugate pair */
567:       VecGetArrayRead(xr,&px);
568:       VecPlaceArray(w,px);
569:       BVInsertVec(pep->V,i,w);
570:       VecResetArray(w);
571:       VecRestoreArrayRead(xr,&px);
572:       VecGetArrayRead(xi,&px);
573:       VecPlaceArray(w,px);
574:       BVInsertVec(pep->V,i+1,w);
575:       VecResetArray(w);
576:       VecRestoreArrayRead(xi,&px);
577:       i++;
578:     } else   /* real eigenvalue */
579: #endif
580:     {
581:       VecGetArrayRead(xr,&px);
582:       VecPlaceArray(w,px);
583:       BVInsertVec(pep->V,i,w);
584:       VecResetArray(w);
585:       VecRestoreArrayRead(xr,&px);
586:     }
587:   }
588:   VecDestroy(&w);
589:   VecDestroy(&xr);
590:   VecDestroy(&xi);
591:   return(0);
592: }

594: /*
595:    PEPLinearExtract_Norm - Auxiliary routine that copies the solution of the
596:    linear eigenproblem to the PEP object. The eigenvector of the generalized
597:    problem is supposed to be
598:                                z = [  x  ]
599:                                    [ l*x ]
600:    If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
601:    Finally, x is normalized so that ||x||_2 = 1.
602: */
603: static PetscErrorCode PEPLinearExtract_Norm(PEP pep,EPS eps)
604: {
605:   PetscErrorCode    ierr;
606:   PetscInt          i,offset;
607:   const PetscScalar *px;
608:   PetscScalar       *er=pep->eigr;
609:   Mat               A;
610:   Vec               xr,xi=NULL,w;
611: #if !defined(PETSC_USE_COMPLEX)
612:   PetscScalar       *ei=pep->eigi;
613: #endif

616:   EPSGetOperators(eps,&A,NULL);
617:   MatCreateVecs(A,&xr,NULL);
618: #if !defined(PETSC_USE_COMPLEX)
619:   VecDuplicate(xr,&xi);
620: #endif
621:   MatCreateVecsEmpty(pep->A[0],&w,NULL);
622:   for (i=0;i<pep->nconv;i++) {
623:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
624:     if (SlepcAbsEigenvalue(er[i],ei[i])>1.0) offset = (pep->nmat-2)*pep->nloc;
625:     else offset = 0;
626: #if !defined(PETSC_USE_COMPLEX)
627:     if (ei[i]!=0.0) {   /* complex conjugate pair */
628:       VecGetArrayRead(xr,&px);
629:       VecPlaceArray(w,px+offset);
630:       BVInsertVec(pep->V,i,w);
631:       VecResetArray(w);
632:       VecRestoreArrayRead(xr,&px);
633:       VecGetArrayRead(xi,&px);
634:       VecPlaceArray(w,px+offset);
635:       BVInsertVec(pep->V,i+1,w);
636:       VecResetArray(w);
637:       VecRestoreArrayRead(xi,&px);
638:       i++;
639:     } else /* real eigenvalue */
640: #endif
641:     {
642:       VecGetArrayRead(xr,&px);
643:       VecPlaceArray(w,px+offset);
644:       BVInsertVec(pep->V,i,w);
645:       VecResetArray(w);
646:       VecRestoreArrayRead(xr,&px);
647:     }
648:   }
649:   VecDestroy(&w);
650:   VecDestroy(&xr);
651: #if !defined(PETSC_USE_COMPLEX)
652:   VecDestroy(&xi);
653: #endif
654:   return(0);
655: }

657: PetscErrorCode PEPExtractVectors_Linear(PEP pep)
658: {
660:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

663:   switch (pep->extract) {
664:   case PEP_EXTRACT_NONE:
665:     PEPLinearExtract_None(pep,ctx->eps);
666:     break;
667:   case PEP_EXTRACT_NORM:
668:     PEPLinearExtract_Norm(pep,ctx->eps);
669:     break;
670:   case PEP_EXTRACT_RESIDUAL:
671:     PEPLinearExtract_Residual(pep,ctx->eps);
672:     break;
673:   case PEP_EXTRACT_STRUCTURED:
674:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
675:   }
676:   return(0);
677: }

679: PetscErrorCode PEPSolve_Linear(PEP pep)
680: {
682:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
683:   PetscScalar    sigma;
684:   PetscBool      flg;
685:   PetscInt       i;

688:   EPSSolve(ctx->eps);
689:   EPSGetConverged(ctx->eps,&pep->nconv);
690:   EPSGetIterationNumber(ctx->eps,&pep->its);
691:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);

693:   /* recover eigenvalues */
694:   for (i=0;i<pep->nconv;i++) {
695:     EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
696:     pep->eigr[i] *= pep->sfactor;
697:     pep->eigi[i] *= pep->sfactor;
698:   }

700:   /* restore target */
701:   EPSGetTarget(ctx->eps,&sigma);
702:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

704:   STGetTransform(pep->st,&flg);
705:   if (flg && pep->ops->backtransform) {
706:     (*pep->ops->backtransform)(pep);
707:   }
708:   if (pep->sfactor!=1.0) {
709:     /* Restore original values */
710:     for (i=0;i<pep->nmat;i++){
711:       pep->pbc[pep->nmat+i] *= pep->sfactor;
712:       pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
713:     }
714:     if (!flg && !ctx->explicitmatrix) {
715:       STScaleShift(pep->st,pep->sfactor);
716:     }
717:   }
718:   if (ctx->explicitmatrix || !flg) {
719:     RGPopScale(pep->rg);
720:   }
721:   return(0);
722: }

724: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
725: {
726:   PEP            pep = (PEP)ctx;

730:   PEPMonitor(pep,its,nconv,eigr,eigi,errest,nest);
731:   return(0);
732: }

734: PetscErrorCode PEPSetFromOptions_Linear(PetscOptionItems *PetscOptionsObject,PEP pep)
735: {
737:   PetscBool      set,val;
738:   PetscInt       i;
739:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

742:   PetscOptionsHead(PetscOptionsObject,"PEP Linear Options");

744:     PetscOptionsInt("-pep_linear_cform","Number of the companion form (1 or 2)","PEPLinearSetCompanionForm",ctx->cform,&i,&set);
745:     if (set) { PEPLinearSetCompanionForm(pep,i); }

747:     PetscOptionsBool("-pep_linear_explicitmatrix","Use explicit matrix in linearization","PEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
748:     if (set) { PEPLinearSetExplicitMatrix(pep,val); }

750:   PetscOptionsTail();

752:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
753:   EPSSetFromOptions(ctx->eps);
754:   return(0);
755: }

757: static PetscErrorCode PEPLinearSetCompanionForm_Linear(PEP pep,PetscInt cform)
758: {
759:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

762:   if (!cform) return(0);
763:   if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
764:   else {
765:     if (cform!=1 && cform!=2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
766:     ctx->cform = cform;
767:   }
768:   return(0);
769: }

771: /*@
772:    PEPLinearSetCompanionForm - Choose between the two companion forms available
773:    for the linearization of a quadratic eigenproblem.

775:    Logically Collective on PEP

777:    Input Parameters:
778: +  pep   - polynomial eigenvalue solver
779: -  cform - 1 or 2 (first or second companion form)

781:    Options Database Key:
782: .  -pep_linear_cform <int> - Choose the companion form

784:    Level: advanced

786: .seealso: PEPLinearGetCompanionForm()
787: @*/
788: PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform)
789: {

795:   PetscTryMethod(pep,"PEPLinearSetCompanionForm_C",(PEP,PetscInt),(pep,cform));
796:   return(0);
797: }

799: static PetscErrorCode PEPLinearGetCompanionForm_Linear(PEP pep,PetscInt *cform)
800: {
801:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

804:   *cform = ctx->cform;
805:   return(0);
806: }

808: /*@
809:    PEPLinearGetCompanionForm - Returns the number of the companion form that
810:    will be used for the linearization of a quadratic eigenproblem.

812:    Not Collective

814:    Input Parameter:
815: .  pep  - polynomial eigenvalue solver

817:    Output Parameter:
818: .  cform - the companion form number (1 or 2)

820:    Level: advanced

822: .seealso: PEPLinearSetCompanionForm()
823: @*/
824: PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform)
825: {

831:   PetscUseMethod(pep,"PEPLinearGetCompanionForm_C",(PEP,PetscInt*),(pep,cform));
832:   return(0);
833: }

835: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
836: {
837:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

840:   if (ctx->explicitmatrix != explicitmatrix) {
841:     ctx->explicitmatrix = explicitmatrix;
842:     pep->state = PEP_STATE_INITIAL;
843:   }
844:   return(0);
845: }

847: /*@
848:    PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
849:    linearization of the problem must be built explicitly.

851:    Logically Collective on PEP

853:    Input Parameters:
854: +  pep      - polynomial eigenvalue solver
855: -  explicit - boolean flag indicating if the matrices are built explicitly

857:    Options Database Key:
858: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

860:    Level: advanced

862: .seealso: PEPLinearGetExplicitMatrix()
863: @*/
864: PetscErrorCode PEPLinearSetExplicitMatrix(PEP pep,PetscBool explicitmatrix)
865: {

871:   PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmatrix));
872:   return(0);
873: }

875: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmatrix)
876: {
877:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

880:   *explicitmatrix = ctx->explicitmatrix;
881:   return(0);
882: }

884: /*@
885:    PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
886:    A and B for the linearization are built explicitly.

888:    Not Collective

890:    Input Parameter:
891: .  pep  - polynomial eigenvalue solver

893:    Output Parameter:
894: .  explicitmatrix - the mode flag

896:    Level: advanced

898: .seealso: PEPLinearSetExplicitMatrix()
899: @*/
900: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmatrix)
901: {

907:   PetscUseMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmatrix));
908:   return(0);
909: }

911: static PetscErrorCode PEPLinearSetEPS_Linear(PEP pep,EPS eps)
912: {
914:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

917:   PetscObjectReference((PetscObject)eps);
918:   EPSDestroy(&ctx->eps);
919:   ctx->eps = eps;
920:   ctx->usereps = PETSC_TRUE;
921:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
922:   pep->state = PEP_STATE_INITIAL;
923:   return(0);
924: }

926: /*@
927:    PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
928:    polynomial eigenvalue solver.

930:    Collective on PEP

932:    Input Parameters:
933: +  pep - polynomial eigenvalue solver
934: -  eps - the eigensolver object

936:    Level: advanced

938: .seealso: PEPLinearGetEPS()
939: @*/
940: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
941: {

948:   PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
949:   return(0);
950: }

952: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
953: {
955:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

958:   if (!ctx->eps) {
959:     EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
960:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
961:     EPSAppendOptionsPrefix(ctx->eps,"pep_linear_");
962:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
963:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
964:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
965:   }
966:   *eps = ctx->eps;
967:   return(0);
968: }

970: /*@
971:    PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
972:    to the polynomial eigenvalue solver.

974:    Not Collective

976:    Input Parameter:
977: .  pep - polynomial eigenvalue solver

979:    Output Parameter:
980: .  eps - the eigensolver object

982:    Level: advanced

984: .seealso: PEPLinearSetEPS()
985: @*/
986: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
987: {

993:   PetscUseMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
994:   return(0);
995: }

997: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
998: {
1000:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
1001:   PetscBool      isascii;

1004:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1005:   if (isascii) {
1006:     if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
1007:     PetscViewerASCIIPrintf(viewer,"  %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
1008:     PetscViewerASCIIPrintf(viewer,"  %s companion form\n",ctx->cform==1? "1st": "2nd");
1009:     PetscViewerASCIIPushTab(viewer);
1010:     EPSView(ctx->eps,viewer);
1011:     PetscViewerASCIIPopTab(viewer);
1012:   }
1013:   return(0);
1014: }

1016: PetscErrorCode PEPReset_Linear(PEP pep)
1017: {
1019:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

1022:   if (!ctx->eps) { EPSReset(ctx->eps); }
1023:   MatDestroy(&ctx->A);
1024:   MatDestroy(&ctx->B);
1025:   VecDestroy(&ctx->w[0]);
1026:   VecDestroy(&ctx->w[1]);
1027:   VecDestroy(&ctx->w[2]);
1028:   VecDestroy(&ctx->w[3]);
1029:   VecDestroy(&ctx->w[4]);
1030:   VecDestroy(&ctx->w[5]);
1031:   return(0);
1032: }

1034: PetscErrorCode PEPDestroy_Linear(PEP pep)
1035: {
1037:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

1040:   EPSDestroy(&ctx->eps);
1041:   PetscFree(pep->data);
1042:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",NULL);
1043:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",NULL);
1044:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
1045:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
1046:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
1047:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
1048:   return(0);
1049: }

1051: PETSC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
1052: {
1054:   PEP_LINEAR     *ctx;

1057:   PetscNewLog(pep,&ctx);
1058:   ctx->explicitmatrix = PETSC_FALSE;
1059:   pep->data = (void*)ctx;

1061:   pep->ops->solve          = PEPSolve_Linear;
1062:   pep->ops->setup          = PEPSetUp_Linear;
1063:   pep->ops->setfromoptions = PEPSetFromOptions_Linear;
1064:   pep->ops->destroy        = PEPDestroy_Linear;
1065:   pep->ops->reset          = PEPReset_Linear;
1066:   pep->ops->view           = PEPView_Linear;
1067:   pep->ops->backtransform  = PEPBackTransform_Default;
1068:   pep->ops->computevectors = PEPComputeVectors_Default;
1069:   pep->ops->extractvectors = PEPExtractVectors_Linear;

1071:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",PEPLinearSetCompanionForm_Linear);
1072:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",PEPLinearGetCompanionForm_Linear);
1073:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
1074:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
1075:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
1076:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
1077:   return(0);
1078: }