Actual source code: linear.c

slepc-3.18.0 2022-10-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>
 15: #include "linear.h"

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

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

 40:   VecSet(y,0.0);
 41:   VecGetArrayRead(x,&px);
 42:   VecGetArray(y,&py);
 43:   a = 1.0;

 45:   /* first block */
 46:   VecPlaceArray(x2,px);
 47:   VecPlaceArray(x3,px+m);
 48:   VecPlaceArray(y1,py);
 49:   VecAXPY(y1,cb[0]-sigma,x2);
 50:   VecAXPY(y1,ca[0],x3);
 51:   VecResetArray(x2);
 52:   VecResetArray(x3);
 53:   VecResetArray(y1);

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

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

 90:   VecRestoreArrayRead(x,&px);
 91:   VecRestoreArray(y,&py);
 92:   return 0;
 93: }

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

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

120:   /* first block */
121:   VecPlaceArray(x1,px);
122:   VecPlaceArray(y1,py+m);
123:   VecCopy(x1,y1);
124:   VecScale(y1,1.0/ca[0]);
125:   VecResetArray(x1);
126:   VecResetArray(y1);

128:   /* second block */
129:   if (deg>2) {
130:     VecPlaceArray(x1,px+m);
131:     VecPlaceArray(y1,py+m);
132:     VecPlaceArray(y2,py+2*m);
133:     VecCopy(x1,y2);
134:     VecAXPY(y2,sigma-cb[1],y1);
135:     VecScale(y2,1.0/ca[1]);
136:     VecResetArray(x1);
137:     VecResetArray(y1);
138:     VecResetArray(y2);
139:   }

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

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

187:   VecCopy(y1,aux);
188:   STMatSolve(pep->st,aux,y1);
189:   VecScale(y1,-1.0);

191:   /* final update */
192:   for (i=1;i<deg;i++) {
193:     VecPlaceArray(y2,py+i*m);
194:     tt = t;
195:     t = ((sigma-cb[i-1])*t-cg[i-1]*tp)/ca[i-1]; /* i-th basis polynomial */
196:     tp = tt;
197:     VecAXPY(y2,t,y1);
198:     VecResetArray(y2);
199:   }
200:   VecResetArray(y1);

202:   VecRestoreArrayRead(x,&px);
203:   VecRestoreArray(y,&py);
204:   return 0;
205: }

