Actual source code: pjd.c

slepc-3.7.0 2016-05-16
Report Typos and Errors
  1: /*

  3:    SLEPc polynomial eigensolver: "jd"

  5:    Method: Jacobi-Davidson

  7:    Algorithm:

  9:        Jacobi-Davidson for polynomial eigenvalue problems.
 10:        Based on code contributed by the authors of [2] below.

 12:    References:

 14:        [1] G.L.G. Sleijpen et al., "Jacobi-Davidson type methods for
 15:            generalized eigenproblems and polynomial eigenproblems", BIT
 16:            36(3):595-633, 1996.

 18:        [2] Feng-Nan Hwang, Zih-Hao Wei, Tsung-Ming Huang, Weichung Wang,
 19:            "A Parallel Additive Schwarz Preconditioned Jacobi-Davidson
 20:            Algorithm for Polynomial Eigenvalue Problems in Quantum Dot
 21:            Simulation", J. Comput. Phys. 229(8):2932-2947, 2010.

 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 25:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

 27:    This file is part of SLEPc.

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

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

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

 43: #include <slepc/private/pepimpl.h>    /*I "slepcpep.h" I*/
 44:  #include pjdp.h
 45: #include <slepcblaslapack.h>

 49: /*
 50:    Duplicate and resize auxiliary basis
 51: */
 52: static PetscErrorCode PEPJDDuplicateBasis(PEP pep,BV *basis)
 53: {
 54:   PetscErrorCode     ierr;
 55:   PetscInt           nloc,m;
 56:   PetscMPIInt        rank,nproc;
 57:   BVType             type;
 58:   BVOrthogType       otype;
 59:   BVOrthogRefineType oref;
 60:   PetscReal          oeta;
 61:   BVOrthogBlockType  oblock;

 64:   if (pep->nev>1) {
 65:     BVCreate(PetscObjectComm((PetscObject)pep),basis);
 66:     MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rank);
 67:     MPI_Comm_size(PetscObjectComm((PetscObject)pep),&nproc);
 68:     BVGetSizes(pep->V,&nloc,NULL,&m);
 69:     if (rank==nproc-1) nloc += pep->nev-1;
 70:     BVSetSizes(*basis,nloc,PETSC_DECIDE,m);
 71:     BVGetType(pep->V,&type);
 72:     BVSetType(*basis,type);
 73:     BVGetOrthogonalization(pep->V,&otype,&oref,&oeta,&oblock);
 74:     BVSetOrthogonalization(*basis,otype,oref,oeta,oblock);
 75:     PetscObjectStateIncrease((PetscObject)*basis);
 76:   } else {
 77:     BVDuplicate(pep->V,basis);
 78:   }
 79:   return(0);
 80: }

 84: PetscErrorCode PEPSetUp_JD(PEP pep)
 85: {
 87:   PEP_JD         *pjd = (PEP_JD*)pep->data;
 88:   PetscBool      isshift,flg;
 89:   PetscInt       i;

 92:   pep->lineariz = PETSC_FALSE;
 93:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
 94:   if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
 95:   if (!pep->which) pep->which = PEP_LARGEST_MAGNITUDE;

 97:   /* Set STSHIFT as the default ST */
 98:   if (!((PetscObject)pep->st)->type_name) {
 99:     STSetType(pep->st,STSHIFT);
100:   }
101:   PetscObjectTypeCompare((PetscObject)pep->st,STSHIFT,&isshift);
102:   if (!isshift) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"JD only works with shift spectral transformation");

104:   if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
105:   STGetTransform(pep->st,&flg);
106:   if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");

108:   if (!pjd->keep) pjd->keep = 0.5;

110:   PEPAllocateSolution(pep,0);
111:   PEPSetWorkVecs(pep,5);
112:   PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
113:   for (i=0;i<pep->nmat;i++) {
114:     PEPJDDuplicateBasis(pep,pjd->TV+i);
115:   }
116:   PEPJDDuplicateBasis(pep,&pjd->W);
117:   if (pep->nev>1) {
118:     PEPJDDuplicateBasis(pep,&pjd->V);
119:     for (i=0;i<pep->nmat;i++) {
120:       BVDuplicateResize(pep->V,pep->nev-1,pjd->AX+i);
121:     }
122:     BVDuplicateResize(pep->V,pep->nev,&pjd->X);
123:     PetscCalloc3((pep->nev)*(pep->nev),&pjd->XpX,pep->nev*pep->nev,&pjd->T,pep->nev*pep->nev*pep->nmat,&pjd->Tj);
124:   } else pjd->V = pep->V;
125:   DSSetType(pep->ds,DSPEP);
126:   DSPEPSetDegree(pep->ds,pep->nmat-1);
127:   DSAllocate(pep->ds,pep->ncv);
128:   return(0);
129: }

