Actual source code: linear.c

slepc-3.17.0 2022-03-31
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:   PetscFunctionReturn(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:   PetscFunctionReturn(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:   PetscFunctionReturn(0);
216: }

218: /*
219:    Dummy backtransform operation
220:  */
221: static PetscErrorCode BackTransform_Skip(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
222: {
223:   PetscFunctionReturn(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:   PetscFunctionReturn(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);
296:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
297:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);

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

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

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

369:     EPSSetDimensions(ctx->eps,pep->nev,pep->ncv,pep->mpd);
370:     EPSSetTolerances(ctx->eps,SlepcDefaultTol(pep->tol),pep->max_it);
371:   }
372:   RGIsTrivial(pep->rg,&istrivial);
373:   if (!istrivial) {
375:     EPSSetRG(ctx->eps,pep->rg);
376:   }
377:   /* Transfer the trackall option from pep to eps */
378:   PEPGetTrackAll(pep,&trackall);
379:   EPSSetTrackAll(ctx->eps,trackall);

381:   /* temporary change of target */
382:   if (pep->sfactor!=1.0) {
383:     EPSGetTarget(ctx->eps,&sigma);
384:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
385:   }

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

411:   EPSSetUp(ctx->eps);
412:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
413:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
414:   PEPAllocateSolution(pep,0);
415:   PetscFunctionReturn(0);
416: }

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

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

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

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

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

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

628: PetscErrorCode PEPExtractVectors_Linear(PEP pep)
629: {
630:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

632:   switch (pep->extract) {
633:   case PEP_EXTRACT_NONE:
634:     PEPLinearExtract_None(pep,ctx->eps);
635:     break;
636:   case PEP_EXTRACT_NORM:
637:     PEPLinearExtract_Norm(pep,ctx->eps);
638:     break;
639:   case PEP_EXTRACT_RESIDUAL:
640:     PEPLinearExtract_Residual(pep,ctx->eps);
641:     break;
642:   case PEP_EXTRACT_STRUCTURED:
643:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
644:   }
645:   PetscFunctionReturn(0);
646: }

648: PetscErrorCode PEPSolve_Linear(PEP pep)
649: {
650:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
651:   PetscScalar    sigma;
652:   PetscBool      flg;
653:   PetscInt       i;

655:   EPSSolve(ctx->eps);
656:   EPSGetConverged(ctx->eps,&pep->nconv);
657:   EPSGetIterationNumber(ctx->eps,&pep->its);
658:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);

660:   /* recover eigenvalues */
661:   for (i=0;i<pep->nconv;i++) {
662:     EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
663:     pep->eigr[i] *= pep->sfactor;
664:     pep->eigi[i] *= pep->sfactor;
665:   }

667:   /* restore target */
668:   EPSGetTarget(ctx->eps,&sigma);
669:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

671:   STGetTransform(pep->st,&flg);
672:   if (flg && pep->ops->backtransform) (*pep->ops->backtransform)(pep);
673:   if (pep->sfactor!=1.0) {
674:     /* Restore original values */
675:     for (i=0;i<pep->nmat;i++) {
676:       pep->pbc[pep->nmat+i] *= pep->sfactor;
677:       pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
678:     }
679:     if (!flg && !ctx->explicitmatrix) STScaleShift(pep->st,pep->sfactor);
680:   }
681:   if (ctx->explicitmatrix || !flg) RGPopScale(pep->rg);
682:   PetscFunctionReturn(0);
683: }

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

689:   PEPMonitor(pep,its,nconv,eigr,eigi,errest,nest);
690:   PetscFunctionReturn(0);
691: }

693: PetscErrorCode PEPSetFromOptions_Linear(PetscOptionItems *PetscOptionsObject,PEP pep)
694: {
695:   PetscBool      set,val;
696:   PetscInt       k;
697:   PetscReal      array[2]={0,0};
698:   PetscBool      flg;
699:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

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

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

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

710:   PetscOptionsTail();

712:   if (!ctx->eps) PEPLinearGetEPS(pep,&ctx->eps);
713:   EPSSetFromOptions(ctx->eps);
714:   PetscFunctionReturn(0);
715: }