207: static PetscErrorCode BackTransform_Linear(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
208: {
209:   PEP_LINEAR     *ctx;
210:   ST             stctx;

212:   STShellGetContext(st,&ctx);
213:   PEPGetST(ctx->pep,&stctx);
214:   STBackTransform(stctx,n,eigr,eigi);
215:   return 0;
216: }

218: /*
219:    Dummy backtransform operation
220:  */
221: static PetscErrorCode BackTransform_Skip(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
222: {
223:   return 0;
224: }

226: static PetscErrorCode Apply_Linear(ST st,Vec x,Vec y)
227: {
228:   PEP_LINEAR     *ctx;

230:   STShellGetContext(st,&ctx);
231:   MatMult(ctx->A,x,y);
232:   return 0;
233: }

235: PetscErrorCode PEPSetUp_Linear(PEP pep)
236: {
237:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
238:   ST             st;
239:   PetscInt       i=0,deg=pep->nmat-1;
240:   EPSWhich       which = EPS_LARGEST_MAGNITUDE;
241:   EPSProblemType ptype;
242:   PetscBool      trackall,istrivial,transf,sinv,ks;
243:   PetscScalar    sigma,*epsarray,*peparray;
244:   Vec            veps,w=NULL;
245:   /* function tables */
246:   PetscErrorCode (*fcreate[][2])(MPI_Comm,PEP_LINEAR*,Mat*) = {
247:     { MatCreateExplicit_Linear_NA, MatCreateExplicit_Linear_NB },
248:     { MatCreateExplicit_Linear_SA, MatCreateExplicit_Linear_SB },
249:     { MatCreateExplicit_Linear_HA, MatCreateExplicit_Linear_HB },
250:   };

252:   PEPCheckShiftSinvert(pep);
253:   PEPCheckUnsupported(pep,PEP_FEATURE_STOPPING);
254:   PEPCheckIgnored(pep,PEP_FEATURE_CONVERGENCE);
255:   STGetTransform(pep->st,&transf);
256:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
257:   if (!pep->which) PEPSetWhichEigenpairs_Default(pep);
259:   STSetUp(pep->st);
260:   if (!ctx->eps) PEPLinearGetEPS(pep,&ctx->eps);
261:   EPSGetST(ctx->eps,&st);
262:   if (!transf && !ctx->usereps) EPSSetTarget(ctx->eps,pep->target);
263:   if (sinv && !transf && !ctx->usereps) STSetDefaultShift(st,pep->target);
264:   /* compute scale factor if not set by user */
265:   PEPComputeScaleFactor(pep);

267:   if (ctx->explicitmatrix) {
268:     PEPCheckQuadraticCondition(pep,PETSC_TRUE," (with explicit matrix)");
269:     PEPCheckUnsupportedCondition(pep,PEP_FEATURE_NONMONOMIAL,PETSC_TRUE," (with explicit matrix)");
272:     if (sinv && !transf) STSetType(st,STSINVERT);
273:     RGPushScale(pep->rg,1.0/pep->sfactor);
274:     STGetMatrixTransformed(pep->st,0,&ctx->K);
275:     STGetMatrixTransformed(pep->st,1,&ctx->C);
276:     STGetMatrixTransformed(pep->st,2,&ctx->M);
277:     ctx->sfactor = pep->sfactor;
278:     ctx->dsfactor = pep->dsfactor;

280:     MatDestroy(&ctx->A);
281:     MatDestroy(&ctx->B);
282:     VecDestroy(&ctx->w[0]);
283:     VecDestroy(&ctx->w[1]);
284:     VecDestroy(&ctx->w[2]);
285:     VecDestroy(&ctx->w[3]);

287:     switch (pep->problem_type) {
288:       case PEP_GENERAL:    i = 0; break;
289:       case PEP_HERMITIAN:
290:       case PEP_HYPERBOLIC: i = 1; break;
291:       case PEP_GYROSCOPIC: i = 2; break;
292:     }

294:     (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
295:     (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);

297:   } else {   /* implicit matrix */
299:     if (!((PetscObject)(ctx->eps))->type_name) EPSSetType(ctx->eps,EPSKRYLOVSCHUR);
300:     else {
301:       PetscObjectTypeCompare((PetscObject)ctx->eps,EPSKRYLOVSCHUR,&ks);
303:     }
305:     STSetType(st,STSHELL);
306:     STShellSetContext(st,ctx);
307:     if (!transf) STShellSetBackTransform(st,BackTransform_Linear);
308:     else STShellSetBackTransform(st,BackTransform_Skip);
309:     MatCreateVecsEmpty(pep->A[0],&ctx->w[0],&ctx->w[1]);
310:     MatCreateVecsEmpty(pep->A[0],&ctx->w[2],&ctx->w[3]);
311:     MatCreateVecs(pep->A[0],&ctx->w[4],&ctx->w[5]);
312:     MatCreateShell(PetscObjectComm((PetscObject)pep),deg*pep->nloc,deg*pep->nloc,deg*pep->n,deg*pep->n,ctx,&ctx->A);
313:     if (sinv && !transf) MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Sinvert);
314:     else MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Shift);
315:     STShellSetApply(st,Apply_Linear);
316:     ctx->pep = pep;

318:     PEPBasisCoefficients(pep,pep->pbc);
319:     if (!transf) {
320:       PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
321:       if (sinv) PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
322:       else {
323:         for (i=0;i<deg;i++) pep->solvematcoeffs[i] = 0.0;
324:         pep->solvematcoeffs[deg] = 1.0;
325:       }
326:       STScaleShift(pep->st,1.0/pep->sfactor);
327:       RGPushScale(pep->rg,1.0/pep->sfactor);
328:     }
329:     if (pep->sfactor!=1.0) {
330:       for (i=0;i<pep->nmat;i++) {
331:         pep->pbc[pep->nmat+i] /= pep->sfactor;
332:         pep->pbc[2*pep->nmat+i] /= pep->sfactor*pep->sfactor;
333:       }
334:     }
335:   }

337:   EPSSetOperators(ctx->eps,ctx->A,ctx->B);
338:   EPSGetProblemType(ctx->eps,&ptype);
339:   if (!ptype) {
340:     if (ctx->explicitmatrix) EPSSetProblemType(ctx->eps,EPS_GNHEP);
341:     else EPSSetProblemType(ctx->eps,EPS_NHEP);
342:   }
343:   if (!ctx->usereps) {
344:     if (transf) which = EPS_LARGEST_MAGNITUDE;
345:     else {
346:       switch (pep->which) {
347:         case PEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
348:         case PEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
349:         case PEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
350:         case PEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
351:         case PEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
352:         case PEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
353:         case PEP_TARGET_MAGNITUDE:   which = EPS_TARGET_MAGNITUDE; break;
354:         case PEP_TARGET_REAL:        which = EPS_TARGET_REAL; break;
355:         case PEP_TARGET_IMAGINARY:   which = EPS_TARGET_IMAGINARY; break;
356:         case PEP_ALL:                which = EPS_ALL; break;
357:         case PEP_WHICH_USER:         which = EPS_WHICH_USER;
358:           EPSSetEigenvalueComparison(ctx->eps,pep->sc->comparison,pep->sc->comparisonctx);
359:           break;
360:       }
361:     }
362:     EPSSetWhichEigenpairs(ctx->eps,which);

364:     EPSSetDimensions(ctx->eps,pep->nev,pep->ncv,pep->mpd);
365:     EPSSetTolerances(ctx->eps,SlepcDefaultTol(pep->tol),pep->max_it);
366:   }
367:   RGIsTrivial(pep->rg,&istrivial);
368:   if (!istrivial) {
370:     EPSSetRG(ctx->eps,pep->rg);
371:   }
372:   /* Transfer the trackall option from pep to eps */
373:   PEPGetTrackAll(pep,&trackall);
374:   EPSSetTrackAll(ctx->eps,trackall);

376:   /* temporary change of target */
377:   if (pep->sfactor!=1.0) {
378:     EPSGetTarget(ctx->eps,&sigma);
379:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
380:   }

382:   /* process initial vector */
383:   if (pep->nini<0) {
384:     VecCreateMPI(PetscObjectComm((PetscObject)ctx->eps),deg*pep->nloc,deg*pep->n,&veps);
385:     VecGetArray(veps,&epsarray);
386:     for (i=0;i<deg;i++) {
387:       if (i<-pep->nini) {
388:         VecGetArray(pep->IS[i],&peparray);
389:         PetscArraycpy(epsarray+i*pep->nloc,peparray,pep->nloc);
390:         VecRestoreArray(pep->IS[i],&peparray);
391:       } else {
392:         if (!w) VecDuplicate(pep->IS[0],&w);
393:         VecSetRandom(w,NULL);
394:         VecGetArray(w,&peparray);
395:         PetscArraycpy(epsarray+i*pep->nloc,peparray,pep->nloc);
396:         VecRestoreArray(w,&peparray);
397:       }
398:     }
399:     VecRestoreArray(veps,&epsarray);
400:     EPSSetInitialSpace(ctx->eps,1,&veps);
401:     VecDestroy(&veps);
402:     VecDestroy(&w);
403:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
404:   }

406:   EPSSetUp(ctx->eps);
407:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
408:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
409:   PEPAllocateSolution(pep,0);
410:   return 0;
411: }

413: /*
414:    PEPLinearExtract_Residual - Auxiliary routine that copies the solution of the
415:    linear eigenproblem to the PEP object. The eigenvector of the generalized
416:    problem is supposed to be
417:                                z = [  x  ]
418:                                    [ l*x ]
419:    The eigenvector is taken from z(1:n) or z(n+1:2*n) depending on the explicitly
420:    computed residual norm.
421:    Finally, x is normalized so that ||x||_2 = 1.
422: */
423: static PetscErrorCode PEPLinearExtract_Residual(PEP pep,EPS eps)
424: {
425:   PetscInt          i,k;
426:   const PetscScalar *px;
427:   PetscScalar       *er=pep->eigr,*ei=pep->eigi;
428:   PetscReal         rn1,rn2;
429:   Vec               xr,xi=NULL,wr;
430:   Mat               A;
431: #if !defined(PETSC_USE_COMPLEX)
432:   Vec               wi;
433:   const PetscScalar *py;
434: #endif

436: #if defined(PETSC_USE_COMPLEX)
437:   PEPSetWorkVecs(pep,2);
438: #else
439:   PEPSetWorkVecs(pep,4);
440: #endif
441:   EPSGetOperators(eps,&A,NULL);
442:   MatCreateVecs(A,&xr,NULL);
443:   MatCreateVecsEmpty(pep->A[0],&wr,NULL);
444: #if !defined(PETSC_USE_COMPLEX)
445:   VecDuplicate(xr,&xi);
446:   VecDuplicateEmpty(wr,&wi);
447: #endif
448:   for (i=0;i<pep->nconv;i++) {
449:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
450: #if !defined(PETSC_USE_COMPLEX)
451:     if (ei[i]!=0.0) {   /* complex conjugate pair */
452:       VecGetArrayRead(xr,&px);
453:       VecGetArrayRead(xi,&py);
454:       VecPlaceArray(wr,px);
455:       VecPlaceArray(wi,py);
456:       VecNormalizeComplex(wr,wi,PETSC_TRUE,NULL);
457:       PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn1);
458:       BVInsertVec(pep->V,i,wr);
459:       BVInsertVec(pep->V,i+1,wi);
460:       for (k=1;k<pep->nmat-1;k++) {
461:         VecResetArray(wr);
462:         VecResetArray(wi);
463:         VecPlaceArray(wr,px+k*pep->nloc);
464:         VecPlaceArray(wi,py+k*pep->nloc);
465:         VecNormalizeComplex(wr,wi,PETSC_TRUE,NULL);
466:         PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,wi,pep->work,&rn2);
467:         if (rn1>rn2) {
468:           BVInsertVec(pep->V,i,wr);
469:           BVInsertVec(pep->V,i+1,wi);
470:           rn1 = rn2;
471:         }
472:       }
473:       VecResetArray(wr);
474:       VecResetArray(wi);
475:       VecRestoreArrayRead(xr,&px);
476:       VecRestoreArrayRead(xi,&py);
477:       i++;
478:     } else   /* real eigenvalue */
479: #endif
480:     {
481:       VecGetArrayRead(xr,&px);
482:       VecPlaceArray(wr,px);
483:       VecNormalizeComplex(wr,NULL,PETSC_FALSE,NULL);
484:       PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn1);
485:       BVInsertVec(pep->V,i,wr);
486:       for (k=1;k<pep->nmat-1;k++) {
487:         VecResetArray(wr);
488:         VecPlaceArray(wr,px+k*pep->nloc);
489:         VecNormalizeComplex(wr,NULL,PETSC_FALSE,NULL);
490:         PEPComputeResidualNorm_Private(pep,er[i],ei[i],wr,NULL,pep->work,&rn2);
491:         if (rn1>rn2) {
492:           BVInsertVec(pep->V,i,wr);
493:           rn1 = rn2;
494:         }
495:       }
496:       VecResetArray(wr);
497:       VecRestoreArrayRead(xr,&px);
498:     }
499:   }
500:   VecDestroy(&wr);
501:   VecDestroy(&xr);
502: #if !defined(PETSC_USE_COMPLEX)
503:   VecDestroy(&wi);
504:   VecDestroy(&xi);
505: #endif
506:   return 0;
507: }