133: /*
134:    Updates columns (low to (high-1)) of TV[i]
135: */
136: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
137: {
139:   PEP_JD         *pjd = (PEP_JD*)pep->data;
140:   PetscInt       pp,col,i,j,nloc,nconv,deg=pep->nmat-1;
141:   Vec            v1,v2,t1,t2;
142:   PetscScalar    *array1,*array2,*x2,*tt,*xx,*y2,zero=0.0,sone=1.0;
143:   PetscMPIInt    rk,np,count;
144:   PetscBLASInt   n,ld,one=1;

147:   nconv = pjd->nconv;
148:   PetscMalloc3(nconv,&tt,nconv,&x2,nconv,&xx);
149:   MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
150:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
151:   BVGetSizes(pep->V,&nloc,NULL,NULL); 
152:   t1 = w[0];
153:   t2 = w[1];
154:   for (col=low;col<high;col++) {
155:     BVGetColumn(pjd->V,col,&v1);
156:     VecGetArray(v1,&array1);
157:     if (nconv>0) {
158:       if (rk==np-1) { for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]; }
159:       PetscMPIIntCast(nconv,&count);
160:       MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
161:     }
162:     VecPlaceArray(t1,array1);
163:     for (pp=0;pp<pep->nmat;pp++) {
164:       BVGetColumn(pjd->TV[pp],col,&v2);
165:       VecGetArray(v2,&array2);
166:       VecPlaceArray(t2,array2);
167:       MatMult(pep->A[pp],t1,t2);
168:       if (nconv) {
169:         PetscBLASIntCast(pjd->nconv,&n);
170:         PetscBLASIntCast(pep->nev,&ld);
171:         for (j=0;j<nconv;j++) tt[j] = x2[j];
172:         for (i=pp+1;i<pep->nmat;i++) {
173:           BVMultVec(pjd->AX[i],1.0,1.0,t2,tt);
174:           if (i!=pep->nmat-1) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,tt,&one));
175:         }
176:         BVDotVec(pjd->X,t1,xx);
177:         if (rk==np-1 && pp<deg) {
178:           y2 = array2+nloc;
179:           for (j=0;j<nconv;j++) { y2[j] = xx[j]; xx[j] = x2[j]; }
180:           PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*pp,&ld,y2,&one));
181:           for (i=pp+1;i<pep->nmat-1;i++) {
182:             PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,xx,&one,&zero,tt,&one));
183:             PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*i,&ld,tt,&one));
184:             for (j=0;j<nconv;j++) y2[j] += tt[j];
185:             if (i<pep->nmat-2) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,xx,&one));
186:           }
187:         }
188:       }
189:       VecResetArray(t2);
190:       VecRestoreArray(v2,&array2);
191:       BVRestoreColumn(pjd->TV[pp],col,&v2);        
192:     }
193:     VecResetArray(t1);
194:     VecRestoreArray(v1,&array1);
195:     BVRestoreColumn(pjd->V,col,&v1);
196:   }
197:   PetscFree3(tt,x2,xx);
198:   return(0);
199: }

203: /*
204:    RRQR of X. Xin*P=Xou*R. Rank of R is rk
205: */
206: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
207: {
208: #if defined(SLEPC_MISSING_LAPACK_GEQP3)
210:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQP3 - Lapack routine is unavailable");
211: #else
213:   PetscInt       i,j,n,r;
214:   PetscBLASInt   row_,col_,ldx_,*p,lwork,info,n_;
215:   PetscScalar    *tau,*work;
216:   PetscReal      tol,*rwork;

219:   PetscBLASIntCast(row,&row_);
220:   PetscBLASIntCast(col,&col_);
221:   PetscBLASIntCast(ldx,&ldx_);
222:   n = PetscMin(row,col);
223:   PetscBLASIntCast(n,&n_);
224:   lwork = 3*col_+1;
225:   PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
226:   for (i=1;i<col;i++) p[i] = 0;
227:   p[0] = 1;

229:   /* rank revealing QR */
230: #if defined(PETSC_USE_COMPLEX)
231:   PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
232: #else
233:   PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
234: #endif
235:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQP3 %d",info);
236:   if (P) for (i=0;i<col;i++) P[i] = p[i];

238:   /* rank computation */
239:   tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
240:   r = 1;
241:   for (i=1;i<n;i++) { 
242:     if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
243:     else break;
244:   }
245:   if (rk) *rk=r;

247:   /* copy upper triangular matrix if requested */
248:   if (R) {
249:      for (i=0;i<r;i++) {
250:        PetscMemzero(R+i*ldr,r*sizeof(PetscScalar));
251:        for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
252:      }
253:   }
254:   PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
255:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
256:   PetscFree4(p,tau,work,rwork);
257:   return(0);
258: #endif
259: }

263: /*
264:    Application of extended preconditioner
265: */
266: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
267: {
268:   PetscInt          i,j,nloc,n,ld;
269:   PetscMPIInt       rk,np,count;
270:   Vec               tx,ty;
271:   PEP_JD_PCSHELL    *ctx;
272:   PetscErrorCode    ierr;
273:   const PetscScalar *array1;
274:   PetscScalar       *x2=NULL,*t=NULL,*ps,*array2;
275:   PetscBLASInt      one=1.0,ld_,n_;

278:   PCShellGetContext(pc,(void**)&ctx);
279:   n  = ctx->n;
280:   ps = ctx->ps;
281:   ld = ctx->ld;
282:   if (n) {
283:     PetscMalloc2(n,&x2,n,&t);
284:     MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rk);
285:     MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
286:     if (rk==np-1) {
287:       VecGetSize(ctx->work[0],&nloc); 
288:       VecGetArrayRead(x,&array1);
289:       for (i=0;i<n;i++) x2[i] = array1[nloc+i];
290:       VecRestoreArrayRead(x,&array1);
291:     }
292:     PetscMPIIntCast(n,&count);
293:     MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pc));
294:   }