717: static PetscErrorCode PEPLinearSetLinearization_Linear(PEP pep,PetscReal alpha,PetscReal beta)
718: {
719:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

722:   ctx->alpha = alpha;
723:   ctx->beta  = beta;
724:   PetscFunctionReturn(0);
725: }

727: /*@
728:    PEPLinearSetLinearization - Set the coefficients that define
729:    the linearization of a quadratic eigenproblem.

731:    Logically Collective on pep

733:    Input Parameters:
734: +  pep   - polynomial eigenvalue solver
735: .  alpha - first parameter of the linearization
736: -  beta  - second parameter of the linearization

738:    Options Database Key:
739: .  -pep_linear_linearization <alpha,beta> - Sets the coefficients

741:    Notes:
742:    Cannot pass zero for both alpha and beta. The default values are
743:    alpha=1 and beta=0.

745:    Level: advanced

747: .seealso: PEPLinearGetLinearization()
748: @*/
749: PetscErrorCode PEPLinearSetLinearization(PEP pep,PetscReal alpha,PetscReal beta)
750: {
754:   PetscTryMethod(pep,"PEPLinearSetLinearization_C",(PEP,PetscReal,PetscReal),(pep,alpha,beta));
755:   PetscFunctionReturn(0);
756: }

758: static PetscErrorCode PEPLinearGetLinearization_Linear(PEP pep,PetscReal *alpha,PetscReal *beta)
759: {
760:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

762:   if (alpha) *alpha = ctx->alpha;
763:   if (beta)  *beta  = ctx->beta;
764:   PetscFunctionReturn(0);
765: }

767: /*@
768:    PEPLinearGetLinearization - Returns the coefficients that define
769:    the linearization of a quadratic eigenproblem.

771:    Not Collective

773:    Input Parameter:
774: .  pep  - polynomial eigenvalue solver

776:    Output Parameters:
777: +  alpha - the first parameter of the linearization
778: -  beta  - the second parameter of the linearization

780:    Level: advanced

782: .seealso: PEPLinearSetLinearization()
783: @*/
784: PetscErrorCode PEPLinearGetLinearization(PEP pep,PetscReal *alpha,PetscReal *beta)
785: {
787:   PetscUseMethod(pep,"PEPLinearGetLinearization_C",(PEP,PetscReal*,PetscReal*),(pep,alpha,beta));
788:   PetscFunctionReturn(0);
789: }

791: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
792: {
793:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

795:   if (ctx->explicitmatrix != explicitmatrix) {
796:     ctx->explicitmatrix = explicitmatrix;
797:     pep->state = PEP_STATE_INITIAL;
798:   }
799:   PetscFunctionReturn(0);
800: }

802: /*@
803:    PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
804:    linearization of the problem must be built explicitly.

806:    Logically Collective on pep

808:    Input Parameters:
809: +  pep         - polynomial eigenvalue solver
810: -  explicitmat - boolean flag indicating if the matrices are built explicitly

812:    Options Database Key:
813: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

815:    Level: advanced

817: .seealso: PEPLinearGetExplicitMatrix()
818: @*/
819: PetscErrorCode PEPLinearSetExplicitMatrix(PEP pep,PetscBool explicitmat)
820: {
823:   PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmat));
824:   PetscFunctionReturn(0);
825: }

827: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmat)
828: {
829:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

831:   *explicitmat = ctx->explicitmatrix;
832:   PetscFunctionReturn(0);
833: }

835: /*@
836:    PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
837:    A and B for the linearization are built explicitly.

839:    Not Collective

841:    Input Parameter:
842: .  pep  - polynomial eigenvalue solver

844:    Output Parameter:
845: .  explicitmat - the mode flag

847:    Level: advanced

849: .seealso: PEPLinearSetExplicitMatrix()
850: @*/
851: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmat)
852: {
855:   PetscUseMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmat));
856:   PetscFunctionReturn(0);
857: }

859: static PetscErrorCode PEPLinearSetEPS_Linear(PEP pep,EPS eps)
860: {
861:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

863:   PetscObjectReference((PetscObject)eps);
864:   EPSDestroy(&ctx->eps);
865:   ctx->eps = eps;
866:   ctx->usereps = PETSC_TRUE;
867:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
868:   pep->state = PEP_STATE_INITIAL;
869:   PetscFunctionReturn(0);
870: }