509: /*
510:    PEPLinearExtract_None - Same as PEPLinearExtract_Norm but always takes
511:    the first block.
512: */
513: static PetscErrorCode PEPLinearExtract_None(PEP pep,EPS eps)
514: {
515:   PetscInt          i;
516:   const PetscScalar *px;
517:   Mat               A;
518:   Vec               xr,xi=NULL,w;
519: #if !defined(PETSC_USE_COMPLEX)
520:   PetscScalar       *ei=pep->eigi;
521: #endif

523:   EPSGetOperators(eps,&A,NULL);
524:   MatCreateVecs(A,&xr,NULL);
525: #if !defined(PETSC_USE_COMPLEX)
526:   VecDuplicate(xr,&xi);
527: #endif
528:   MatCreateVecsEmpty(pep->A[0],&w,NULL);
529:   for (i=0;i<pep->nconv;i++) {
530:     EPSGetEigenvector(eps,i,xr,xi);
531: #if !defined(PETSC_USE_COMPLEX)
532:     if (ei[i]!=0.0) {   /* complex conjugate pair */
533:       VecGetArrayRead(xr,&px);
534:       VecPlaceArray(w,px);
535:       BVInsertVec(pep->V,i,w);
536:       VecResetArray(w);
537:       VecRestoreArrayRead(xr,&px);
538:       VecGetArrayRead(xi,&px);
539:       VecPlaceArray(w,px);
540:       BVInsertVec(pep->V,i+1,w);
541:       VecResetArray(w);
542:       VecRestoreArrayRead(xi,&px);
543:       i++;
544:     } else   /* real eigenvalue */
545: #endif
546:     {
547:       VecGetArrayRead(xr,&px);
548:       VecPlaceArray(w,px);
549:       BVInsertVec(pep->V,i,w);
550:       VecResetArray(w);
551:       VecRestoreArrayRead(xr,&px);
552:     }
553:   }
554:   VecDestroy(&w);
555:   VecDestroy(&xr);
556: #if !defined(PETSC_USE_COMPLEX)
557:   VecDestroy(&xi);
558: #endif
559:   return 0;
560: }

