Actual source code: linear.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*
  2:    Explicit linearization for polynomial eigenproblems.

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

  8:    This file is part of SLEPc.

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

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

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

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

 29: static PetscErrorCode MatMult_Linear_Shift(Mat M,Vec x,Vec y)
 30: {
 31:   PetscErrorCode    ierr;
 32:   PEP_LINEAR        *ctx;
 33:   PEP               pep;
 34:   const PetscScalar *px;
 35:   PetscScalar       *py,a,sigma=0.0;
 36:   PetscInt          nmat,deg,i,m;
 37:   Vec               x1,x2,x3,y1,aux;
 38:   PetscReal         *ca,*cb,*cg;
 39:   PetscBool         flg;

 42:   MatShellGetContext(M,(void**)&ctx);
 43:   pep = ctx->pep;
 44:   STGetTransform(pep->st,&flg);
 45:   if (!flg) {
 46:     STGetShift(pep->st,&sigma);
 47:   }
 48:   nmat = pep->nmat;
 49:   deg = nmat-1;
 50:   m = pep->nloc;
 51:   ca = pep->pbc;
 52:   cb = pep->pbc+nmat;
 53:   cg = pep->pbc+2*nmat;
 54:   x1=ctx->w[0];x2=ctx->w[1];x3=ctx->w[2];y1=ctx->w[3];aux=ctx->w[4];
 55:   
 56:   VecSet(y,0.0);
 57:   VecGetArrayRead(x,&px);
 58:   VecGetArray(y,&py);
 59:   a = 1.0;

 61:   /* first block */
 62:   VecPlaceArray(x2,px);
 63:   VecPlaceArray(x3,px+m);
 64:   VecPlaceArray(y1,py);
 65:   VecAXPY(y1,cb[0]-sigma,x2);
 66:   VecAXPY(y1,ca[0],x3);
 67:   VecResetArray(x2);
 68:   VecResetArray(x3);
 69:   VecResetArray(y1);

 71:   /* inner blocks */
 72:   for (i=1;i<deg-1;i++) {
 73:     VecPlaceArray(x1,px+(i-1)*m);
 74:     VecPlaceArray(x2,px+i*m);
 75:     VecPlaceArray(x3,px+(i+1)*m);
 76:     VecPlaceArray(y1,py+i*m);
 77:     VecAXPY(y1,cg[i],x1);
 78:     VecAXPY(y1,cb[i]-sigma,x2);
 79:     VecAXPY(y1,ca[i],x3);
 80:     VecResetArray(x1);
 81:     VecResetArray(x2);
 82:     VecResetArray(x3);
 83:     VecResetArray(y1);
 84:   }

 86:   /* last block */
 87:   VecPlaceArray(y1,py+(deg-1)*m);
 88:   for (i=0;i<deg;i++) {
 89:     VecPlaceArray(x1,px+i*m);
 90:     STMatMult(pep->st,i,x1,aux);
 91:     VecAXPY(y1,a,aux);
 92:     VecResetArray(x1);
 93:     a *= pep->sfactor;
 94:   }
 95:   VecCopy(y1,aux);
 96:   STMatSolve(pep->st,aux,y1);
 97:   VecScale(y1,-ca[deg-1]/a);
 98:   VecPlaceArray(x1,px+(deg-2)*m);
 99:   VecPlaceArray(x2,px+(deg-1)*m);
100:   VecAXPY(y1,cg[deg-1],x1);
101:   VecAXPY(y1,cb[deg-1]-sigma,x2);
102:   VecResetArray(x1);
103:   VecResetArray(x2);
104:   VecResetArray(y1);

106:   VecRestoreArrayRead(x,&px);
107:   VecRestoreArray(y,&py);
108:   return(0);
109: }