872: /*@
873:    PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
874:    polynomial eigenvalue solver.

876:    Collective on pep

878:    Input Parameters:
879: +  pep - polynomial eigenvalue solver
880: -  eps - the eigensolver object

882:    Level: advanced

884: .seealso: PEPLinearGetEPS()
885: @*/
886: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
887: {
891:   PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
892:   PetscFunctionReturn(0);
893: }

895: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
896: {
897:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

899:   if (!ctx->eps) {
900:     EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
901:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
902:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
903:     EPSAppendOptionsPrefix(ctx->eps,"pep_linear_");
904:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
905:     PetscObjectSetOptions((PetscObject)ctx->eps,((PetscObject)pep)->options);
906:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
907:   }
908:   *eps = ctx->eps;
909:   PetscFunctionReturn(0);
910: }

912: /*@
913:    PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
914:    to the polynomial eigenvalue solver.

916:    Not Collective

918:    Input Parameter:
919: .  pep - polynomial eigenvalue solver

921:    Output Parameter:
922: .  eps - the eigensolver object

924:    Level: advanced

926: .seealso: PEPLinearSetEPS()
927: @*/
928: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
929: {
932:   PetscUseMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
933:   PetscFunctionReturn(0);
934: }

936: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
937: {
938:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
939:   PetscBool      isascii;

941:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
942:   if (isascii) {
943:     if (!ctx->eps) PEPLinearGetEPS(pep,&ctx->eps);
944:     PetscViewerASCIIPrintf(viewer,"  %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
945:     PetscViewerASCIIPrintf(viewer,"  linearization parameters: alpha=%g beta=%g\n",(double)ctx->alpha,(double)ctx->beta);
946:     PetscViewerASCIIPushTab(viewer);
947:     EPSView(ctx->eps,viewer);
948:     PetscViewerASCIIPopTab(viewer);
949:   }
950:   PetscFunctionReturn(0);
951: }

953: PetscErrorCode PEPReset_Linear(PEP pep)
954: {
955:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

957:   if (!ctx->eps) EPSReset(ctx->eps);
958:   MatDestroy(&ctx->A);
959:   MatDestroy(&ctx->B);
960:   VecDestroy(&ctx->w[0]);
961:   VecDestroy(&ctx->w[1]);
962:   VecDestroy(&ctx->w[2]);
963:   VecDestroy(&ctx->w[3]);
964:   VecDestroy(&ctx->w[4]);
965:   VecDestroy(&ctx->w[5]);
966:   PetscFunctionReturn(0);
967: }

969: PetscErrorCode PEPDestroy_Linear(PEP pep)
970: {
971:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

973:   EPSDestroy(&ctx->eps);
974:   PetscFree(pep->data);
975:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetLinearization_C",NULL);
976:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetLinearization_C",NULL);
977:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
978:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
979:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
980:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
981:   PetscFunctionReturn(0);
982: }

984: SLEPC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
985: {
986:   PEP_LINEAR     *ctx;

988:   PetscNewLog(pep,&ctx);
989:   pep->data = (void*)ctx;

991:   pep->lineariz       = PETSC_TRUE;
992:   ctx->explicitmatrix = PETSC_FALSE;
993:   ctx->alpha          = 1.0;
994:   ctx->beta           = 0.0;

996:   pep->ops->solve          = PEPSolve_Linear;
997:   pep->ops->setup          = PEPSetUp_Linear;
998:   pep->ops->setfromoptions = PEPSetFromOptions_Linear;
999:   pep->ops->destroy        = PEPDestroy_Linear;
1000:   pep->ops->reset          = PEPReset_Linear;
1001:   pep->ops->view           = PEPView_Linear;
1002:   pep->ops->backtransform  = PEPBackTransform_Default;
1003:   pep->ops->computevectors = PEPComputeVectors_Default;
1004:   pep->ops->extractvectors = PEPExtractVectors_Linear;

1006:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetLinearization_C",PEPLinearSetLinearization_Linear);
1007:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetLinearization_C",PEPLinearGetLinearization_Linear);
1008:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
1009:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
1010:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
1011:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
1012:   PetscFunctionReturn(0);
1013: }