562: /*
563:    PEPLinearExtract_Norm - Auxiliary routine that copies the solution of the
564:    linear eigenproblem to the PEP object. The eigenvector of the generalized
565:    problem is supposed to be
566:                                z = [  x  ]
567:                                    [ l*x ]
568:    If |l|<1.0, the eigenvector is taken from z(1:n), otherwise from z(n+1:2*n).
569:    Finally, x is normalized so that ||x||_2 = 1.
570: */
571: static PetscErrorCode PEPLinearExtract_Norm(PEP pep,EPS eps)
572: {
573:   PetscInt          i,offset;
574:   const PetscScalar *px;
575:   PetscScalar       *er=pep->eigr;
576:   Mat               A;
577:   Vec               xr,xi=NULL,w;
578: #if !defined(PETSC_USE_COMPLEX)
579:   PetscScalar       *ei=pep->eigi;
580: #endif

582:   EPSGetOperators(eps,&A,NULL);
583:   MatCreateVecs(A,&xr,NULL);
584: #if !defined(PETSC_USE_COMPLEX)
585:   VecDuplicate(xr,&xi);
586: #endif
587:   MatCreateVecsEmpty(pep->A[0],&w,NULL);
588:   for (i=0;i<pep->nconv;i++) {
589:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
590:     if (SlepcAbsEigenvalue(er[i],ei[i])>1.0) offset = (pep->nmat-2)*pep->nloc;
591:     else offset = 0;
592: #if !defined(PETSC_USE_COMPLEX)
593:     if (ei[i]!=0.0) {   /* complex conjugate pair */
594:       VecGetArrayRead(xr,&px);
595:       VecPlaceArray(w,px+offset);
596:       BVInsertVec(pep->V,i,w);
597:       VecResetArray(w);
598:       VecRestoreArrayRead(xr,&px);
599:       VecGetArrayRead(xi,&px);
600:       VecPlaceArray(w,px+offset);
601:       BVInsertVec(pep->V,i+1,w);
602:       VecResetArray(w);
603:       VecRestoreArrayRead(xi,&px);
604:       i++;
605:     } else /* real eigenvalue */
606: #endif
607:     {
608:       VecGetArrayRead(xr,&px);
609:       VecPlaceArray(w,px+offset);
610:       BVInsertVec(pep->V,i,w);
611:       VecResetArray(w);
612:       VecRestoreArrayRead(xr,&px);
613:     }
614:   }
615:   VecDestroy(&w);
616:   VecDestroy(&xr);
617: #if !defined(PETSC_USE_COMPLEX)
618:   VecDestroy(&xi);
619: #endif
620:   return 0;
621: }