113: static PetscErrorCode MatMult_Linear_Sinvert(Mat M,Vec x,Vec y)
114: {
115:   PetscErrorCode    ierr;
116:   PEP_LINEAR        *ctx;
117:   PEP               pep;
118:   const PetscScalar *px;
119:   PetscScalar       *py,a,sigma,t=1.0,tp=0.0,tt;
120:   PetscInt          nmat,deg,i,m;
121:   Vec               x1,y1,y2,y3,aux,aux2;
122:   PetscReal         *ca,*cb,*cg;

125:   MatShellGetContext(M,(void**)&ctx);
126:   pep = ctx->pep;
127:   nmat = pep->nmat;
128:   deg = nmat-1;
129:   m = pep->nloc;
130:   ca = pep->pbc;
131:   cb = pep->pbc+nmat;
132:   cg = pep->pbc+2*nmat;
133:   x1=ctx->w[0];y1=ctx->w[1];y2=ctx->w[2];y3=ctx->w[3];aux=ctx->w[4];aux2=ctx->w[5];
134:   EPSGetTarget(ctx->eps,&sigma);
135:   VecSet(y,0.0);
136:   VecGetArrayRead(x,&px);
137:   VecGetArray(y,&py);
138:   a = pep->sfactor;

140:   /* first block */
141:   VecPlaceArray(x1,px);
142:   VecPlaceArray(y1,py+m);
143:   VecCopy(x1,y1);
144:   VecScale(y1,1.0/ca[0]);
145:   VecResetArray(x1);
146:   VecResetArray(y1);

148:   /* second block */
149:   if (deg>2) {
150:     VecPlaceArray(x1,px+m);
151:     VecPlaceArray(y1,py+m);
152:     VecPlaceArray(y2,py+2*m);
153:     VecCopy(x1,y2);
154:     VecAXPY(y2,sigma-cb[1],y1);
155:     VecScale(y2,1.0/ca[1]);
156:     VecResetArray(x1);
157:     VecResetArray(y1);
158:     VecResetArray(y2);
159:   }

161:   /* inner blocks */
162:   for (i=2;i<deg-1;i++) {
163:     VecPlaceArray(x1,px+i*m);
164:     VecPlaceArray(y1,py+(i-1)*m);
165:     VecPlaceArray(y2,py+i*m);
166:     VecPlaceArray(y3,py+(i+1)*m);
167:     VecCopy(x1,y3);
168:     VecAXPY(y3,sigma-cb[i],y2);
169:     VecAXPY(y3,-cg[i],y1);
170:     VecScale(y3,1.0/ca[i]);
171:     VecResetArray(x1);
172:     VecResetArray(y1);
173:     VecResetArray(y2);
174:     VecResetArray(y3);
175:   }

177:   /* last block */
178:   VecPlaceArray(y1,py);
179:   for (i=0;i<deg-2;i++) {
180:     VecPlaceArray(y2,py+(i+1)*m);
181:     STMatMult(pep->st,i+1,y2,aux);
182:     VecAXPY(y1,a,aux);
183:     VecResetArray(y2);
184:     a *= pep->sfactor;
185:   }
186:   i = deg-2;
187:   VecPlaceArray(y2,py+(i+1)*m);
188:   VecPlaceArray(y3,py+i*m);
189:   VecCopy(y2,aux2);
190:   VecAXPY(aux2,cg[i+1]/ca[i+1],y3);
191:   STMatMult(pep->st,i+1,aux2,aux);
192:   VecAXPY(y1,a,aux);
193:   VecResetArray(y2);
194:   VecResetArray(y3);
195:   a *= pep->sfactor;
196:   i = deg-1;
197:   VecPlaceArray(x1,px+i*m);
198:   VecPlaceArray(y3,py+i*m);
199:   VecCopy(x1,aux2);
200:   VecAXPY(aux2,sigma-cb[i],y3);
201:   VecScale(aux2,1.0/ca[i]);
202:   STMatMult(pep->st,i+1,aux2,aux);
203:   VecAXPY(y1,a,aux);
204:   VecResetArray(x1);
205:   VecResetArray(y3);

207:   VecCopy(y1,aux);
208:   STMatSolve(pep->st,aux,y1);
209:   VecScale(y1,-1.0);

211:   /* final update */
212:   for (i=1;i<deg;i++) {
213:     VecPlaceArray(y2,py+i*m);
214:     tt = t;
215:     t = ((sigma-cb[i-1])*t-cg[i-1]*tp)/ca[i-1]; /* i-th basis polynomial */
216:     tp = tt;
217:     VecAXPY(y2,t,y1);
218:     VecResetArray(y2);
219:   }
220:   VecResetArray(y1);

222:   VecRestoreArrayRead(x,&px);
223:   VecRestoreArray(y,&py);
224:   return(0);
225: }

229: static PetscErrorCode BackTransform_Linear(ST st,PetscInt n,PetscScalar *eigr,PetscScalar *eigi)
230: {
232:   PEP_LINEAR     *ctx;
233:   ST             stctx;

236:   STShellGetContext(st,(void**)&ctx);
237:   PEPGetST(ctx->pep,&stctx);
238:   STBackTransform(stctx,n,eigr,eigi);
239:   return(0);
240: }

244: static PetscErrorCode Apply_Linear(ST st,Vec x,Vec y)
245: {
247:   PEP_LINEAR     *ctx;

250:   STShellGetContext(st,(void**)&ctx);
251:   MatMult(ctx->A,x,y);
252:   return(0);
253: }