296:   /* y = B\x apply PC */
297:   tx = ctx->work[0];
298:   ty = ctx->work[1];
299:   VecGetArrayRead(x,&array1);
300:   VecPlaceArray(tx,array1);
301:   VecGetArray(y,&array2);
302:   VecPlaceArray(ty,array2);
303:   PCApply(ctx->pc,tx,ty);
304:   if (n) {
305:     for (j=0;j<n;j++) {
306:       t[j] = 0.0;
307:       for (i=0;i<n;i++) t[j] += ctx->M[i+j*ld]*x2[i];
308:     }
309:     if (rk==np-1) for (i=0;i<n;i++) array2[nloc+i] = t[i];
310:     PetscBLASIntCast(ld,&ld_);
311:     PetscBLASIntCast(n,&n_);
312:     PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n_,ps,&ld_,t,&one));
313:     BVMultVec(ctx->X,-1.0,1.0,ty,t);
314:     PetscFree2(x2,t);
315:   }
316:   VecResetArray(tx);
317:   VecResetArray(ty);
318:   VecRestoreArrayRead(x,&array1);
319:   VecRestoreArray(y,&array2);
320:   return(0);
321: }

325: /*
326:    Application of shell preconditioner:
327:       y = B\x - eta*B\p,  with eta = (u'*B\x)/(u'*B\p)
328: */
329: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
330: {
332:   PetscScalar    eta;
333:   PEP_JD_PCSHELL *ctx;

336:   PCShellGetContext(pc,(void**)&ctx);

338:   /* y = B\x apply extended PC */
339:   PEPJDExtendedPCApply(pc,x,y);

341:   /* Compute eta = u'*y / u'*Bp */
342:   VecDot(y,ctx->u,&eta);
343:   eta /= ctx->gamma;
344:   
345:   /* y = y - eta*Bp */
346:   VecAXPY(y,-eta,ctx->Bp); 
347:   return(0);
348: }

352: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
353: {
355:   PetscMPIInt    np,rk,count;
356:   PetscScalar    *array1,*array2;
357:   PetscInt       nloc;

360:   MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
361:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
362:   BVGetSizes(pep->V,&nloc,NULL,NULL);
363:   if (v) {
364:     VecGetArray(v,&array1);
365:     VecGetArray(vex,&array2);
366:     if (back) {
367:       PetscMemcpy(array1,array2,nloc*sizeof(PetscScalar));
368:     } else {
369:       PetscMemcpy(array2,array1,nloc*sizeof(PetscScalar));
370:     }
371:     VecRestoreArray(v,&array1);
372:     VecRestoreArray(vex,&array2);
373:   }
374:   if (a) {
375:     if (rk==np-1) {
376:       VecGetArray(vex,&array2);
377:       if (back) {
378:         PetscMemcpy(a,array2+nloc+off,na*sizeof(PetscScalar));
379:       } else {
380:         PetscMemcpy(array2+nloc+off,a,na*sizeof(PetscScalar));
381:       }
382:       VecRestoreArray(vex,&array2);
383:     }
384:     if (back) {
385:       PetscMPIIntCast(na,&count);
386:       MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
387:     }
388:   }
389:   return(0);
390: }

394: static PetscErrorCode PEPJDComputePResidual(PEP pep,Vec u,PetscScalar theta,Vec p,Vec *work)
395: {
396:   PEP_JD         *pjd = (PEP_JD*)pep->data;
398:   PetscMPIInt    rk,np,count;
399:   Vec            tu,tp,w;
400:   PetscScalar    *array1,*array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,sone=1.0,zero=0.0;
401:   PetscInt       i,j,nconv=pjd->nconv,nloc,deg=pep->nmat-1;
402:   PetscBLASInt   n,ld,one=1;

405:   if (nconv>0) {
406:     PetscMalloc4(nconv,&xx,nconv,&tt,nconv,&x2,nconv,&q);
407:     MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
408:     MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
409:     if (rk==np-1) {
410:       BVGetSizes(pep->V,&nloc,NULL,NULL); 
411:       VecGetArray(u,&array1);
412:       for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
413:       VecRestoreArray(u,&array1);
414:     }
415:     PetscMPIIntCast(nconv,&count);
416:     MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
417:   }
418:   tu = work[0];
419:   tp = work[1];
420:   w  = work[2];
421:   VecGetArray(u,&array1);
422:   VecPlaceArray(tu,array1);
423:   VecGetArray(p,&array2);
424:   VecPlaceArray(tp,array2);
425:   VecSet(tp,0.0);
426:   for (i=1;i<pep->nmat;i++) {
427:     MatMult(pep->A[i],tu,w);
428:     VecAXPY(tp,fact*(PetscReal)i,w);
429:     fact *= theta;
430:   }
431:   if (nconv) {
432:     PetscBLASIntCast(nconv,&n);
433:     PetscBLASIntCast(pep->nev,&ld);
434:     for (j=0;j<nconv;j++) q[j] = x2[j];
435:     fact = theta;
436:     for (i=2;i<pep->nmat;i++) {
437:       BVMultVec(pjd->AX[i],1.0,1.0,tp,q);
438:       PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
439:       for (j=0;j<nconv;j++) q[j] += (PetscReal)i*fact*x2[j];
440:       fact *= theta;
441:     }
442:     BVSetActiveColumns(pjd->X,0,nconv);
443:     BVDotVec(pjd->X,tu,xx);
444:     if (rk==np-1) {
445:       y2 = array2+nloc;
446:       for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
447:       PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld,&ld,y2,&one));
448:       fact = theta;
449:       for (j=2;j<deg;j++) {
450:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
451:         for (i=0;i<nconv;i++) tt[i] += (PetscReal)j*fact*xx[i];
452:         PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
453:         for (i=0;i<nconv;i++) y2[i] += tt[i];
454:         PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
455:         for (i=0;i<nconv;i++) q[i] += (PetscReal)j*fact*x2[i];
456:         fact *= theta;
457:       }
458:     }
459:     PetscFree4(xx,x2,q,tt);
460:   }
461:   VecResetArray(tu);
462:   VecRestoreArray(u,&array1);
463:   VecResetArray(tp);
464:   VecRestoreArray(p,&array2);
465:   return(0);
466: }