623: PetscErrorCode PEPExtractVectors_Linear(PEP pep)
624: {
625:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

627:   switch (pep->extract) {
628:   case PEP_EXTRACT_NONE:
629:     PEPLinearExtract_None(pep,ctx->eps);
630:     break;
631:   case PEP_EXTRACT_NORM:
632:     PEPLinearExtract_Norm(pep,ctx->eps);
633:     break;
634:   case PEP_EXTRACT_RESIDUAL:
635:     PEPLinearExtract_Residual(pep,ctx->eps);
636:     break;
637:   case PEP_EXTRACT_STRUCTURED:
638:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
639:   }
640:   return 0;
641: }

643: PetscErrorCode PEPSolve_Linear(PEP pep)
644: {
645:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
646:   PetscScalar    sigma;
647:   PetscBool      flg;
648:   PetscInt       i;

650:   EPSSolve(ctx->eps);
651:   EPSGetConverged(ctx->eps,&pep->nconv);
652:   EPSGetIterationNumber(ctx->eps,&pep->its);
653:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);

655:   /* recover eigenvalues */
656:   for (i=0;i<pep->nconv;i++) {
657:     EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
658:     pep->eigr[i] *= pep->sfactor;
659:     pep->eigi[i] *= pep->sfactor;
660:   }