257: PetscErrorCode PEPSetUp_Linear(PEP pep)
258: {
260:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
261:   ST             st;
262:   PetscInt       i=0,deg=pep->nmat-1;
263:   EPSWhich       which;
264:   EPSProblemType ptype;
265:   PetscBool      trackall,istrivial,transf,sinv,ks;
266:   PetscScalar    sigma,*epsarray,*peparray;
267:   Vec            veps;
268:   /* function tables */
269:   PetscErrorCode (*fcreate[][2])(MPI_Comm,PEP_LINEAR*,Mat*) = {
270:     { MatCreateExplicit_Linear_N1A, MatCreateExplicit_Linear_N1B },   /* N1 */
271:     { MatCreateExplicit_Linear_N2A, MatCreateExplicit_Linear_N2B },   /* N2 */
272:     { MatCreateExplicit_Linear_S1A, MatCreateExplicit_Linear_S1B },   /* S1 */
273:     { MatCreateExplicit_Linear_S2A, MatCreateExplicit_Linear_S2B },   /* S2 */
274:     { MatCreateExplicit_Linear_H1A, MatCreateExplicit_Linear_H1B },   /* H1 */
275:     { MatCreateExplicit_Linear_H2A, MatCreateExplicit_Linear_H2B }    /* H2 */
276:   };

279:   if (pep->stopping!=PEPStoppingBasic) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"User-defined stopping test not supported");
280:   pep->lineariz = PETSC_TRUE;
281:   if (!ctx->cform) ctx->cform = 1;
282:   STGetTransform(pep->st,&transf);
283:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
284:   if (!pep->which) {
285:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
286:     else pep->which = PEP_LARGEST_MAGNITUDE;
287:   }
288:   STSetUp(pep->st);
289:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
290:   EPSGetST(ctx->eps,&st);
291:   if (!transf) { EPSSetTarget(ctx->eps,pep->target); }
292:   if (sinv && !transf) { STSetDefaultShift(st,pep->target); }
293:   /* compute scale factor if not set by user */
294:   PEPComputeScaleFactor(pep);

296:   if (ctx->explicitmatrix) {
297:     if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option is not implemented with st-transform flag active");
298:     if (pep->nmat!=3) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option only available for quadratic problems");
299:     if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Explicit matrix option not implemented for non-monomial bases");
300:     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");
301:     if (sinv && !transf) { STSetType(st,STSINVERT); }
302:     RGPushScale(pep->rg,1.0/pep->sfactor);
303:     STGetTOperators(pep->st,0,&ctx->K);
304:     STGetTOperators(pep->st,1,&ctx->C);
305:     STGetTOperators(pep->st,2,&ctx->M);
306:     ctx->sfactor = pep->sfactor;
307:     ctx->dsfactor = pep->dsfactor;
308:   
309:     MatDestroy(&ctx->A);
310:     MatDestroy(&ctx->B);
311:     VecDestroy(&ctx->w[0]);
312:     VecDestroy(&ctx->w[1]);
313:     VecDestroy(&ctx->w[2]);
314:     VecDestroy(&ctx->w[3]);
315:   
316:     switch (pep->problem_type) {
317:       case PEP_GENERAL:    i = 0; break;
318:       case PEP_HERMITIAN:  i = 2; break;
319:       case PEP_GYROSCOPIC: i = 4; break;
320:       default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of pep->problem_type");
321:     }
322:     i += ctx->cform-1;

324:     (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
325:     (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);
326:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
327:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);

329:   } else {   /* implicit matrix */
330:     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");
331:     if (!((PetscObject)(ctx->eps))->type_name) {
332:       EPSSetType(ctx->eps,EPSKRYLOVSCHUR);
333:     } else {
334:       PetscObjectTypeCompare((PetscObject)ctx->eps,EPSKRYLOVSCHUR,&ks);
335:       if (!ks) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option only implemented for Krylov-Schur");
336:     }
337:     if (ctx->cform!=1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Implicit matrix option not available for 2nd companion form");
338:     STSetType(st,STSHELL);
339:     STShellSetContext(st,(PetscObject)ctx);
340:     if (!transf) { STShellSetBackTransform(st,BackTransform_Linear); }
341:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[0]);
342:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[1]);
343:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[2]);
344:     VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&ctx->w[3]);
345:     MatCreateVecs(pep->A[0],&ctx->w[4],NULL);
346:     MatCreateVecs(pep->A[0],&ctx->w[5],NULL);
347:     PetscLogObjectParents(pep,6,ctx->w);
348:     MatCreateShell(PetscObjectComm((PetscObject)pep),deg*pep->nloc,deg*pep->nloc,deg*pep->n,deg*pep->n,ctx,&ctx->A);
349:     PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
350:     if (sinv && !transf) {
351:       MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Sinvert);
352:     } else {
353:       MatShellSetOperation(ctx->A,MATOP_MULT,(void(*)(void))MatMult_Linear_Shift);
354:     }
355:     STShellSetApply(st,Apply_Linear);
356:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
357:     ctx->pep = pep;

