Actual source code: linear.c

slepc-3.7.1 2016-05-27
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,shift,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:   /* Set STSHIFT as the default ST */
284:   if (!((PetscObject)pep->st)->type_name) {
285:     STSetType(pep->st,STSHIFT);
286:   }
287:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&shift);
288:   PetscObjectTypeCompare((PetscObject)pep->st,STSINVERT,&sinv);
289:   if (!shift && !sinv) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used");
290:   if (!pep->which) {
291:     if (sinv) pep->which = PEP_TARGET_MAGNITUDE;
292:     else pep->which = PEP_LARGEST_MAGNITUDE;
293:   }
294:   STSetUp(pep->st);
295:   if (!ctx->eps) { PEPLinearGetEPS(pep,&ctx->eps); }
296:   EPSGetST(ctx->eps,&st);
297:   if (!transf) { EPSSetTarget(ctx->eps,pep->target); }
298:   if (sinv && !transf) { STSetDefaultShift(st,pep->target); }
299:   /* compute scale factor if not set by user */
300:   PEPComputeScaleFactor(pep);

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

330:     (*fcreate[i][0])(PetscObjectComm((PetscObject)pep),ctx,&ctx->A);
331:     (*fcreate[i][1])(PetscObjectComm((PetscObject)pep),ctx,&ctx->B);
332:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->A);
333:     PetscLogObjectParent((PetscObject)pep,(PetscObject)ctx->B);

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

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

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

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

424:   /* temporary change of target */
425:   if (pep->sfactor!=1.0) {
426:     EPSGetTarget(ctx->eps,&sigma);
427:     EPSSetTarget(ctx->eps,sigma/pep->sfactor);
428:   }

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

447:   EPSSetUp(ctx->eps);
448:   EPSGetDimensions(ctx->eps,NULL,&pep->ncv,&pep->mpd);
449:   EPSGetTolerances(ctx->eps,NULL,&pep->max_it);
450:   if (pep->nini>0) { PetscInfo(pep,"Ignoring initial vectors\n"); }
451:   PEPAllocateSolution(pep,0);
452:   return(0);
453: }

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

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

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

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

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

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

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

699: PetscErrorCode PEPSolve_Linear(PEP pep)
700: {
702:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
703:   PetscScalar    sigma;
704:   PetscBool      flg;
705:   PetscInt       i;

708:   EPSSolve(ctx->eps);
709:   EPSGetConverged(ctx->eps,&pep->nconv);
710:   EPSGetIterationNumber(ctx->eps,&pep->its);
711:   EPSGetConvergedReason(ctx->eps,(EPSConvergedReason*)&pep->reason);

713:   /* recover eigenvalues */
714:   for (i=0;i<pep->nconv;i++) {
715:     EPSGetEigenpair(ctx->eps,i,&pep->eigr[i],&pep->eigi[i],NULL,NULL);
716:     pep->eigr[i] *= pep->sfactor;
717:     pep->eigi[i] *= pep->sfactor;
718:   }

720:   /* restore target */
721:   EPSGetTarget(ctx->eps,&sigma);
722:   EPSSetTarget(ctx->eps,sigma*pep->sfactor);

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

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

752:   PEPMonitor(pep,its,nconv,eigr,eigi,errest,nest);
753:   return(0);
754: }

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

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

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

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

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

803:    Logically Collective on PEP

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

809:    Options Database Key:
810: .  -pep_linear_cform <int> - Choose the companion form

812:    Level: advanced

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

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

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

834:   *cform = ctx->cform;
835:   return(0);
836: }

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

844:    Not Collective

846:    Input Parameter:
847: .  pep  - polynomial eigenvalue solver

849:    Output Parameter:
850: .  cform - the companion form number (1 or 2)

852:    Level: advanced

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

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

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

874:   ctx->explicitmatrix = explicitmatrix;
875:   return(0);
876: }

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

884:    Logically Collective on PEP

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

890:    Options Database Key:
891: .  -pep_linear_explicitmatrix <boolean> - Indicates the boolean flag

893:    Level: advanced

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

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

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

915:   *explicitmatrix = ctx->explicitmatrix;
916:   return(0);
917: }

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

925:    Not Collective

927:    Input Parameter:
928: .  pep  - polynomial eigenvalue solver

930:    Output Parameter:
931: .  explicitmatrix - the mode flag

933:    Level: advanced

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

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

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

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

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

970:    Collective on PEP

972:    Input Parameters:
973: +  pep - polynomial eigenvalue solver
974: -  eps - the eigensolver object

976:    Level: advanced

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

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

994: static PetscErrorCode PEPLinearGetEPS_Linear(PEP pep,EPS *eps)
995: {
997:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
998:   ST             st;

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

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

1021:    Not Collective

1023:    Input Parameter:
1024: .  pep - polynomial eigenvalue solver

1026:    Output Parameter:
1027: .  eps - the eigensolver object

1029:    Level: advanced

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

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

1046: PetscErrorCode PEPView_Linear(PEP pep,PetscViewer viewer)
1047: {
1049:   PEP_LINEAR     *ctx = (PEP_LINEAR*)pep->data;
1050:   PetscBool      isascii;

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

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

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

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

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

1106: PETSC_EXTERN PetscErrorCode PEPCreate_Linear(PEP pep)
1107: {
1109:   PEP_LINEAR     *ctx;

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

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