662:   /* restore target */
663:   EPSGetTarget(ctx->eps,&sigma);
664:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

666:   STGetTransform(pep->st,&flg);
667:   if (flg) PetscTryTypeMethod(pep,backtransform);
668:   if (pep->sfactor!=1.0) {
669:     /* Restore original values */
670:     for (i=0;i<pep->nmat;i++) {
671:       pep->pbc[pep->nmat+i] *= pep->sfactor;
672:       pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
673:     }
674:     if (!flg && !ctx->explicitmatrix) STScaleShift(pep->st,pep->sfactor);
675:   }
676:   if (ctx->explicitmatrix || !flg) RGPopScale(pep->rg);
677:   return 0;
678: }

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

684:   PEPMonitor(pep,its,nconv,eigr,eigi,errest,nest);
685:   return 0;
686: }

688: PetscErrorCode PEPSetFromOptions_Linear(PEP pep,PetscOptionItems *PetscOptionsObject)
689: {
690:   PetscBool      set,val;
691:   PetscInt       k;
692:   PetscReal      array[2]={0,0};
693:   PetscBool      flg;
694:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

696:   PetscOptionsHeadBegin(PetscOptionsObject,"PEP Linear Options");

698:     k = 2;
699:     PetscOptionsRealArray("-pep_linear_linearization","Parameters of the linearization","PEPLinearSetLinearization",array,&k,&flg);
700:     if (flg) PEPLinearSetLinearization(pep,array[0],array[1]);

702:     PetscOptionsBool("-pep_linear_explicitmatrix","Use explicit matrix in linearization","PEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
703:     if (set) PEPLinearSetExplicitMatrix(pep,val);

705:   PetscOptionsHeadEnd();

707:   if (!ctx->eps) PEPLinearGetEPS(pep,&ctx->eps);
708:   EPSSetFromOptions(ctx->eps);
709:   return 0;
710: }

712: static PetscErrorCode PEPLinearSetLinearization_Linear(PEP pep,PetscReal alpha,PetscReal beta)
713: {
714:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

717:   ctx->alpha = alpha;
718:   ctx->beta  = beta;
719:   return 0;
720: }

722: /*@
723:    PEPLinearSetLinearization - Set the coefficients that define
724:    the linearization of a quadratic eigenproblem.

726:    Logically Collective on pep

728:    Input Parameters:
729: +  pep   - polynomial eigenvalue solver
730: .  alpha - first parameter of the linearization
731: -  beta  - second parameter of the linearization

733:    Options Database Key:
734: .  -pep_linear_linearization <alpha,beta> - Sets the coefficients

736:    Notes:
737:    Cannot pass zero for both alpha and beta. The default values are
738:    alpha=1 and beta=0.

740:    Level: advanced

742: .seealso: PEPLinearGetLinearization()
743: @*/
744: PetscErrorCode PEPLinearSetLinearization(PEP pep,PetscReal alpha,PetscReal beta)
745: {
749:   PetscTryMethod(pep,"PEPLinearSetLinearization_C",(PEP,PetscReal,PetscReal),(pep,alpha,beta));
750:   return 0;
751: }

753: static PetscErrorCode PEPLinearGetLinearization_Linear(PEP pep,PetscReal *alpha,PetscReal *beta)
754: {
755:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

757:   if (alpha) *alpha = ctx->alpha;
758:   if (beta)  *beta  = ctx->beta;
759:   return 0;
760: }

762: /*@
763:    PEPLinearGetLinearization - Returns the coefficients that define
764:    the linearization of a quadratic eigenproblem.

766:    Not Collective

768:    Input Parameter:
769: .  pep  - polynomial eigenvalue solver

771:    Output Parameters:
772: +  alpha - the first parameter of the linearization
773: -  beta  - the second parameter of the linearization

775:    Level: advanced

777: .seealso: PEPLinearSetLinearization()
778: @*/
779: PetscErrorCode PEPLinearGetLinearization(PEP pep,PetscReal *alpha,PetscReal *beta)
780: {
782:   PetscUseMethod(pep,"PEPLinearGetLinearization_C",(PEP,PetscReal*,PetscReal*),(pep,alpha,beta));
783:   return 0;
784: }

786: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
787: {
788:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

790:   if (ctx->explicitmatrix != explicitmatrix) {
791:     ctx->explicitmatrix = explicitmatrix;
792:     pep->state = PEP_STATE_INITIAL;
793:   }
794:   return 0;
795: }

797: /*@
798:    PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
799:    linearization of the problem must be built explicitly.

801:    Logically Collective on pep

803:    Input Parameters:
804: +  pep         - polynomial eigenvalue solver
805: -  explicitmat - boolean flag indicating if the matrices are built explicitly

807:    Options Database Key:
808: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

810:    Level: advanced

812: .seealso: PEPLinearGetExplicitMatrix()
813: @*/
814: PetscErrorCode PEPLinearSetExplicitMatrix(PEP pep,PetscBool explicitmat)
815: {
818:   PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmat));
819:   return 0;
820: }