359:     PEPBasisCoefficients(pep,pep->pbc);
360:     if (!transf) {
361:       PetscMalloc1(pep->nmat,&pep->solvematcoeffs);
362:       if (sinv) {
363:         PEPEvaluateBasis(pep,pep->target,0,pep->solvematcoeffs,NULL);
364:       } else {
365:         for (i=0;i<deg;i++) pep->solvematcoeffs[i] = 0.0;
366:         pep->solvematcoeffs[deg] = 1.0;
367:       }
368:       STScaleShift(pep->st,1.0/pep->sfactor);
369:       RGPushScale(pep->rg,1.0/pep->sfactor);
370:     }
371:     if (pep->sfactor!=1.0) {
372:       for (i=0;i<pep->nmat;i++) {
373:         pep->pbc[pep->nmat+i] /= pep->sfactor;
374:         pep->pbc[2*pep->nmat+i] /= pep->sfactor*pep->sfactor; 
375:       }
376:     }
377:   }

379:   EPSSetOperators(ctx->eps,ctx->A,ctx->B);
380:   EPSGetProblemType(ctx->eps,&ptype);
381:   if (!ptype) {
382:     if (ctx->explicitmatrix) {
383:       EPSSetProblemType(ctx->eps,EPS_GNHEP);
384:     } else {
385:       EPSSetProblemType(ctx->eps,EPS_NHEP);
386:     }
387:   }
388:   if (transf) which = EPS_LARGEST_MAGNITUDE;
389:   else {
390:     switch (pep->which) {
391:         case PEP_LARGEST_MAGNITUDE:  which = EPS_LARGEST_MAGNITUDE; break;
392:         case PEP_SMALLEST_MAGNITUDE: which = EPS_SMALLEST_MAGNITUDE; break;
393:         case PEP_LARGEST_REAL:       which = EPS_LARGEST_REAL; break;
394:         case PEP_SMALLEST_REAL:      which = EPS_SMALLEST_REAL; break;
395:         case PEP_LARGEST_IMAGINARY:  which = EPS_LARGEST_IMAGINARY; break;
396:         case PEP_SMALLEST_IMAGINARY: which = EPS_SMALLEST_IMAGINARY; break;
397:         case PEP_TARGET_MAGNITUDE:   which = EPS_TARGET_MAGNITUDE; break;
398:         case PEP_TARGET_REAL:        which = EPS_TARGET_REAL; break;
399:         case PEP_TARGET_IMAGINARY:   which = EPS_TARGET_IMAGINARY; break;
400:         case PEP_WHICH_USER:         which = EPS_WHICH_USER;
401:           EPSSetEigenvalueComparison(ctx->eps,pep->sc->comparison,pep->sc->comparisonctx);
402:           break;
403:         default: SETERRQ(PetscObjectComm((PetscObject)pep),1,"Wrong value of which");
404:     }
405:   }
406:   EPSSetWhichEigenpairs(ctx->eps,which);

408:   EPSSetDimensions(ctx->eps,pep->nev,pep->ncv?pep->ncv:PETSC_DEFAULT,pep->mpd?pep->mpd:PETSC_DEFAULT);
409:   EPSSetTolerances(ctx->eps,pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,pep->max_it?pep->max_it:PETSC_DEFAULT);
410:   RGIsTrivial(pep->rg,&istrivial);
411:   if (!istrivial) {
412:     if (transf) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEPLINEAR does not support a nontrivial region with st-transform");
413:     EPSSetRG(ctx->eps,pep->rg);
414:   }
415:   /* Transfer the trackall option from pep to eps */
416:   PEPGetTrackAll(pep,&trackall);
417:   EPSSetTrackAll(ctx->eps,trackall);

419:   /* temporary change of target */
420:   if (pep->sfactor!=1.0) {
421:     EPSGetTarget(ctx->eps,&sigma);
422:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
423:   }