470: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
471: {
472:   PEP_JD         *pjd = (PEP_JD*)pep->data;
474:   PetscScalar    *tt;
475:   Vec            vg,wg;
476:   PetscInt       i;
477:   PetscReal      norm;

480:   PetscMalloc1(pep->nev-1,&tt);
481:   if (pep->nini==0) {
482:     BVSetRandomColumn(pjd->V,0);
483:     for (i=0;i<pep->nev-1;i++) tt[i] = 0.0;
484:     BVGetColumn(pjd->V,0,&vg);
485:     PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,vg,PETSC_FALSE);
486:     BVRestoreColumn(pjd->V,0,&vg);
487:     BVNormColumn(pjd->V,0,NORM_2,&norm);
488:     BVScaleColumn(pjd->V,0,1.0/norm);
489:     BVGetColumn(pjd->V,0,&vg);
490:     BVGetColumn(pjd->W,0,&wg);
491:     VecSet(wg,0.0);
492:     PEPJDComputePResidual(pep,vg,pep->target,wg,w);
493:     BVRestoreColumn(pjd->W,0,&wg);
494:     BVRestoreColumn(pjd->V,0,&vg);
495:     BVNormColumn(pjd->W,0,NORM_2,&norm);
496:     BVScaleColumn(pjd->W,0,1.0/norm);
497:   } else {
498:    SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TO DO");
499:   }
500:   PetscFree(tt);
501:   return(0);
502: }

506: static PetscErrorCode PEPJDShellMatMult(Mat P,Vec x,Vec y)
507: {
508:   PetscErrorCode    ierr;
509:   PEP_JD_MATSHELL   *matctx;
510:   PEP_JD            *pjd;
511:   PetscMPIInt       rk,np,count;
512:   PetscInt          i,j,nconv,nloc,nmat,ldt,deg;
513:   Vec               tx,ty;
514:   PetscScalar       *array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,theta,*yy=NULL,sone=1.0,zero=0.0;
515:   PetscBLASInt      n,ld,one=1;
516:   const PetscScalar *array1;

519:   MatShellGetContext(P,(void**)&matctx);
520:   pjd   = (PEP_JD*)(matctx->pep->data);
521:   nconv = pjd->nconv;
522:   theta = matctx->theta;
523:   nmat  = matctx->pep->nmat;
524:   deg   = nmat-1;
525:   ldt   = matctx->pep->nev;
526:   if (nconv>0) {
527:     PetscMalloc5(nconv,&tt,nconv,&x2,nconv,&q,nconv,&xx,nconv,&yy);
528:     MPI_Comm_rank(PetscObjectComm((PetscObject)P),&rk);
529:     MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
530:     if (rk==np-1) {
531:       BVGetSizes(matctx->pep->V,&nloc,NULL,NULL); 
532:       VecGetArrayRead(x,&array1);
533:       for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
534:       VecRestoreArrayRead(x,&array1);
535:     }
536:     PetscMPIIntCast(nconv,&count);
537:     MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)P));
538:   }
539:   tx = matctx->work[0];
540:   ty = matctx->work[1];
541:   VecGetArrayRead(x,&array1);
542:   VecPlaceArray(tx,array1);
543:   VecGetArray(y,&array2);
544:   VecPlaceArray(ty,array2);
545:   VecSet(ty,0.0);
546:   MatMult(matctx->P,tx,ty);
547:   if (nconv) {
548:     PetscBLASIntCast(pjd->nconv,&n);
549:     PetscBLASIntCast(ldt,&ld);
550:     for (j=0;j<nconv;j++) q[j] = x2[j];
551:     fact = theta;
552:     for (i=1;i<nmat;i++) {
553:       BVMultVec(pjd->AX[i],1.0,1.0,ty,q);
554:       PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
555:       for (j=0;j<nconv;j++) q[j] += fact*x2[j];
556:       fact *= theta;
557:     }
558:     BVSetActiveColumns(pjd->X,0,nconv);    
559:     BVDotVec(pjd->X,tx,xx);
560:     if (rk==np-1) {
561:       y2 = array2+nloc;
562:       for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
563:       fact = theta;
564:       for (j=1;j<deg;j++) {
565:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
566:         for (i=0;i<nconv;i++) tt[i] += fact*xx[i];
567:         PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
568:         for (i=0;i<nconv;i++) y2[i] += tt[i];
569:         PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
570:         for (i=0;i<nconv;i++) q[i] += fact*x2[i];
571:         fact *= theta;
572:       }
573:     }
574:     PetscFree5(tt,x2,q,xx,yy);
575:   }
576:   VecResetArray(tx);
577:   VecRestoreArrayRead(x,&array1);
578:   VecResetArray(ty);
579:   VecRestoreArray(y,&array2);
580:   return(0);
581: }