822: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmat)
823: {
824:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

826:   *explicitmat = ctx->explicitmatrix;
827:   return 0;
828: }

830: /*@
831:    PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
832:    A and B for the linearization are built explicitly.

834:    Not Collective

836:    Input Parameter:
837: .  pep  - polynomial eigenvalue solver

839:    Output Parameter:
840: .  explicitmat - the mode flag

842:    Level: advanced

844: .seealso: PEPLinearSetExplicitMatrix()
845: @*/
846: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmat)
847: {
850:   PetscUseMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmat));
851:   return 0;
852: }

854: static PetscErrorCode PEPLinearSetEPS_Linear(PEP pep,EPS eps)
855: {
856:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

858:   PetscObjectReference((PetscObject)eps);
859:   EPSDestroy(&ctx->eps);
860:   ctx->eps     = eps;
861:   ctx->usereps = PETSC_TRUE;
862:   pep->state   = PEP_STATE_INITIAL;
863:   return 0;
864: }

866: /*@
867:    PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
868:    polynomial eigenvalue solver.

870:    Collective on pep

872:    Input Parameters:
873: +  pep - polynomial eigenvalue solver
874: -  eps - the eigensolver object

876:    Level: advanced

878: .seealso: PEPLinearGetEPS()
879: @*/
880: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
881: {
885:   PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
886:   return 0;
887: }