425:   /* process initial vector */
426:   if (pep->nini<=-deg) {
427:     VecCreateMPI(PetscObjectComm((PetscObject)ctx->eps),deg*pep->nloc,deg*pep->n,&veps);
428:     VecGetArray(veps,&epsarray);
429:     for (i=0;i<deg;i++) {
430:       VecGetArray(pep->IS[i],&peparray);
431:       PetscMemcpy(epsarray+i*pep->nloc,peparray,pep->nloc*sizeof(PetscScalar));
432:       VecRestoreArray(pep->IS[i],&peparray);
433:     }
434:     VecRestoreArray(veps,&epsarray);
435:     EPSSetInitialSpace(ctx->eps,1,&veps);
436:     VecDestroy(&veps);
437:   }
438:   if (pep->nini<0) {
439:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
440:   }

442:   EPSSetUp(ctx->eps);
443:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
444:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
445:   if (pep->nini>0) { PetscInfo(pep,"Ignoring initial vectors\n"); }
446:   PEPAllocateSolution(pep,0);
447:   return(0);
448: }

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

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

552: /*
553:    PEPLinearExtract_None - Same as PEPLinearExtract_Norm but always takes
554:    the first block.
555: */
556: static PetscErrorCode PEPLinearExtract_None(PEP pep,EPS eps)
557: {
558:   PetscErrorCode    ierr;
559:   PetscInt          i;
560:   const PetscScalar *px;
561:   Mat               A;
562:   Vec               xr,xi,w;
563: #if !defined(PETSC_USE_COMPLEX)
564:   PetscScalar       *ei=pep->eigi;
565: #endif

568:   EPSGetOperators(eps,&A,NULL);
569:   MatCreateVecs(A,&xr,NULL);
570:   VecDuplicate(xr,&xi);
571:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&w);
572:   for (i=0;i<pep->nconv;i++) {
573:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
574: #if !defined(PETSC_USE_COMPLEX)
575:     if (ei[i]!=0.0) {   /* complex conjugate pair */
576:       VecGetArrayRead(xr,&px);
577:       VecPlaceArray(w,px);
578:       BVInsertVec(pep->V,i,w);
579:       VecResetArray(w);
580:       VecRestoreArrayRead(xr,&px);
581:       VecGetArrayRead(xi,&px);
582:       VecPlaceArray(w,px);
583:       BVInsertVec(pep->V,i+1,w);
584:       VecResetArray(w);
585:       VecRestoreArrayRead(xi,&px);
586:       i++;
587:     } else   /* real eigenvalue */
588: #endif
589:     {
590:       VecGetArrayRead(xr,&px);
591:       VecPlaceArray(w,px);
592:       BVInsertVec(pep->V,i,w);
593:       VecResetArray(w);
594:       VecRestoreArrayRead(xr,&px);
595:     }
596:   }
597:   VecDestroy(&w);
598:   VecDestroy(&xr);
599:   VecDestroy(&xi);
600:   return(0);
601: }

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

627:   EPSGetOperators(eps,&A,NULL);
628:   MatCreateVecs(A,&xr,NULL);
629: #if !defined(PETSC_USE_COMPLEX)
630:   VecDuplicate(xr,&xi);
631: #endif
632:   VecCreateMPIWithArray(PetscObjectComm((PetscObject)pep),1,pep->nloc,pep->n,NULL,&w);
633:   for (i=0;i<pep->nconv;i++) {
634:     EPSGetEigenpair(eps,i,NULL,NULL,xr,xi);
635:     if (SlepcAbsEigenvalue(er[i],ei[i])>1.0) offset = (pep->nmat-2)*pep->nloc;
636:     else offset = 0;
637: #if !defined(PETSC_USE_COMPLEX)
638:     if (ei[i]!=0.0) {   /* complex conjugate pair */
639:       VecGetArrayRead(xr,&px);
640:       VecPlaceArray(w,px+offset);
641:       BVInsertVec(pep->V,i,w);
642:       VecResetArray(w);
643:       VecRestoreArrayRead(xr,&px);
644:       VecGetArrayRead(xi,&px);
645:       VecPlaceArray(w,px+offset);
646:       BVInsertVec(pep->V,i+1,w);
647:       VecResetArray(w);
648:       VecRestoreArrayRead(xi,&px);
649:       i++;
650:     } else /* real eigenvalue */
651: #endif
652:     {
653:       VecGetArrayRead(xr,&px);
654:       VecPlaceArray(w,px+offset);
655:       BVInsertVec(pep->V,i,w);
656:       VecResetArray(w);
657:       VecRestoreArrayRead(xr,&px);
658:     }
659:   }
660:   VecDestroy(&w);
661:   VecDestroy(&xr);
662: #if !defined(PETSC_USE_COMPLEX)
663:   VecDestroy(&xi);
664: #endif
665:   return(0);
666: }