585: static PetscErrorCode PEPJDCreateShellPC(PEP pep)
586: {
587:   PEP_JD          *pjd = (PEP_JD*)pep->data;
588:   PEP_JD_PCSHELL  *pcctx;
589:   PEP_JD_MATSHELL *matctx;
590:   KSP             ksp;
591:   PetscInt        nloc,mloc;
592:   PetscMPIInt     np,rk;
593:   PetscErrorCode  ierr;

596:   PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
597:   PCSetType(pjd->pcshell,PCSHELL);
598:   PCShellSetName(pjd->pcshell,"PCPEPJD");
599:   PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
600:   PetscNew(&pcctx);
601:   PCShellSetContext(pjd->pcshell,pcctx);
602:   STGetKSP(pep->st,&ksp);
603:   BVCreateVec(pjd->V,&pcctx->Bp);
604:   KSPGetPC(ksp,&pcctx->pc);
605:   PetscObjectReference((PetscObject)pcctx->pc);
606:   MatGetLocalSize(pep->A[0],&mloc,&nloc);
607:   if (pep->nev>1) {
608:     MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
609:     MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
610:     if (rk==np-1) { nloc += pep->nev-1; mloc += pep->nev-1; }
611:   }
612:   PetscNew(&matctx);
613:   MatCreateShell(PetscObjectComm((PetscObject)pep),nloc,mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
614:   MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)())PEPJDShellMatMult);
615:   matctx->pep = pep;
616:   MatDuplicate(pep->A[0],MAT_DO_NOT_COPY_VALUES,&matctx->P);
617:   PCSetOperators(pcctx->pc,matctx->P,matctx->P);
618:   KSPSetPC(ksp,pjd->pcshell);
619:   KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
620:   if (pep->nev>1) {
621:     PetscMalloc2(pep->nev*pep->nev,&pcctx->M,pep->nev*pep->nev,&pcctx->ps);
622:     pcctx->X  = pjd->X;
623:     pcctx->ld = pep->nev;
624:   }
625:   return(0);
626: }

630: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
631: {
633:   PEP_JD         *pjd = (PEP_JD*)pep->data;
634:   PEP_JD_PCSHELL *pcctx;  
635:   PetscInt       i,j,k,n=pjd->nconv,ld=pep->nev,deg=pep->nmat-1;
636:   PetscScalar    fact,*M,*ps,*work,*U,*V,*S,sone=1.0,zero=0.0;
637:   PetscReal      tol,maxeig=0.0,*sg,*rwork;
638:   PetscBLASInt   n_,info,ld_,*p,lw_,rk=0;

641: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GETRI) || defined(PETSC_MISSING_LAPACK_GETRF)
643:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GETRI/GETRF - Lapack routine is unavailable");
644: #else
645:   if (n) { 
646:     PCShellGetContext(pjd->pcshell,(void**)&pcctx);
647:     pcctx->n = n;
648:     M  = pcctx->M;
649:     ps = pcctx->ps;
650:                       /* h, and q are vectors containing diagonal matrices */
651:     PetscCalloc7(n*n,&U,n*n,&V,n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p);
652:     /* pseudo-inverse */
653:     for (j=0;j<n;j++) {
654:       for (i=0;i<j;i++) S[n*j+i] = -pjd->T[pep->nev*j+i];
655:       S[n*i+i] = theta-pjd->T[pep->nev*i+i];
656:     }
657:     PetscBLASIntCast(n,&n_);
658:     PetscBLASIntCast(ld,&ld_);
659:     lw_ = 10*n_;
660: #if !defined (PETSC_USE_COMPLEX)
661:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
662: #else
663:     PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
664: #endif
665:     for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
666:     tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
667:     for (j=0;j<n;j++) {
668:       if (sg[j]>tol) {
669:         for (i=0;i<n;i++) U[j*n+i] /= sg[j];
670:         rk++;
671:       } else break;
672:     }
673:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));

675:     /* compute M */
676:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
677:     fact = theta;
678:     PetscMemzero(S,n*n*sizeof(PetscScalar));
679:     for (j=0;j<n;j++) S[j*(n+1)] = 1.0; /* q=S */
680:     for (k=0;k<deg;k++) {
681:       for (j=0;j<n;j++) for (i=0;i<n;i++) V[j*n+i] = S[j*n+i] + M[j*ld+i]*fact;
682:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
683:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
684:       PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n_,&n_,&sone,pjd->T,&ld_,S,&n_));
685:       for (j=0;j<n;j++) S[j*(n+1)] += fact;
686:       fact *=theta;
687:     }
688:     /* inverse */
689:     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
690:     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
691:     PetscFree7(U,V,S,sg,work,rwork,p);
692:   }
693:   return(0);
694: #endif
695: }

699: static PetscErrorCode PEPJDPCMatSetUp(PEP pep,PetscScalar theta)
700: {
701:   PetscErrorCode  ierr;
702:   PEP_JD          *pjd = (PEP_JD*)pep->data;
703:   PEP_JD_MATSHELL *matctx;
704:   PEP_JD_PCSHELL  *pcctx;  
705:   MatStructure    str;
706:   PetscScalar     t;
707:   PetscInt        i;

710:   MatShellGetContext(pjd->Pshell,(void**)&matctx);
711:   PCShellGetContext(pjd->pcshell,(void**)&pcctx);
712:   STGetMatStructure(pep->st,&str);
713:   MatCopy(pep->A[0],matctx->P,str);
714:   t = theta;
715:   for (i=1;i<pep->nmat;i++) {
716:     if (t!=0.0) { MatAXPY(matctx->P,t,pep->A[i],str); }
717:     t *= theta;
718:   }
719:   PCSetOperators(pcctx->pc,matctx->P,matctx->P);
720:   PCSetUp(pcctx->pc);
721:   matctx->theta = theta;
722:   return(0);
723: }