889: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
890: {
891:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

893:   if (!ctx->eps) {
894:     EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
895:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
896:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
897:     EPSAppendOptionsPrefix(ctx->eps,"pep_linear_");
898:     PetscObjectSetOptions((PetscObject)ctx->eps,((PetscObject)pep)->options);
899:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
900:   }
901:   *eps = ctx->eps;
902:   return 0;
903: }

905: /*@
906:    PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
907:    to the polynomial eigenvalue solver.

909:    Not Collective

911:    Input Parameter:
912: .  pep - polynomial eigenvalue solver

914:    Output Parameter:
915: .  eps - the eigensolver object

917:    Level: advanced

919: .seealso: PEPLinearSetEPS()
920: @*/
921: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
922: {
925:   PetscUseMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
926:   return 0;
927: }

929: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
930: {
931:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
932:   PetscBool      isascii;

934:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
935:   if (isascii) {
936:     if (!ctx->eps) PEPLinearGetEPS(pep,&ctx->eps);
937:     PetscViewerASCIIPrintf(viewer,"  %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
938:     PetscViewerASCIIPrintf(viewer,"  linearization parameters: alpha=%g beta=%g\n",(double)ctx->alpha,(double)ctx->beta);
939:     PetscViewerASCIIPushTab(viewer);
940:     EPSView(ctx->eps,viewer);
941:     PetscViewerASCIIPopTab(viewer);
942:   }
943:   return 0;
944: }

946: PetscErrorCode PEPReset_Linear(PEP pep)
947: {
948:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

950:   if (!ctx->eps) EPSReset(ctx->eps);
951:   MatDestroy(&ctx->A);
952:   MatDestroy(&ctx->B);
953:   VecDestroy(&ctx->w[0]);
954:   VecDestroy(&ctx->w[1]);
955:   VecDestroy(&ctx->w[2]);
956:   VecDestroy(&ctx->w[3]);
957:   VecDestroy(&ctx->w[4]);
958:   VecDestroy(&ctx->w[5]);
959:   return 0;
960: }

962: PetscErrorCode PEPDestroy_Linear(PEP pep)
963: {
964:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

966:   EPSDestroy(&ctx->eps);
967:   PetscFree(pep->data);
968:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetLinearization_C",NULL);
969:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetLinearization_C",NULL);
970:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
971:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
972:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
973:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
974:   return 0;
975: }

977: SLEPC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
978: {
979:   PEP_LINEAR     *ctx;

981:   PetscNew(&ctx);
982:   pep->data = (void*)ctx;

984:   pep->lineariz       = PETSC_TRUE;
985:   ctx->explicitmatrix = PETSC_FALSE;
986:   ctx->alpha          = 1.0;
987:   ctx->beta           = 0.0;

989:   pep->ops->solve          = PEPSolve_Linear;
990:   pep->ops->setup          = PEPSetUp_Linear;
991:   pep->ops->setfromoptions = PEPSetFromOptions_Linear;
992:   pep->ops->destroy        = PEPDestroy_Linear;
993:   pep->ops->reset          = PEPReset_Linear;
994:   pep->ops->view           = PEPView_Linear;
995:   pep->ops->backtransform  = PEPBackTransform_Default;
996:   pep->ops->computevectors = PEPComputeVectors_Default;
997:   pep->ops->extractvectors = PEPExtractVectors_Linear;

999:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetLinearization_C",PEPLinearSetLinearization_Linear);
1000:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetLinearization_C",PEPLinearGetLinearization_Linear);
1001:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
1002:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
1003:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
1004:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
1005:   return 0;
1006: }