670: PetscErrorCode PEPExtractVectors_Linear(PEP pep)
671: {
673:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
674:   
676:   switch (pep->extract) {
677:   case PEP_EXTRACT_NONE:
678:     PEPLinearExtract_None(pep,ctx->eps);
679:     break;
680:   case PEP_EXTRACT_NORM:
681:     PEPLinearExtract_Norm(pep,ctx->eps);
682:     break;
683:   case PEP_EXTRACT_RESIDUAL:
684:     PEPLinearExtract_Residual(pep,ctx->eps);
685:     break;
686:   default:
687:     SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Extraction not implemented in this solver");
688:   }
689:   return(0);
690: }

694: PetscErrorCode PEPSolve_Linear(PEP pep)
695: {
697:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
698:   PetscScalar    sigma;
699:   PetscBool      flg;
700:   PetscInt       i;

703:   EPSSolve(ctx->eps);
704:   EPSGetConverged(ctx->eps,&pep->nconv);
705:   EPSGetIterationNumber(ctx->eps,&pep->its);
706:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);

708:   /* recover eigenvalues */
709:   for (i=0;i<pep->nconv;i++) {
710:     EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
711:     pep->eigr[i] *= pep->sfactor;
712:     pep->eigi[i] *= pep->sfactor;
713:   }

715:   /* restore target */
716:   EPSGetTarget(ctx->eps,&sigma);
717:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

719:   STGetTransform(pep->st,&flg);
720:   if (flg && pep->ops->backtransform) {
721:     (*pep->ops->backtransform)(pep);
722:   }
723:   if (pep->sfactor!=1.0) {
724:     /* Restore original values */
725:     for (i=0;i<pep->nmat;i++){
726:       pep->pbc[pep->nmat+i] *= pep->sfactor;
727:       pep->pbc[2*pep->nmat+i] *= pep->sfactor*pep->sfactor;
728:     }
729:     if (!flg && !ctx->explicitmatrix) {
730:       STScaleShift(pep->st,pep->sfactor);
731:     } 
732:   }
733:   RGPopScale(pep->rg);
734:   return(0);
735: }

739: static PetscErrorCode EPSMonitor_Linear(EPS eps,PetscInt its,PetscInt nconv,PetscScalar *eigr,PetscScalar *eigi,PetscReal *errest,PetscInt nest,void *ctx)
740: {
741:   PetscInt       i;
742:   PEP            pep = (PEP)ctx;
743:   ST             st;

747:   for (i=0;i<PetscMin(nest,pep->ncv);i++) {
748:     pep->eigr[i] = eigr[i];
749:     pep->eigi[i] = eigi[i];
750:     pep->errest[i] = errest[i];
751:   }
752:   EPSGetST(eps,&st);
753:   STBackTransform(st,nest,pep->eigr,pep->eigi);
754:   PEPMonitor(pep,its,nconv,pep->eigr,pep->eigi,pep->errest,nest);
755:   return(0);
756: }

760: PetscErrorCode PEPSetFromOptions_Linear(PetscOptionItems *PetscOptionsObject,PEP pep)
761: {
763:   PetscBool      set,val;
764:   PetscInt       i;
765:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

768:   PetscOptionsHead(PetscOptionsObject,"PEP Linear Options");
769:   PetscOptionsInt("-pep_linear_cform","Number of the companion form","PEPLinearSetCompanionForm",ctx->cform,&i,&set);
770:   if (set) {
771:     PEPLinearSetCompanionForm(pep,i);
772:   }
773:   PetscOptionsBool("-pep_linear_explicitmatrix","Use explicit matrix in linearization","PEPLinearSetExplicitMatrix",ctx->explicitmatrix,&val,&set);
774:   if (set) {
775:     PEPLinearSetExplicitMatrix(pep,val);
776:   }
777:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
778:   EPSSetFromOptions(ctx->eps);
779:   PetscOptionsTail();
780:   return(0);
781: }