727: static PetscErrorCode PEPJDEigenvectors(PEP pep)
728: {
730:   PEP_JD         *pjd = (PEP_JD*)pep->data;
731:   PetscBLASInt   ld,nconv,info,nc;
732:   PetscScalar    *Z,*w;
733:   PetscReal      *wr,norm;
734:   PetscInt       i;
735:   Mat            U;
736:  
738: #if defined(SLEPC_MISSING_LAPACK_TREVC)
740:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
741: #else
742:   PetscMalloc3(pjd->nconv*pjd->nconv,&Z,3*pep->nev,&wr,2*pep->nev,&w);
743:   PetscBLASIntCast(pep->nev,&ld);
744:   PetscBLASIntCast(pjd->nconv,&nconv);
745: #if !defined(PETSC_USE_COMPLEX)
746:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
747: #else
748:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
749: #endif
750:   MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
751:   BVSetActiveColumns(pjd->X,0,pjd->nconv);
752:   BVMultInPlace(pjd->X,U,0,pjd->nconv);
753:   for (i=0;i<pjd->nconv;i++) {
754:     BVNormColumn(pjd->X,i,NORM_2,&norm);
755:     BVScaleColumn(pjd->X,i,1.0/norm);  
756:   }
757:   MatDestroy(&U);
758:   PetscFree3(Z,wr,w);
759:   return(0);
760: #endif
761: }

765: PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv,Vec u,Vec *ww)
766: {
767:   PetscErrorCode    ierr;
768:   PEP_JD            *pjd = (PEP_JD*)pep->data;
769:   PetscInt          j,i,ldds,rk,*P,nvv=*nv;
770:   Vec               v;
771:   PetscBLASInt      n,ld,rk_,nv_,info,one=1;
772:   PetscScalar       sone=1.0,*Tj,*R,*r,*tt,*pX;
773:   Mat               X;
774:   const PetscScalar *array;

777: #if defined(SLEPC_MISSING_LAPACK_TRTRI)
779:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TRTRI - Lapack routine is unavailable");
780: #else
781:   /* update AX and XpX */
782:   VecGetArrayRead(u,&array);
783:   VecPlaceArray(ww[0],array);
784:   for (j=0;j<pep->nmat;j++) {
785:     BVGetColumn(pjd->AX[j],pjd->nconv-1,&v);
786:     MatMult(pep->A[j],ww[0],v);
787:     BVRestoreColumn(pjd->AX[j],pjd->nconv-1,&v);
788:     BVSetActiveColumns(pjd->AX[j],0,pjd->nconv);
789:   }
790:   BVDotVec(pjd->X,ww[0],pjd->XpX+(pjd->nconv-1)*(pep->nev));
791:   for (j=0;j<pjd->nconv-1;j++) pjd->XpX[j*(pep->nev)+pjd->nconv-1] = PetscConj(pjd->XpX[(pjd->nconv-1)*(pep->nev)+j]);
792:   VecResetArray(ww[0]);
793:   VecRestoreArrayRead(u,&array);
794:   
795:   /* Compute powers of T */
796:   PetscBLASIntCast(pjd->nconv,&n);
797:   PetscBLASIntCast(pep->nev,&ld);
798:   PetscMemzero(pjd->Tj,pep->nev*pep->nev*pep->nmat*sizeof(PetscScalar));
799:   Tj = pjd->Tj;
800:   for (j=0;j<pep->nmat;j++) Tj[(pep->nev+1)*j] = 1.0;
801:   Tj = pjd->Tj+pep->nev*pep->nev;
802:   PetscMemcpy(Tj,pjd->T,pep->nev*pjd->nconv*sizeof(PetscScalar));
803:   for (j=2;j<pep->nmat;j++) {
804:     PetscMemcpy(Tj+pep->nev*pep->nev,Tj,pep->nev*pjd->nconv*sizeof(PetscScalar));
805:     Tj += pep->nev*pep->nev;
806:     PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n,&n,&sone,pjd->T,&ld,Tj,&ld));
807:   }

809:   /* Extend search space */
810:   PetscCalloc4(nvv,&P,nvv*nvv,&R,nvv,&r,pep->nev-1,&tt);
811:   DSGetLeadingDimension(pep->ds,&ldds);
812:   DSGetArray(pep->ds,DS_MAT_X,&pX);
813:   PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
814:   for (i=0;i<rk-1;i++) r[i] = PetscConj(R[nvv*i]*pep->eigr[P[i+1]]); /* first row scaled with permuted diagonal */
815:   PetscBLASIntCast(rk,&rk_);
816:   PetscBLASIntCast(nvv,&nv_);
817:   PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
818:   if (info) SETERRQ1(PETSC_COMM_SELF,1,"Error in xTRTRI, info=%D",(PetscInt)info);
819:   PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r,&one));
820:   for (i=0;i<rk;i++) r[i] = PetscConj(r[i]); /* revert */
821:   BVSetActiveColumns(pjd->V,0,nvv);
822:   for (j=0;j<rk-1;j++) {
823:     PetscMemcpy(R+j*nvv,pX+(j+1)*ldds,nvv*sizeof(PetscScalar));
824:   } 
825:   DSRestoreArray(pep->ds,DS_MAT_X,&pX);
826:   MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk-1,R,&X);
827:   BVMultInPlace(pjd->V,X,0,rk-1);
828:   MatDestroy(&X);
829:   BVSetActiveColumns(pjd->V,0,rk-1);
830:   for (j=0;j<rk-1;j++) {
831:     BVGetColumn(pjd->V,j,&v);
832:     PEPJDCopyToExtendedVec(pep,NULL,r+j,1,pjd->nconv-1,v,PETSC_FALSE);
833:     BVRestoreColumn(pjd->V,j,&v);
834:   }
835:   BVOrthogonalize(pjd->V,NULL); 
836:   for (j=0;j<rk-1;j++) {
837:     BVGetColumn(pjd->W,j,&v);
838:     PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,v,PETSC_FALSE);
839:     BVRestoreColumn(pjd->W,j,&v);
840:   }
841:   *nv = rk-1;
842:   PetscFree4(P,R,r,tt);
843: #endif
844:   return(0);
845: }