785: static PetscErrorCode PEPLinearSetCompanionForm_Linear(PEP pep,PetscInt cform)
786: {
787:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

790:   if (!cform) return(0);
791:   if (cform==PETSC_DECIDE || cform==PETSC_DEFAULT) ctx->cform = 1;
792:   else {
793:     if (cform!=1 && cform!=2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'cform'");
794:     ctx->cform = cform;
795:   }
796:   return(0);
797: }

801: /*@
802:    PEPLinearSetCompanionForm - Choose between the two companion forms available
803:    for the linearization of a quadratic eigenproblem.

805:    Logically Collective on PEP

807:    Input Parameters:
808: +  pep   - polynomial eigenvalue solver
809: -  cform - 1 or 2 (first or second companion form)

811:    Options Database Key:
812: .  -pep_linear_cform <int> - Choose the companion form

814:    Level: advanced

816: .seealso: PEPLinearGetCompanionForm()
817: @*/
818: PetscErrorCode PEPLinearSetCompanionForm(PEP pep,PetscInt cform)
819: {

825:   PetscTryMethod(pep,"PEPLinearSetCompanionForm_C",(PEP,PetscInt),(pep,cform));
826:   return(0);
827: }

831: static PetscErrorCode PEPLinearGetCompanionForm_Linear(PEP pep,PetscInt *cform)
832: {
833:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

836:   *cform = ctx->cform;
837:   return(0);
838: }

842: /*@
843:    PEPLinearGetCompanionForm - Returns the number of the companion form that
844:    will be used for the linearization of a quadratic eigenproblem.

846:    Not Collective

848:    Input Parameter:
849: .  pep  - polynomial eigenvalue solver

851:    Output Parameter:
852: .  cform - the companion form number (1 or 2)

854:    Level: advanced

856: .seealso: PEPLinearSetCompanionForm()
857: @*/
858: PetscErrorCode PEPLinearGetCompanionForm(PEP pep,PetscInt *cform)
859: {

865:   PetscUseMethod(pep,"PEPLinearGetCompanionForm_C",(PEP,PetscInt*),(pep,cform));
866:   return(0);
867: }

871: static PetscErrorCode PEPLinearSetExplicitMatrix_Linear(PEP pep,PetscBool explicitmatrix)
872: {
873:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

876:   ctx->explicitmatrix = explicitmatrix;
877:   return(0);
878: }

882: /*@
883:    PEPLinearSetExplicitMatrix - Indicate if the matrices A and B for the
884:    linearization of the problem must be built explicitly.

886:    Logically Collective on PEP

888:    Input Parameters:
889: +  pep      - polynomial eigenvalue solver
890: -  explicit - boolean flag indicating if the matrices are built explicitly

892:    Options Database Key:
893: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

895:    Level: advanced

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

906:   PetscTryMethod(pep,"PEPLinearSetExplicitMatrix_C",(PEP,PetscBool),(pep,explicitmatrix));
907:   return(0);
908: }

912: static PetscErrorCode PEPLinearGetExplicitMatrix_Linear(PEP pep,PetscBool *explicitmatrix)
913: {
914:   PEP_LINEAR *ctx = (PEP_LINEAR*)pep->data;

917:   *explicitmatrix = ctx->explicitmatrix;
918:   return(0);
919: }

923: /*@
924:    PEPLinearGetExplicitMatrix - Returns the flag indicating if the matrices
925:    A and B for the linearization are built explicitly.

927:    Not Collective

929:    Input Parameter:
930: .  pep  - polynomial eigenvalue solver

932:    Output Parameter:
933: .  explicitmatrix - the mode flag

935:    Level: advanced

937: .seealso: PEPLinearSetExplicitMatrix()
938: @*/
939: PetscErrorCode PEPLinearGetExplicitMatrix(PEP pep,PetscBool *explicitmatrix)
940: {

946:   PetscUseMethod(pep,"PEPLinearGetExplicitMatrix_C",(PEP,PetscBool*),(pep,explicitmatrix));
947:   return(0);
948: }

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

958:   PetscObjectReference((PetscObject)eps);
959:   EPSDestroy(&ctx->eps);
960:   ctx->eps = eps;
961:   PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
962:   pep->state = PEP_STATE_INITIAL;
963:   return(0);
964: }

968: /*@
969:    PEPLinearSetEPS - Associate an eigensolver object (EPS) to the
970:    polynomial eigenvalue solver.

972:    Collective on PEP

974:    Input Parameters:
975: +  pep - polynomial eigenvalue solver
976: -  eps - the eigensolver object

978:    Level: advanced

980: .seealso: PEPLinearGetEPS()
981: @*/
982: PetscErrorCode PEPLinearSetEPS(PEP pep,EPS eps)
983: {

990:   PetscTryMethod(pep,"PEPLinearSetEPS_C",(PEP,EPS),(pep,eps));
991:   return(0);
992: }

996: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
997: {
999:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
1000:   ST             st;

1003:   if (!ctx->eps) {
1004:     EPSCreate(PetscObjectComm((PetscObject)pep),&ctx->eps);
1005:     EPSSetOptionsPrefix(ctx->eps,((PetscObject)pep)->prefix);
1006:     EPSAppendOptionsPrefix(ctx->eps,"pep_linear_");
1007:     EPSGetST(ctx->eps,&st);
1008:     STSetOptionsPrefix(st,((PetscObject)ctx->eps)->prefix);
1009:     PetscObjectIncrementTabLevel((PetscObject)ctx->eps,(PetscObject)pep,1);
1010:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->eps);
1011:     EPSMonitorSet(ctx->eps,EPSMonitor_Linear,pep,NULL);
1012:   }
1013:   *eps = ctx->eps;
1014:   return(0);
1015: }

1019: /*@
1020:    PEPLinearGetEPS - Retrieve the eigensolver object (EPS) associated
1021:    to the polynomial eigenvalue solver.

1023:    Not Collective

1025:    Input Parameter:
1026: .  pep - polynomial eigenvalue solver

1028:    Output Parameter:
1029: .  eps - the eigensolver object

1031:    Level: advanced

1033: .seealso: PEPLinearSetEPS()
1034: @*/
1035: PetscErrorCode PEPLinearGetEPS(PEP pep,EPS *eps)
1036: {

1042:   PetscUseMethod(pep,"PEPLinearGetEPS_C",(PEP,EPS*),(pep,eps));
1043:   return(0);
1044: }

1048: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
1049: {
1051:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
1052:   PetscBool      isascii;

1055:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1056:   if (isascii) {
1057:     if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
1058:     PetscViewerASCIIPrintf(viewer,"  Linear: %s matrices\n",ctx->explicitmatrix? "explicit": "implicit");
1059:     PetscViewerASCIIPrintf(viewer,"  Linear: %s companion form\n",ctx->cform==1? "1st": "2nd");
1060:     PetscViewerASCIIPushTab(viewer);
1061:     EPSView(ctx->eps,viewer);
1062:     PetscViewerASCIIPopTab(viewer);
1063:   }
1064:   return(0);
1065: }

1069: PetscErrorCode PEPReset_Linear(PEP pep)
1070: {
1072:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

1075:   if (!ctx->eps) { EPSReset(ctx->eps); }
1076:   MatDestroy(&ctx->A);
1077:   MatDestroy(&ctx->B);
1078:   VecDestroy(&ctx->w[0]);
1079:   VecDestroy(&ctx->w[1]);
1080:   VecDestroy(&ctx->w[2]);
1081:   VecDestroy(&ctx->w[3]);
1082:   VecDestroy(&ctx->w[4]);
1083:   VecDestroy(&ctx->w[5]);
1084:   return(0);
1085: }

1089: PetscErrorCode PEPDestroy_Linear(PEP pep)
1090: {
1092:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;

1095:   EPSDestroy(&ctx->eps);
1096:   PetscFree(pep->data);
1097:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",NULL);
1098:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",NULL);
1099:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",NULL);
1100:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",NULL);
1101:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",NULL);
1102:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",NULL);
1103:   return(0);
1104: }

1108: PETSC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
1109: {
1111:   PEP_LINEAR     *ctx;

1114:   PetscNewLog(pep,&ctx);
1115:   ctx->explicitmatrix = PETSC_FALSE;
1116:   pep->data = (void*)ctx;

1118:   pep->ops->solve          = PEPSolve_Linear;
1119:   pep->ops->setup          = PEPSetUp_Linear;
1120:   pep->ops->setfromoptions = PEPSetFromOptions_Linear;
1121:   pep->ops->destroy        = PEPDestroy_Linear;
1122:   pep->ops->reset          = PEPReset_Linear;
1123:   pep->ops->view           = PEPView_Linear;
1124:   pep->ops->backtransform  = PEPBackTransform_Default;
1125:   pep->ops->computevectors = PEPComputeVectors_Default;
1126:   pep->ops->extractvectors = PEPExtractVectors_Linear;
1127:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetCompanionForm_C",PEPLinearSetCompanionForm_Linear);
1128:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetCompanionForm_C",PEPLinearGetCompanionForm_Linear);
1129:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetEPS_C",PEPLinearSetEPS_Linear);
1130:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetEPS_C",PEPLinearGetEPS_Linear);
1131:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearSetExplicitMatrix_C",PEPLinearSetExplicitMatrix_Linear);
1132:   PetscObjectComposeFunction((PetscObject)pep,"PEPLinearGetExplicitMatrix_C",PEPLinearGetExplicitMatrix_Linear);
1133:   return(0);
1134: }