849: PetscErrorCode PEPSolve_JD(PEP pep)
850: {
851:   PetscErrorCode  ierr;
852:   PEP_JD          *pjd = (PEP_JD*)pep->data;
853:   PetscInt        k,nv,ld,minv,low,high,*P,dim;
854:   PetscScalar     theta=0.0,*pX,*stt,*exu,*exr,*exp,*R,*eig;
855:   PetscReal       norm,*res;
856:   PetscBool       lindep,initial=PETSC_FALSE,flglk=PETSC_FALSE,flgre=PETSC_FALSE;
857:   Vec             t,u,p,r,*ww=pep->work,v;
858:   Mat             G,X,Y;
859:   KSP             ksp;
860:   PEP_JD_PCSHELL  *pcctx;
861:   PEP_JD_MATSHELL *matctx;

864:   DSGetLeadingDimension(pep->ds,&ld);
865:   PetscMalloc5(ld,&P,ld,&stt,pep->nev-1,&exu,pep->nev-1,&exr,pep->nev-1,&exp);
866:   PetscMalloc3(ld*ld,&R,pep->ncv,&eig,pep->ncv,&res);
867:   BVCreateVec(pjd->V,&u);
868:   VecDuplicate(u,&p);
869:   VecDuplicate(u,&r);
870:   STGetKSP(pep->st,&ksp);

872:   if (pep->nini) {
873:     nv = pep->nini; initial = PETSC_TRUE;
874:   } else {
875:     theta = pep->target;
876:     nv = 1;
877:   }
878:   PEPJDProcessInitialSpace(pep,ww);
879:   BVCopyVec(pjd->V,0,u);

881:   /* Restart loop */
882:   while (pep->reason == PEP_CONVERGED_ITERATING) {
883:     pep->its++;

885:     low = (flglk || flgre)? 0: nv-1;
886:     high = nv;
887:     DSSetDimensions(pep->ds,nv,0,0,0);
888:     BVSetActiveColumns(pjd->V,low,high);
889:     PEPJDUpdateTV(pep,low,high,ww);
890:     BVSetActiveColumns(pjd->W,low,high);
891:     for (k=0;k<pep->nmat;k++) {
892:       BVSetActiveColumns(pjd->TV[k],low,high);
893:       DSGetMat(pep->ds,DSMatExtra[k],&G);
894:       BVMatProject(pjd->TV[k],NULL,pjd->W,G);
895:       DSRestoreMat(pep->ds,DSMatExtra[k],&G);
896:     }
897:     BVSetActiveColumns(pjd->V,0,nv);
898:     BVSetActiveColumns(pjd->W,0,nv);

900:     /* Solve projected problem */
901:     if (nv>1 || initial ) {
902:       DSSetState(pep->ds,DS_STATE_RAW);
903:       DSSolve(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv);
904:       DSSort(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv,NULL,NULL,NULL);
905:       DSSort(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv,NULL,NULL,NULL);
906:       theta = pep->eigr[0];
907: #if !defined(PETSC_USE_COMPLEX)
908:       if (PetscAbsScalar(pep->eigi[pep->nconv])!=0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PJD solver not implemented for complex Ritz values in real arithmetic");
909: #endif

911:       /* Compute Ritz vector u=V*X(:,1) */
912:       DSGetArray(pep->ds,DS_MAT_X,&pX);
913:       BVSetActiveColumns(pjd->V,0,nv);
914:       BVMultVec(pjd->V,1.0,0.0,u,pX);
915:       DSRestoreArray(pep->ds,DS_MAT_X,&pX);
916:     }
917:     PEPJDUpdateExtendedPC(pep,theta);

919:     /* Replace preconditioner with one containing projectors */
920:     if (!pjd->pcshell) {
921:       PEPJDCreateShellPC(pep);
922:       PCShellGetContext(pjd->pcshell,(void**)&pcctx);
923:       MatShellGetContext(pjd->Pshell,(void**)&matctx);
924:       matctx->work = ww;
925:       pcctx->work  = ww;
926:     }
927:     PEPJDPCMatSetUp(pep,theta);
928:     
929:     /* Compute r and r' */
930:     MatMult(pjd->Pshell,u,r);
931:     PEPJDComputePResidual(pep,u,theta,p,ww);
932:     pcctx->u = u;

934:     /* Check convergence */
935:     VecNorm(r,NORM_2,&norm);
936:     (*pep->converged)(pep,theta,0,norm,&pep->errest[pep->nconv],pep->convergedctx);
937:     if (pep->its >= pep->max_it) pep->reason = PEP_DIVERGED_ITS;

939:     if (pep->errest[pep->nconv]<pep->tol) {

941:       /* Ritz pair converged */
942:       minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
943:       if (pep->nev>1) {
944:         BVGetColumn(pjd->X,pjd->nconv,&v);
945:         PEPJDCopyToExtendedVec(pep,v,pjd->T+pep->nev*pjd->nconv,pep->nev-1,0,u,PETSC_TRUE);
946:         BVRestoreColumn(pjd->X,pjd->nconv,&v);
947:         BVSetActiveColumns(pjd->X,0,pjd->nconv+1);
948:         BVNormColumn(pjd->X,pjd->nconv,NORM_2,&norm);
949:         BVScaleColumn(pjd->X,pjd->nconv,1.0/norm);
950:         for (k=0;k<pjd->nconv;k++) pjd->T[pep->nev*pjd->nconv+k] /= norm;
951:         pjd->T[(pep->nev+1)*pjd->nconv] = pep->eigr[0];
952:       } else {
953:         BVInsertVec(pep->V,pep->nconv,u);
954:       }
955:       pjd->nconv++;
956:       if (pjd->nconv >= pep->nev) pep->reason = PEP_CONVERGED_TOL;

958:       if (pep->reason==PEP_CONVERGED_ITERATING) {
959:         PEPJDLockConverged(pep,&nv,u,ww);
960:         BVCopyVec(pjd->V,nv,u);
961:         if (nv==1) theta = pep->target;
962:       }
963:       flglk = PETSC_TRUE;
964:     } else if (nv==pep->ncv-1) {

966:       /* Basis full, force restart */
967:       minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
968:       DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
969:       DSGetArray(pep->ds,DS_MAT_X,&pX);
970:       PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
971:       DSRestoreArray(pep->ds,DS_MAT_X,&pX);
972:       DSGetArray(pep->ds,DS_MAT_Y,&pX);
973:       PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
974:       DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
975:       DSGetMat(pep->ds,DS_MAT_X,&X);
976:       BVMultInPlace(pjd->V,X,pep->nconv,minv);
977:       DSRestoreMat(pep->ds,DS_MAT_X,&X);
978:       DSGetMat(pep->ds,DS_MAT_Y,&Y);
979:       BVMultInPlace(pjd->W,Y,pep->nconv,minv);
980:       DSRestoreMat(pep->ds,DS_MAT_Y,&Y);
981:       nv = minv;
982:       flgre = PETSC_TRUE;
983:     } else {
984:       /* Solve correction equation to expand basis */
985:       PEPJDExtendedPCApply(pjd->pcshell,p,pcctx->Bp);
986:       VecDot(pcctx->Bp,u,&pcctx->gamma);
987:       BVGetColumn(pjd->V,nv,&t);
988:       KSPSolve(ksp,r,t);
989:       BVRestoreColumn(pjd->V,nv,&t);
990:       BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
991:       if (lindep) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
992:       BVScaleColumn(pjd->V,nv,1.0/norm);
993:       BVInsertVec(pjd->W,nv,r);
994:       BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
995:       if (lindep) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
996:       BVScaleColumn(pjd->W,nv,1.0/norm);
997:       nv++;
998:       flglk = PETSC_FALSE;
999:       flgre = PETSC_FALSE;
1000:     }
1001:     for (k=pjd->nconv;k<nv;k++) {
1002:       eig[k] = pep->eigr[k-pjd->nconv];
1003:       res[k] = pep->errest[k-pjd->nconv];
1004:     }
1005:     PEPMonitor(pep,pep->its,pjd->nconv,eig,pep->eigi,res,nv);
1006:   }
1007:   if (pep->nev>1) {
1008:     PEPJDEigenvectors(pep);
1009:     for (k=0;k<pjd->nconv;k++) {
1010:       BVGetColumn(pjd->X,k,&v);
1011:       BVInsertVec(pep->V,k,v);
1012:       BVRestoreColumn(pjd->X,k,&v);
1013:       pep->eigr[k] = pjd->T[(pep->nev+1)*k]; 
1014:     }
1015:     PetscFree2(pcctx->M,pcctx->ps); 
1016:   }
1017:   pep->nconv = pjd->nconv; 
1018:   KSPSetPC(ksp,pcctx->pc);
1019:   MatDestroy(&matctx->P);
1020:   VecDestroy(&pcctx->Bp);
1021:   MatDestroy(&pjd->Pshell);
1022:   PCDestroy(&pcctx->pc);
1023:   PetscFree(pcctx);
1024:   PetscFree(matctx);
1025:   PCDestroy(&pjd->pcshell);
1026:   PetscFree5(P,stt,exu,exr,exp);
1027:   PetscFree3(R,eig,res);
1028:   VecDestroy(&u);
1029:   VecDestroy(&r);
1030:   VecDestroy(&p);
1031:   return(0);
1032: }

1036: PetscErrorCode PEPReset_JD(PEP pep)
1037: {
1039:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1040:   PetscInt       i;

1043:   for (i=0;i<pep->nmat;i++) {
1044:     BVDestroy(pjd->TV+i);
1045:   }
1046:   BVDestroy(&pjd->W);
1047:   if (pep->nev>1) {
1048:     BVDestroy(&pjd->V);
1049:     for (i=0;i<pep->nmat;i++) {
1050:       BVDestroy(pjd->AX+i);
1051:     }
1052:     BVDestroy(&pjd->X);
1053:     PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
1054:   }
1055:   PetscFree2(pjd->TV,pjd->AX);
1056:   return(0);
1057: }

1061: PetscErrorCode PEPDestroy_JD(PEP pep)
1062: {

1066:   PetscFree(pep->data);
1067:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
1068:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
1069:   return(0);
1070: }

1074: PETSC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
1075: {
1076:   PEP_JD         *pjd;

1080:   PetscNewLog(pep,&pjd);
1081:   pep->data = (void*)pjd;

1083:   pjd->keep = 0;
1084:   pep->ops->solve          = PEPSolve_JD;
1085:   pep->ops->setup          = PEPSetUp_JD;
1086:   pep->ops->setfromoptions = PEPSetFromOptions_JD;
1087:   pep->ops->reset          = PEPReset_JD;
1088:   pep->ops->destroy        = PEPDestroy_JD;
1089:   pep->ops->view           = PEPView_JD;
1090:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
1091:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
1092:   return(0);
1093: }