Actual source code: pjd.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    SLEPc polynomial eigensolver: "jd"

 13:    Method: Jacobi-Davidson

 15:    Algorithm:

 17:        Jacobi-Davidson for polynomial eigenvalue problems.

 19:    References:

 21:        [1] C. Campos and J.E. Roman, "A polynomial Jacobi-Davidson solver
 22:            with support for non-monomial bases and deflation", BIT Numer.
 23:            Math. 60:295-318, 2020.

 25:        [2] G.L.G. Sleijpen et al., "Jacobi-Davidson type methods for
 26:            generalized eigenproblems and polynomial eigenproblems", BIT
 27:            36(3):595-633, 1996.

 29:        [3] Feng-Nan Hwang, Zih-Hao Wei, Tsung-Ming Huang, Weichung Wang,
 30:            "A Parallel Additive Schwarz Preconditioned Jacobi-Davidson
 31:            Algorithm for Polynomial Eigenvalue Problems in Quantum Dot
 32:            Simulation", J. Comput. Phys. 229(8):2932-2947, 2010.
 33: */

 35: #include <slepc/private/pepimpl.h>
 36: #include <slepcblaslapack.h>

 38: static PetscBool  cited = PETSC_FALSE;
 39: static const char citation[] =
 40:   "@Article{slepc-slice-qep,\n"
 41:   "   author = \"C. Campos and J. E. Roman\",\n"
 42:   "   title = \"A polynomial {Jacobi-Davidson} solver with support for non-monomial bases and deflation\",\n"
 43:   "   journal = \"{BIT} Numer. Math.\",\n"
 44:   "   volume = \"60\",\n"
 45:   "   pages = \"295--318\",\n"
 46:   "   year = \"2020,\"\n"
 47:   "   doi = \"https://doi.org/10.1007/s10543-019-00778-z\"\n"
 48:   "}\n";

 50: typedef struct {
 51:   PetscReal   keep;          /* restart parameter */
 52:   PetscReal   fix;           /* fix parameter */
 53:   PetscBool   reusepc;       /* flag indicating whether pc is rebuilt or not */
 54:   BV          V;             /* work basis vectors to store the search space */
 55:   BV          W;             /* work basis vectors to store the test space */
 56:   BV          *TV;           /* work basis vectors to store T*V (each TV[i] is the coefficient for \lambda^i of T*V for the extended T) */
 57:   BV          *AX;           /* work basis vectors to store A_i*X for locked eigenvectors */
 58:   BV          N[2];          /* auxiliary work BVs */
 59:   BV          X;             /* locked eigenvectors */
 60:   PetscScalar *T;            /* matrix of the invariant pair */
 61:   PetscScalar *Tj;           /* matrix containing the powers of the invariant pair matrix */
 62:   PetscScalar *XpX;          /* X^H*X */
 63:   PetscInt    ld;            /* leading dimension for Tj and XpX */
 64:   PC          pcshell;       /* preconditioner including basic precond+projector */
 65:   Mat         Pshell;        /* auxiliary shell matrix */
 66:   PetscInt    nlock;         /* number of locked vectors in the invariant pair */
 67:   Vec         vtempl;        /* reference nested vector */
 68:   PetscInt    midx;          /* minimality index */
 69:   PetscInt    mmidx;         /* maximum allowed minimality index */
 70:   PEPJDProjection proj;      /* projection type (orthogonal, harmonic) */
 71: } PEP_JD;

 73: typedef struct {
 74:   PEP         pep;
 75:   PC          pc;            /* basic preconditioner */
 76:   Vec         Bp[2];         /* preconditioned residual of derivative polynomial, B\p */
 77:   Vec         u[2];          /* Ritz vector */
 78:   PetscScalar gamma[2];      /* precomputed scalar u'*B\p */
 79:   PetscScalar theta;
 80:   PetscScalar *M;
 81:   PetscScalar *ps;
 82:   PetscInt    ld;
 83:   Vec         *work;
 84:   Mat         PPr;
 85:   BV          X;
 86:   PetscInt    n;
 87: } PEP_JD_PCSHELL;

 89: typedef struct {
 90:   Mat         Pr,Pi;         /* matrix polynomial evaluated at theta */
 91:   PEP         pep;
 92:   Vec         *work;
 93:   PetscScalar theta[2];
 94: } PEP_JD_MATSHELL;

 96: /*
 97:    Duplicate and resize auxiliary basis
 98: */
 99: static PetscErrorCode PEPJDDuplicateBasis(PEP pep,BV *basis)
100: {
101:   PetscErrorCode     ierr;
102:   PEP_JD             *pjd = (PEP_JD*)pep->data;
103:   PetscInt           nloc,m;
104:   BVType             type;
105:   BVOrthogType       otype;
106:   BVOrthogRefineType oref;
107:   PetscReal          oeta;
108:   BVOrthogBlockType  oblock;

111:   if (pjd->ld>1) {
112:     BVCreate(PetscObjectComm((PetscObject)pep),basis);
113:     BVGetSizes(pep->V,&nloc,NULL,&m);
114:     nloc += pjd->ld-1;
115:     BVSetSizes(*basis,nloc,PETSC_DECIDE,m);
116:     BVGetType(pep->V,&type);
117:     BVSetType(*basis,type);
118:     BVGetOrthogonalization(pep->V,&otype,&oref,&oeta,&oblock);
119:     BVSetOrthogonalization(*basis,otype,oref,oeta,oblock);
120:     PetscObjectStateIncrease((PetscObject)*basis);
121:   } else {
122:     BVDuplicate(pep->V,basis);
123:   }
124:   return(0);
125: }

127: PetscErrorCode PEPSetUp_JD(PEP pep)
128: {
130:   PEP_JD         *pjd = (PEP_JD*)pep->data;
131:   PetscBool      isprecond,flg;
132:   PetscInt       i;

135:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
136:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
137:   if (!pep->which) pep->which = PEP_TARGET_MAGNITUDE;
138:   if (pep->which!=PEP_TARGET_MAGNITUDE && pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"The JD solver supports only target which, see PEPSetWhichEigenpairs()");

140:   PetscObjectTypeCompare((PetscObject)pep->st,STPRECOND,&isprecond);
141:   if (!isprecond) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"The JD solver only works with PRECOND spectral transformation");

143:   STGetTransform(pep->st,&flg);
144:   if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"The JD solver requires the ST transform flag unset, see STSetTransform()");
145:   PEPCheckIgnored(pep,PEP_FEATURE_EXTRACT);

147:   if (!pjd->mmidx) pjd->mmidx = pep->nmat-1;
148:   pjd->mmidx = PetscMin(pjd->mmidx,pep->nmat-1);
149:   if (!pjd->keep) pjd->keep = 0.5;
150:   PEPBasisCoefficients(pep,pep->pbc);
151:   PEPAllocateSolution(pep,0);
152:   PEPSetWorkVecs(pep,5);
153:   pjd->ld = pep->nev;
154: #if !defined (PETSC_USE_COMPLEX)
155:   pjd->ld++;
156: #endif
157:   PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
158:   for (i=0;i<pep->nmat;i++) {
159:     PEPJDDuplicateBasis(pep,pjd->TV+i);
160:   }
161:   if (pjd->ld>1) {
162:     PEPJDDuplicateBasis(pep,&pjd->V);
163:     BVSetFromOptions(pjd->V);
164:     for (i=0;i<pep->nmat;i++) {
165:       BVDuplicateResize(pep->V,pjd->ld-1,pjd->AX+i);
166:     }
167:     BVDuplicateResize(pep->V,pjd->ld-1,pjd->N);
168:     BVDuplicateResize(pep->V,pjd->ld-1,pjd->N+1);
169:     pjd->X = pep->V;
170:     PetscCalloc3((pjd->ld)*(pjd->ld),&pjd->XpX,pep->ncv*pep->ncv,&pjd->T,pjd->ld*pjd->ld*pep->nmat,&pjd->Tj);
171:   } else pjd->V = pep->V;
172:   if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { PEPJDDuplicateBasis(pep,&pjd->W); }
173:   else pjd->W = pjd->V;
174:   DSSetType(pep->ds,DSPEP);
175:   DSPEPSetDegree(pep->ds,pep->nmat-1);
176:   if (pep->basis!=PEP_BASIS_MONOMIAL) {
177:     DSPEPSetCoefficients(pep->ds,pep->pbc);
178:   }
179:   DSAllocate(pep->ds,pep->ncv);
180:   return(0);
181: }

183: /*
184:    Updates columns (low to (high-1)) of TV[i]
185: */
186: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
187: {
189:   PEP_JD         *pjd = (PEP_JD*)pep->data;
190:   PetscInt       pp,col,i,nloc,nconv;
191:   Vec            v1,v2,t1,t2;
192:   PetscScalar    *array1,*array2,*x2,*xx,*N,*Np,*y2=NULL,zero=0.0,sone=1.0,*pT,fact,*psc;
193:   PetscReal      *cg,*ca,*cb;
194:   PetscMPIInt    rk,np;
195:   PetscBLASInt   n_,ld_,one=1;
196:   Mat            T;
197:   BV             pbv;

200:   ca = pep->pbc; cb = ca+pep->nmat; cg = cb + pep->nmat;
201:   nconv = pjd->nlock;
202:   PetscMalloc5(nconv,&x2,nconv,&xx,nconv*nconv,&pT,nconv*nconv,&N,nconv*nconv,&Np);
203:   MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
204:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
205:   BVGetSizes(pep->V,&nloc,NULL,NULL);
206:   t1 = w[0];
207:   t2 = w[1];
208:   PetscBLASIntCast(pjd->nlock,&n_);
209:   PetscBLASIntCast(pjd->ld,&ld_);
210:   if (nconv){
211:     for (i=0;i<nconv;i++) {
212:       PetscArraycpy(pT+i*nconv,pjd->T+i*pep->ncv,nconv);
213:     }
214:     MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,pT,&T);
215:   }
216:   for (col=low;col<high;col++) {
217:     BVGetColumn(pjd->V,col,&v1);
218:     VecGetArray(v1,&array1);
219:     if (nconv>0) {
220:       for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
221:     }
222:     VecPlaceArray(t1,array1);
223:     if (nconv) {
224:       BVSetActiveColumns(pjd->N[0],0,nconv);
225:       BVSetActiveColumns(pjd->N[1],0,nconv);
226:       BVDotVec(pjd->X,t1,xx);
227:     }
228:     for (pp=pep->nmat-1;pp>=0;pp--) {
229:       BVGetColumn(pjd->TV[pp],col,&v2);
230:       VecGetArray(v2,&array2);
231:       VecPlaceArray(t2,array2);
232:       MatMult(pep->A[pp],t1,t2);
233:       if (nconv) {
234:         if (pp<pep->nmat-3) {
235:           BVMult(pjd->N[0],1.0,-cg[pp+2],pjd->AX[pp+1],NULL);
236:           MatShift(T,-cb[pp+1]);
237:           BVMult(pjd->N[0],1.0/ca[pp],1.0/ca[pp],pjd->N[1],T);
238:           pbv = pjd->N[0]; pjd->N[0] = pjd->N[1]; pjd->N[1] = pbv;
239:           BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
240:           MatShift(T,cb[pp+1]);
241:         } else if (pp==pep->nmat-3) {
242:           BVCopy(pjd->AX[pp+2],pjd->N[0]);
243:           BVScale(pjd->N[0],1/ca[pp+1]);
244:           BVCopy(pjd->AX[pp+1],pjd->N[1]);
245:           MatShift(T,-cb[pp+1]);
246:           BVMult(pjd->N[1],1.0/ca[pp],1.0/ca[pp],pjd->N[0],T);
247:           BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
248:           MatShift(T,cb[pp+1]);
249:         } else if (pp==pep->nmat-2) {
250:           BVMultVec(pjd->AX[pp+1],1.0/ca[pp],1.0,t2,x2);
251:         }
252:         if (pp<pjd->midx) {
253:           y2 = array2+nloc;
254:           PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&n_,&sone,pjd->Tj+pjd->ld*pjd->ld*pp,&ld_,xx,&one,&zero,y2,&one));
255:           if (pp<pjd->midx-2) {
256:             fact = -cg[pp+2];
257:             PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&fact,Np,&n_));
258:             fact = 1/ca[pp];
259:             MatShift(T,-cb[pp+1]);
260:             PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&fact,N,&n_,pT,&n_,&fact,Np,&n_));
261:             MatShift(T,cb[pp+1]);
262:             psc = Np; Np = N; N = psc;
263:             PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
264:           } else if (pp==pjd->midx-2) {
265:             fact = 1/ca[pp];
266:             PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&fact,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&zero,N,&n_));
267:             PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
268:           } else if (pp==pjd->midx-1) {
269:             PetscArrayzero(Np,nconv*nconv);
270:           }
271:         }
272:         for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
273:       }
274:       VecResetArray(t2);
275:       VecRestoreArray(v2,&array2);
276:       BVRestoreColumn(pjd->TV[pp],col,&v2);
277:     }
278:     VecResetArray(t1);
279:     VecRestoreArray(v1,&array1);
280:     BVRestoreColumn(pjd->V,col,&v1);
281:   }
282:   if (nconv) {MatDestroy(&T);}
283:   PetscFree5(x2,xx,pT,N,Np);
284:   return(0);
285: }

287: /*
288:    RRQR of X. Xin*P=Xou*R. Rank of R is rk
289: */
290: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
291: {
293:   PetscInt       i,j,n,r;
294:   PetscBLASInt   row_,col_,ldx_,*p,lwork,info,n_;
295:   PetscScalar    *tau,*work;
296:   PetscReal      tol,*rwork;

299:   PetscBLASIntCast(row,&row_);
300:   PetscBLASIntCast(col,&col_);
301:   PetscBLASIntCast(ldx,&ldx_);
302:   n = PetscMin(row,col);
303:   PetscBLASIntCast(n,&n_);
304:   lwork = 3*col_+1;
305:   PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
306:   for (i=1;i<col;i++) p[i] = 0;
307:   p[0] = 1;

309:   /* rank revealing QR */
310: #if defined(PETSC_USE_COMPLEX)
311:   PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
312: #else
313:   PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
314: #endif
315:   SlepcCheckLapackInfo("geqp3",info);
316:   if (P) for (i=0;i<col;i++) P[i] = p[i]-1;

318:   /* rank computation */
319:   tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
320:   r = 1;
321:   for (i=1;i<n;i++) {
322:     if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
323:     else break;
324:   }
325:   if (rk) *rk=r;

327:   /* copy upper triangular matrix if requested */
328:   if (R) {
329:      for (i=0;i<r;i++) {
330:        PetscArrayzero(R+i*ldr,r);
331:        for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
332:      }
333:   }
334:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
335:   SlepcCheckLapackInfo("orgqr",info);
336:   PetscFree4(p,tau,work,rwork);
337:   return(0);
338: }

340: /*
341:    Application of extended preconditioner
342: */
343: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
344: {
345:   PetscInt          i,j,nloc,n,ld=0;
346:   PetscMPIInt       np;
347:   Vec               tx,ty;
348:   PEP_JD_PCSHELL    *ctx;
349:   PetscErrorCode    ierr;
350:   const PetscScalar *array1;
351:   PetscScalar       *x2=NULL,*t=NULL,*ps=NULL,*array2,zero=0.0,sone=1.0;
352:   PetscBLASInt      one=1,ld_,n_,ncv_;
353:   PEP_JD            *pjd=NULL;

356:   MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
357:   PCShellGetContext(pc,(void**)&ctx);
358:   n  = ctx->n;
359:   if (n) {
360:     pjd = (PEP_JD*)ctx->pep->data;
361:     ps = ctx->ps;
362:     ld = pjd->ld;
363:     PetscMalloc2(n,&x2,n,&t);
364:     VecGetLocalSize(ctx->work[0],&nloc);
365:     VecGetArrayRead(x,&array1);
366:     for (i=0;i<n;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
367:     VecRestoreArrayRead(x,&array1);
368:   }

370:   /* y = B\x apply PC */
371:   tx = ctx->work[0];
372:   ty = ctx->work[1];
373:   VecGetArrayRead(x,&array1);
374:   VecPlaceArray(tx,array1);
375:   VecGetArray(y,&array2);
376:   VecPlaceArray(ty,array2);
377:   PCApply(ctx->pc,tx,ty);
378:   if (n) {
379:     PetscBLASIntCast(ld,&ld_);
380:     PetscBLASIntCast(n,&n_);
381:     for (i=0;i<n;i++) {
382:       t[i] = 0.0;
383:       for (j=0;j<n;j++) t[i] += ctx->M[i+j*ld]*x2[j];
384:     }
385:     if (pjd->midx==1) {
386:       PetscBLASIntCast(ctx->pep->ncv,&ncv_);
387:       for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] -= ctx->theta;
388:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,pjd->T,&ncv_,t,&one,&zero,x2,&one));
389:       for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] += ctx->theta;
390:       for (i=0;i<n;i++) array2[nloc+i] = x2[i];
391:       for (i=0;i<n;i++) x2[i] = -t[i];
392:     } else {
393:       for (i=0;i<n;i++) array2[nloc+i] = t[i];
394:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,ps,&ld_,t,&one,&zero,x2,&one));
395:     }
396:     for (i=0;i<n;i++) array2[nloc+i] /= PetscSqrtReal(np);
397:     BVSetActiveColumns(pjd->X,0,n);
398:     BVMultVec(pjd->X,-1.0,1.0,ty,x2);
399:     PetscFree2(x2,t);
400:   }
401:   VecResetArray(tx);
402:   VecResetArray(ty);
403:   VecRestoreArrayRead(x,&array1);
404:   VecRestoreArray(y,&array2);
405:   return(0);
406: }

408: /*
409:    Application of shell preconditioner:
410:       y = B\x - eta*B\p,  with eta = (u'*B\x)/(u'*B\p)
411: */
412: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
413: {
415:   PetscScalar    rr,eta;
416:   PEP_JD_PCSHELL *ctx;
417:   PetscInt       sz;
418:   const Vec      *xs,*ys;
419: #if !defined(PETSC_USE_COMPLEX)
420:   PetscScalar    rx,xr,xx;
421: #endif

424:   PCShellGetContext(pc,(void**)&ctx);
425:   VecCompGetSubVecs(x,&sz,&xs);
426:   VecCompGetSubVecs(y,NULL,&ys);
427:   /* y = B\x apply extended PC */
428:   PEPJDExtendedPCApply(pc,xs[0],ys[0]);
429: #if !defined(PETSC_USE_COMPLEX)
430:   if (sz==2) {
431:     PEPJDExtendedPCApply(pc,xs[1],ys[1]);
432:   }
433: #endif

435:   /* Compute eta = u'*y / u'*Bp */
436:   VecDot(ys[0],ctx->u[0],&rr);
437:   eta  = -rr*ctx->gamma[0];
438: #if !defined(PETSC_USE_COMPLEX)
439:   if (sz==2) {
440:     VecDot(ys[0],ctx->u[1],&xr);
441:     VecDot(ys[1],ctx->u[0],&rx);
442:     VecDot(ys[1],ctx->u[1],&xx);
443:     eta += -ctx->gamma[0]*xx-ctx->gamma[1]*(-xr+rx);
444:   }
445: #endif
446:   eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];

448:   /* y = y - eta*Bp */
449:   VecAXPY(ys[0],eta,ctx->Bp[0]);
450: #if !defined(PETSC_USE_COMPLEX)
451:   if (sz==2) {
452:     VecAXPY(ys[1],eta,ctx->Bp[1]);
453:     eta = -ctx->gamma[1]*(rr+xx)+ctx->gamma[0]*(-xr+rx);
454:     eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];
455:     VecAXPY(ys[0],eta,ctx->Bp[1]);
456:     VecAXPY(ys[1],-eta,ctx->Bp[0]);
457:   }
458: #endif
459:   return(0);
460: }

462: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
463: {
465:   PetscMPIInt    np,rk,count;
466:   PetscScalar    *array1,*array2;
467:   PetscInt       nloc;

470:   MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
471:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
472:   BVGetSizes(pep->V,&nloc,NULL,NULL);
473:   if (v) {
474:     VecGetArray(v,&array1);
475:     VecGetArray(vex,&array2);
476:     if (back) {
477:       PetscArraycpy(array1,array2,nloc);
478:     } else {
479:       PetscArraycpy(array2,array1,nloc);
480:     }
481:     VecRestoreArray(v,&array1);
482:     VecRestoreArray(vex,&array2);
483:   }
484:   if (a) {
485:     VecGetArray(vex,&array2);
486:     if (back) {
487:       PetscArraycpy(a,array2+nloc+off,na);
488:       PetscMPIIntCast(na,&count);
489:       MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
490:     } else {
491:       PetscArraycpy(array2+nloc+off,a,na);
492:       PetscMPIIntCast(na,&count);
493:       MPI_Bcast(array2+nloc+off,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
494:     }
495:     VecRestoreArray(vex,&array2);
496:   }
497:   return(0);
498: }

500: /* Computes Phi^hat(lambda) times a vector or its derivative (depends on beval)
501:      if no vector is provided returns a matrix
502:  */
503: static PetscErrorCode PEPJDEvaluateHatBasis(PEP pep,PetscInt n,PetscScalar *H,PetscInt ldh,PetscScalar *beval,PetscScalar *t,PetscInt idx,PetscScalar *qpp,PetscScalar *qp,PetscScalar *q)
504: {
506:   PetscInt       j,i;
507:   PetscBLASInt   n_,ldh_,one=1;
508:   PetscReal      *a,*b,*g;
509:   PetscScalar    sone=1.0,zero=0.0;

512:   a = pep->pbc; b=a+pep->nmat; g=b+pep->nmat;
513:   PetscBLASIntCast(n,&n_);
514:   PetscBLASIntCast(ldh,&ldh_);
515:   if (idx<1) {
516:     PetscArrayzero(q,t?n:n*n);
517:   } else if (idx==1) {
518:     if (t) {for (j=0;j<n;j++) q[j] = t[j]*beval[idx-1]/a[0];}
519:     else {
520:       PetscArrayzero(q,n*n);
521:       for (j=0;j<n;j++) q[(j+1)*n] = beval[idx-1]/a[0];
522:     }
523:   } else {
524:     if (t) {
525:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,H,&ldh_,qp,&one,&zero,q,&one));
526:       for (j=0;j<n;j++) {
527:         q[j] += beval[idx-1]*t[j]-b[idx-1]*qp[j]-g[idx-1]*qpp[j];
528:         q[j] /= a[idx-1];
529:       }
530:     } else {
531:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,H,&ldh_,qp,&n_,&zero,q,&n_));
532:       for (j=0;j<n;j++) {
533:         q[j+n*j] += beval[idx-1];
534:         for (i=0;i<n;i++) {
535:           q[i+n*j] += -b[idx-1]*qp[j*n+i]-g[idx-1]*qpp[j*n+i];
536:           q[i+n*j] /= a[idx-1];
537:         }
538:       }
539:     }
540:   }
541:   return(0);
542: }

544: static PetscErrorCode PEPJDComputeResidual(PEP pep,PetscBool derivative,PetscInt sz,Vec *u,PetscScalar *theta,Vec *p,Vec *work)
545: {
546:   PEP_JD         *pjd = (PEP_JD*)pep->data;
548:   PetscMPIInt    rk,np,count;
549:   Vec            tu,tp,w;
550:   PetscScalar    *dval,*dvali,*array1,*array2,*x2=NULL,*y2,*qj=NULL,*tt=NULL,*xx=NULL,*xxi=NULL,sone=1.0;
551:   PetscInt       i,j,nconv,nloc;
552:   PetscBLASInt   n,ld,one=1;
553: #if !defined(PETSC_USE_COMPLEX)
554:   Vec            tui=NULL,tpi=NULL;
555:   PetscScalar    *x2i=NULL,*qji=NULL,*qq,*y2i,*arrayi1,*arrayi2;
556: #endif

559:   nconv = pjd->nlock;
560:   if (!nconv) {
561:     PetscMalloc1(2*sz*pep->nmat,&dval);
562:   } else {
563:     PetscMalloc5(2*pep->nmat,&dval,2*nconv,&xx,nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*pep->nmat,&qj);
564:     MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
565:     MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
566:     BVGetSizes(pep->V,&nloc,NULL,NULL);
567:     VecGetArray(u[0],&array1);
568:     for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]*PetscSqrtReal(np);
569:     VecRestoreArray(u[0],&array1);
570: #if !defined(PETSC_USE_COMPLEX)
571:     if (sz==2) {
572:       x2i = x2+nconv;
573:       VecGetArray(u[1],&arrayi1);
574:       for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]*PetscSqrtReal(np);
575:       VecRestoreArray(u[1],&arrayi1);
576:     }
577: #endif
578:   }
579:   dvali = dval+pep->nmat;
580:   tu = work[0];
581:   tp = work[1];
582:   w  = work[2];
583:   VecGetArray(u[0],&array1);
584:   VecPlaceArray(tu,array1);
585:   VecGetArray(p[0],&array2);
586:   VecPlaceArray(tp,array2);
587:   VecSet(tp,0.0);
588: #if !defined(PETSC_USE_COMPLEX)
589:   if (sz==2) {
590:     tui = work[3];
591:     tpi = work[4];
592:     VecGetArray(u[1],&arrayi1);
593:     VecPlaceArray(tui,arrayi1);
594:     VecGetArray(p[1],&arrayi2);
595:     VecPlaceArray(tpi,arrayi2);
596:     VecSet(tpi,0.0);
597:   }
598: #endif
599:   if (derivative) {
600:     PEPEvaluateBasisDerivative(pep,theta[0],theta[1],dval,dvali);
601:   } else {
602:     PEPEvaluateBasis(pep,theta[0],theta[1],dval,dvali);
603:   }
604:   for (i=derivative?1:0;i<pep->nmat;i++) {
605:     MatMult(pep->A[i],tu,w);
606:     VecAXPY(tp,dval[i],w);
607: #if !defined(PETSC_USE_COMPLEX)
608:     if (sz==2) {
609:       VecAXPY(tpi,dvali[i],w);
610:       MatMult(pep->A[i],tui,w);
611:       VecAXPY(tpi,dval[i],w);
612:       VecAXPY(tp,-dvali[i],w);
613:     }
614: #endif
615:   }
616:   if (nconv) {
617:     for (i=0;i<pep->nmat;i++) {
618:       PEPJDEvaluateHatBasis(pep,nconv,pjd->T,pep->ncv,dval,x2,i,i>1?qj+(i-2)*nconv:NULL,i>0?qj+(i-1)*nconv:NULL,qj+i*nconv);
619:     }
620: #if !defined(PETSC_USE_COMPLEX)
621:     if (sz==2) {
622:       qji = qj+nconv*pep->nmat;
623:       qq = qji+nconv*pep->nmat;
624:       for (i=0;i<pep->nmat;i++) {
625:         PEPJDEvaluateHatBasis(pep,nconv,pjd->T,pep->ncv,dvali,x2i,i,i>1?qji+(i-2)*nconv:NULL,i>0?qji+(i-1)*nconv:NULL,qji+i*nconv);
626:       }
627:       for (i=0;i<nconv*pep->nmat;i++) qj[i] -= qji[i];
628:       for (i=0;i<pep->nmat;i++) {
629:         PEPJDEvaluateHatBasis(pep,nconv,pjd->T,pep->ncv,dval,x2i,i,i>1?qji+(i-2)*nconv:NULL,i>0?qji+(i-1)*nconv:NULL,qji+i*nconv);
630:         PEPJDEvaluateHatBasis(pep,nconv,pjd->T,pep->ncv,dvali,x2,i,i>1?qq+(i-2)*nconv:NULL,i>0?qq+(i-1)*nconv:NULL,qq+i*nconv);
631:       }
632:       for (i=0;i<nconv*pep->nmat;i++) qji[i] += qq[i];
633:       for (i=derivative?2:1;i<pep->nmat;i++) {
634:         BVMultVec(pjd->AX[i],1.0,1.0,tpi,qji+i*nconv);
635:       }
636:     }
637: #endif
638:     for (i=derivative?2:1;i<pep->nmat;i++) {
639:       BVMultVec(pjd->AX[i],1.0,1.0,tp,qj+i*nconv);
640:     }

642:     /* extended vector part */
643:     BVSetActiveColumns(pjd->X,0,nconv);
644:     BVDotVec(pjd->X,tu,xx);
645:     xxi = xx+nconv;
646: #if !defined(PETSC_USE_COMPLEX)
647:     if (sz==2) { BVDotVec(pjd->X,tui,xxi); }
648: #endif
649:     if (sz==1) { PetscArrayzero(xxi,nconv); }
650:     if (rk==np-1) {
651:       PetscBLASIntCast(nconv,&n);
652:       PetscBLASIntCast(pjd->ld,&ld);
653:       y2  = array2+nloc;
654:       PetscArrayzero(y2,nconv);
655:       for (j=derivative?1:0;j<pjd->midx;j++) {
656:         for (i=0;i<nconv;i++) tt[i] = dval[j]*xx[i]-dvali[j]*xxi[i];
657:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
658:         PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
659:       }
660:       for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
661: #if !defined(PETSC_USE_COMPLEX)
662:       if (sz==2) {
663:         y2i = arrayi2+nloc;
664:         PetscArrayzero(y2i,nconv);
665:         for (j=derivative?1:0;j<pjd->midx;j++) {
666:           for (i=0;i<nconv;i++) tt[i] = dval[j]*xxi[i]+dvali[j]*xx[i];
667:           PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
668:           PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
669:         }
670:         for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
671:       }
672: #endif
673:     }
674:     PetscMPIIntCast(nconv,&count);
675:     MPI_Bcast(array2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
676: #if !defined(PETSC_USE_COMPLEX)
677:     if (sz==2) {
678:       MPI_Bcast(arrayi2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
679:     }
680: #endif
681:   }
682:   if (nconv) {
683:     PetscFree5(dval,xx,tt,x2,qj);
684:   } else {
685:     PetscFree(dval);
686:   }
687:   VecResetArray(tu);
688:   VecRestoreArray(u[0],&array1);
689:   VecResetArray(tp);
690:   VecRestoreArray(p[0],&array2);
691: #if !defined(PETSC_USE_COMPLEX)
692:   if (sz==2) {
693:     VecResetArray(tui);
694:     VecRestoreArray(u[1],&arrayi1);
695:     VecResetArray(tpi);
696:     VecRestoreArray(p[1],&arrayi2);
697:   }
698: #endif
699:   return(0);
700: }

702: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
703: {
704:   PEP_JD         *pjd = (PEP_JD*)pep->data;
706:   PetscScalar    *tt,target[2];
707:   Vec            vg,wg;
708:   PetscInt       i;
709:   PetscReal      norm;

712:   PetscMalloc1(pjd->ld-1,&tt);
713:   if (pep->nini==0) {
714:     BVSetRandomColumn(pjd->V,0);
715:     for (i=0;i<pjd->ld-1;i++) tt[i] = 0.0;
716:     BVGetColumn(pjd->V,0,&vg);
717:     PEPJDCopyToExtendedVec(pep,NULL,tt,pjd->ld-1,0,vg,PETSC_FALSE);
718:     BVRestoreColumn(pjd->V,0,&vg);
719:     BVNormColumn(pjd->V,0,NORM_2,&norm);
720:     BVScaleColumn(pjd->V,0,1.0/norm);
721:     if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
722:       BVGetColumn(pjd->V,0,&vg);
723:       BVGetColumn(pjd->W,0,&wg);
724:       VecSet(wg,0.0);
725:       target[0] = pep->target; target[1] = 0.0;
726:       PEPJDComputeResidual(pep,PETSC_TRUE,1,&vg,target,&wg,w);
727:       BVRestoreColumn(pjd->W,0,&wg);
728:       BVRestoreColumn(pjd->V,0,&vg);
729:       BVNormColumn(pjd->W,0,NORM_2,&norm);
730:       BVScaleColumn(pjd->W,0,1.0/norm);
731:     }
732:   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for initial vectors not implemented yet");
733:   PetscFree(tt);
734:   return(0);
735: }

737: static PetscErrorCode MatMult_PEPJD(Mat P,Vec x,Vec y)
738: {
739:   PetscErrorCode  ierr;
740:   PEP_JD_MATSHELL *matctx;
741:   PEP_JD          *pjd;
742:   PetscInt        i,j,nconv,nloc,nmat,ldt,ncv,sz;
743:   Vec             tx,ty;
744:   const Vec       *xs,*ys;
745:   PetscScalar     *array1,*array2,*x2=NULL,*y2,*tt=NULL,*xx=NULL,*xxi,theta[2],sone=1.0,*qj,*val,*vali=NULL;
746:   PetscBLASInt    n,ld,one=1;
747:   PetscMPIInt     np;
748: #if !defined(PETSC_USE_COMPLEX)
749:   Vec             txi=NULL,tyi=NULL;
750:   PetscScalar     *x2i=NULL,*qji=NULL,*qq,*y2i,*arrayi1,*arrayi2;
751: #endif

754:   MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
755:   MatShellGetContext(P,(void**)&matctx);
756:   pjd   = (PEP_JD*)(matctx->pep->data);
757:   nconv = pjd->nlock;
758:   nmat  = matctx->pep->nmat;
759:   ncv   = matctx->pep->ncv;
760:   ldt   = pjd->ld;
761:   VecCompGetSubVecs(x,&sz,&xs);
762:   VecCompGetSubVecs(y,NULL,&ys);
763:   theta[0] = matctx->theta[0];
764:   theta[1] = (sz==2)?matctx->theta[1]:0.0;
765:   if (nconv>0) {
766:     PetscMalloc5(nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*nmat,&qj,2*nconv,&xx,2*nmat,&val);
767:     BVGetSizes(matctx->pep->V,&nloc,NULL,NULL);
768:     VecGetArray(xs[0],&array1);
769:     for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
770:     VecRestoreArray(xs[0],&array1);
771: #if !defined(PETSC_USE_COMPLEX)
772:     if (sz==2) {
773:       x2i = x2+nconv;
774:       VecGetArray(xs[1],&arrayi1);
775:       for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]* PetscSqrtReal(np);
776:       VecRestoreArray(xs[1],&arrayi1);
777:     }
778: #endif
779:     vali = val+nmat;
780:   }
781:   tx = matctx->work[0];
782:   ty = matctx->work[1];
783:   VecGetArray(xs[0],&array1);
784:   VecPlaceArray(tx,array1);
785:   VecGetArray(ys[0],&array2);
786:   VecPlaceArray(ty,array2);
787:   MatMult(matctx->Pr,tx,ty);
788: #if !defined(PETSC_USE_COMPLEX)
789:   if (sz==2) {
790:     txi = matctx->work[2];
791:     tyi = matctx->work[3];
792:     VecGetArray(xs[1],&arrayi1);
793:     VecPlaceArray(txi,arrayi1);
794:     VecGetArray(ys[1],&arrayi2);
795:     VecPlaceArray(tyi,arrayi2);
796:     MatMult(matctx->Pr,txi,tyi);
797:     if (theta[1]!=0.0) {
798:       MatMult(matctx->Pi,txi,matctx->work[4]);
799:       VecAXPY(ty,-1.0,matctx->work[4]);
800:       MatMult(matctx->Pi,tx,matctx->work[4]);
801:       VecAXPY(tyi,1.0,matctx->work[4]);
802:     }
803:   }
804: #endif
805:   if (nconv>0) {
806:     PEPEvaluateBasis(matctx->pep,theta[0],theta[1],val,vali);
807:     for (i=0;i<nmat;i++) {
808:       PEPJDEvaluateHatBasis(matctx->pep,nconv,pjd->T,ncv,val,x2,i,i>1?qj+(i-2)*nconv:NULL,i>0?qj+(i-1)*nconv:NULL,qj+i*nconv);
809:     }
810: #if !defined(PETSC_USE_COMPLEX)
811:     if (sz==2) {
812:       qji = qj+nconv*nmat;
813:       qq = qji+nconv*nmat;
814:       for (i=0;i<nmat;i++) {
815:         PEPJDEvaluateHatBasis(matctx->pep,nconv,pjd->T,matctx->pep->ncv,vali,x2i,i,i>1?qji+(i-2)*nconv:NULL,i>0?qji+(i-1)*nconv:NULL,qji+i*nconv);
816:       }
817:       for (i=0;i<nconv*nmat;i++) qj[i] -= qji[i];
818:       for (i=0;i<nmat;i++) {
819:         PEPJDEvaluateHatBasis(matctx->pep,nconv,pjd->T,matctx->pep->ncv,val,x2i,i,i>1?qji+(i-2)*nconv:NULL,i>0?qji+(i-1)*nconv:NULL,qji+i*nconv);
820:         PEPJDEvaluateHatBasis(matctx->pep,nconv,pjd->T,matctx->pep->ncv,vali,x2,i,i>1?qq+(i-2)*nconv:NULL,i>0?qq+(i-1)*nconv:NULL,qq+i*nconv);
821:       }
822:       for (i=0;i<nconv*nmat;i++) qji[i] += qq[i];
823:       for (i=1;i<matctx->pep->nmat;i++) {
824:         BVMultVec(pjd->AX[i],1.0,1.0,tyi,qji+i*nconv);
825:       }
826:     }
827: #endif
828:     for (i=1;i<nmat;i++) {
829:       BVMultVec(pjd->AX[i],1.0,1.0,ty,qj+i*nconv);
830:     }

832:     /* extended vector part */
833:     BVSetActiveColumns(pjd->X,0,nconv);
834:     BVDotVec(pjd->X,tx,xx);
835:     xxi = xx+nconv;
836: #if !defined(PETSC_USE_COMPLEX)
837:     if (sz==2) { BVDotVec(pjd->X,txi,xxi); }
838: #endif
839:     if (sz==1) { PetscArrayzero(xxi,nconv); }
840:     PetscBLASIntCast(pjd->nlock,&n);
841:     PetscBLASIntCast(ldt,&ld);
842:     y2 = array2+nloc;
843:     PetscArrayzero(y2,nconv);
844:     for (j=0;j<pjd->midx;j++) {
845:       for (i=0;i<nconv;i++) tt[i] = val[j]*xx[i]-vali[j]*xxi[i];
846:       PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
847:       PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
848:     }
849: #if !defined(PETSC_USE_COMPLEX)
850:     if (sz==2) {
851:       y2i = arrayi2+nloc;
852:       PetscArrayzero(y2i,nconv);
853:       for (j=0;j<pjd->midx;j++) {
854:         for (i=0;i<nconv;i++) tt[i] = val[j]*xxi[i]+vali[j]*xx[i];
855:         PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
856:         PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
857:       }
858:       for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
859:     }
860: #endif
861:     for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
862:     PetscFree5(tt,x2,qj,xx,val);
863:   }
864:   VecResetArray(tx);
865:   VecRestoreArray(xs[0],&array1);
866:   VecResetArray(ty);
867:   VecRestoreArray(ys[0],&array2);
868: #if !defined(PETSC_USE_COMPLEX)
869:   if (sz==2) {
870:     VecResetArray(txi);
871:     VecRestoreArray(xs[1],&arrayi1);
872:     VecResetArray(tyi);
873:     VecRestoreArray(ys[1],&arrayi2);
874:   }
875: #endif
876:   return(0);
877: }

879: static PetscErrorCode MatCreateVecs_PEPJD(Mat A,Vec *right,Vec *left)
880: {
881:   PetscErrorCode  ierr;
882:   PEP_JD_MATSHELL *matctx;
883:   PEP_JD          *pjd;
884:   PetscInt        kspsf=1,i;
885:   Vec             v[2];

888:   MatShellGetContext(A,(void**)&matctx);
889:   pjd   = (PEP_JD*)(matctx->pep->data);
890: #if !defined (PETSC_USE_COMPLEX)
891:   kspsf = 2;
892: #endif
893:   for (i=0;i<kspsf;i++){
894:     BVCreateVec(pjd->V,v+i);
895:   }
896:   if (right) {
897:     VecCreateCompWithVecs(v,kspsf,pjd->vtempl,right);
898:   }
899:   if (left) {
900:     VecCreateCompWithVecs(v,kspsf,pjd->vtempl,left);
901:   }
902:   for (i=0;i<kspsf;i++) {
903:     VecDestroy(&v[i]);
904:   }
905:   return(0);
906: }

908: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
909: {
911:   PEP_JD         *pjd = (PEP_JD*)pep->data;
912:   PEP_JD_PCSHELL *pcctx;
913:   PetscInt       i,j,k,n=pjd->nlock,ld=pjd->ld,deg=pep->nmat-1;
914:   PetscScalar    *M,*ps,*work,*U,*V,*S,*Sp,*Spp,snone=-1.0,sone=1.0,zero=0.0,*val;
915:   PetscReal      tol,maxeig=0.0,*sg,*rwork;
916:   PetscBLASInt   n_,info,ld_,*p,lw_,rk=0;

919:   if (n) {
920:     PCShellGetContext(pjd->pcshell,(void**)&pcctx);
921:     pcctx->theta = theta;
922:     pcctx->n = n;
923:     M  = pcctx->M;
924:     PetscBLASIntCast(n,&n_);
925:     PetscBLASIntCast(ld,&ld_);
926:     if (pjd->midx==1) {
927:       PetscArraycpy(M,pjd->XpX,ld*ld);
928:       PetscCalloc2(10*n,&work,n,&p);
929:     } else {
930:       ps = pcctx->ps;
931:       PetscCalloc7(2*n*n,&U,3*n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p,deg+1,&val);
932:       V = U+n*n;
933:       /* pseudo-inverse */
934:       for (j=0;j<n;j++) {
935:         for (i=0;i<n;i++) S[n*j+i] = -pjd->T[pep->ncv*j+i];
936:         S[n*j+j] += theta;
937:       }
938:       lw_ = 10*n_;
939: #if !defined (PETSC_USE_COMPLEX)
940:       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
941: #else
942:       PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
943: #endif
944:       SlepcCheckLapackInfo("gesvd",info);
945:       for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
946:       tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
947:       for (j=0;j<n;j++) {
948:         if (sg[j]>tol) {
949:           for (i=0;i<n;i++) U[j*n+i] /= sg[j];
950:           rk++;
951:         } else break;
952:       }
953:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));

955:       /* compute M */
956:       PEPEvaluateBasis(pep,theta,0.0,val,NULL);
957:       PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&snone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
958:       PetscArrayzero(S,2*n*n);
959:       Sp = S+n*n;
960:       for (j=0;j<n;j++) S[j*(n+1)] = 1.0;
961:       for (k=1;k<pjd->midx;k++) {
962:         for (j=0;j<n;j++) for (i=0;i<n;i++) V[j*n+i] = S[j*n+i] - ps[j*ld+i]*val[k];
963:         PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
964:         PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
965:         Spp = Sp; Sp = S;
966:         PEPJDEvaluateHatBasis(pep,n,pjd->T,pep->ncv,val,NULL,k+1,Spp,Sp,S);
967:       }
968:     }
969:     /* inverse */
970:     PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
971:     SlepcCheckLapackInfo("getrf",info);
972:     PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
973:     SlepcCheckLapackInfo("getri",info);
974:     if (pjd->midx==1) {
975:       PetscFree2(work,p);
976:     } else {
977:       PetscFree7(U,S,sg,work,rwork,p,val);
978:     }
979:   }
980:   return(0);
981: }

983: static PetscErrorCode PEPJDMatSetUp(PEP pep,PetscInt sz,PetscScalar *theta)
984: {
985:   PetscErrorCode  ierr;
986:   PEP_JD          *pjd = (PEP_JD*)pep->data;
987:   PEP_JD_MATSHELL *matctx;
988:   PEP_JD_PCSHELL  *pcctx;
989:   MatStructure    str;
990:   PetscScalar     *vals,*valsi;
991:   PetscBool       skipmat=PETSC_FALSE;
992:   PetscInt        i;
993:   Mat             Pr=NULL;

996:   if (sz==2 && theta[1]==0.0) sz = 1;
997:   MatShellGetContext(pjd->Pshell,(void**)&matctx);
998:   PCShellGetContext(pjd->pcshell,(void**)&pcctx);
999:   if (matctx->Pr && matctx->theta[0]==theta[0] && ((!matctx->Pi && sz==1) || (sz==2 && matctx->theta[1]==theta[1]))) {
1000:     if (pcctx->n == pjd->nlock) return(0);
1001:     skipmat = PETSC_TRUE;
1002:   }
1003:   if (!skipmat) {
1004:     PetscMalloc2(pep->nmat,&vals,pep->nmat,&valsi);
1005:     STGetMatStructure(pep->st,&str);
1006:     PEPEvaluateBasis(pep,theta[0],theta[1],vals,valsi);
1007:     if (!matctx->Pr) {
1008:       MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pr);
1009:     } else {
1010:       MatCopy(pep->A[0],matctx->Pr,str);
1011:     }
1012:     for (i=1;i<pep->nmat;i++) {
1013:       MatAXPY(matctx->Pr,vals[i],pep->A[i],str);
1014:     }
1015:     if (!pjd->reusepc) {
1016:       if (pcctx->PPr && sz==2) {
1017:         MatCopy(matctx->Pr,pcctx->PPr,str);
1018:         Pr = pcctx->PPr;
1019:       } else Pr = matctx->Pr;
1020:     }
1021:     matctx->theta[0] = theta[0];
1022: #if !defined(PETSC_USE_COMPLEX)
1023:     if (sz==2) {
1024:       if (!matctx->Pi) {
1025:         MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pi);
1026:       } else {
1027:         MatCopy(pep->A[1],matctx->Pi,str);
1028:       }
1029:       MatScale(matctx->Pi,valsi[1]);
1030:       for (i=2;i<pep->nmat;i++) {
1031:         MatAXPY(matctx->Pi,valsi[i],pep->A[i],str);
1032:       }
1033:       matctx->theta[1] = theta[1];
1034:     }
1035: #endif
1036:     PetscFree2(vals,valsi);
1037:   }
1038:   if (!pjd->reusepc) {
1039:     if (!skipmat) {
1040:       PCSetOperators(pcctx->pc,Pr,Pr);
1041:       PCSetUp(pcctx->pc);
1042:     }
1043:     PEPJDUpdateExtendedPC(pep,theta[0]);
1044:   }
1045:   return(0);
1046: }

1048: static PetscErrorCode PEPJDCreateShellPC(PEP pep,Vec *ww)
1049: {
1050:   PEP_JD          *pjd = (PEP_JD*)pep->data;
1051:   PEP_JD_PCSHELL  *pcctx;
1052:   PEP_JD_MATSHELL *matctx;
1053:   KSP             ksp;
1054:   PetscInt        nloc,mloc,kspsf=1;
1055:   Vec             v[2];
1056:   PetscScalar     target[2];
1057:   PetscErrorCode  ierr;
1058:   Mat             Pr;

1061:   /* Create the reference vector */
1062:   BVGetColumn(pjd->V,0,&v[0]);
1063:   v[1] = v[0];
1064: #if !defined (PETSC_USE_COMPLEX)
1065:   kspsf = 2;
1066: #endif
1067:   VecCreateCompWithVecs(v,kspsf,NULL,&pjd->vtempl);
1068:   BVRestoreColumn(pjd->V,0,&v[0]);
1069:   PetscLogObjectParent((PetscObject)pep,(PetscObject)pjd->vtempl);

1071:   /* Replace preconditioner with one containing projectors */
1072:   PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
1073:   PCSetType(pjd->pcshell,PCSHELL);
1074:   PCShellSetName(pjd->pcshell,"PCPEPJD");
1075:   PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
1076:   PetscNew(&pcctx);
1077:   PCShellSetContext(pjd->pcshell,pcctx);
1078:   STGetKSP(pep->st,&ksp);
1079:   BVCreateVec(pjd->V,&pcctx->Bp[0]);
1080:   VecDuplicate(pcctx->Bp[0],&pcctx->Bp[1]);
1081:   KSPGetPC(ksp,&pcctx->pc);
1082:   PetscObjectReference((PetscObject)pcctx->pc);
1083:   MatGetLocalSize(pep->A[0],&mloc,&nloc);
1084:   if (pjd->ld>1) {
1085:     nloc += pjd->ld-1; mloc += pjd->ld-1;
1086:   }
1087:   PetscNew(&matctx);
1088:   MatCreateShell(PetscObjectComm((PetscObject)pep),kspsf*nloc,kspsf*mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
1089:   MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)(void))MatMult_PEPJD);
1090:   MatShellSetOperation(pjd->Pshell,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_PEPJD);
1091:   matctx->pep = pep;
1092:   target[0] = pep->target; target[1] = 0.0;
1093:   PEPJDMatSetUp(pep,1,target);
1094:   Pr = matctx->Pr;
1095:   pcctx->PPr = NULL;
1096: #if !defined(PETSC_USE_COMPLEX)
1097:   if (!pjd->reusepc) {
1098:     MatDuplicate(matctx->Pr,MAT_COPY_VALUES,&pcctx->PPr);
1099:     Pr = pcctx->PPr;
1100:   }
1101: #endif
1102:   PCSetOperators(pcctx->pc,Pr,Pr);
1103:   PCSetErrorIfFailure(pcctx->pc,PETSC_TRUE);
1104:   KSPSetPC(ksp,pjd->pcshell);
1105:   if (pjd->reusepc) {
1106:     PCSetReusePreconditioner(pcctx->pc,PETSC_TRUE);
1107:     KSPSetReusePreconditioner(ksp,PETSC_TRUE);
1108:   }
1109:   KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
1110:   KSPSetUp(ksp);
1111:   if (pjd->ld>1) {
1112:     PetscMalloc2(pjd->ld*pjd->ld,&pcctx->M,pjd->ld*pjd->ld,&pcctx->ps);
1113:     pcctx->pep = pep;
1114:   }
1115:   matctx->work = ww;
1116:   pcctx->work  = ww;
1117:   return(0);
1118: }

1120: static PetscErrorCode PEPJDEigenvectors(PEP pep)
1121: {
1123:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1124:   PetscBLASInt   ld,nconv,info,nc;
1125:   PetscScalar    *Z;
1126:   PetscReal      *wr;
1127:   Mat            U;
1128: #if defined(PETSC_USE_COMPLEX)
1129:   PetscScalar    *w;
1130: #endif

1133:   PetscBLASIntCast(pep->ncv,&ld);
1134:   PetscBLASIntCast(pep->nconv,&nconv);
1135: #if !defined(PETSC_USE_COMPLEX)
1136:   PetscMalloc2(pep->nconv*pep->nconv,&Z,3*pep->ncv,&wr);
1137:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
1138: #else
1139:   PetscMalloc3(pep->nconv*pep->nconv,&Z,3*pep->ncv,&wr,2*pep->ncv,&w);
1140:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
1141: #endif
1142:   SlepcCheckLapackInfo("trevc",info);
1143:   MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
1144:   BVSetActiveColumns(pjd->X,0,pep->nconv);
1145:   BVMultInPlace(pjd->X,U,0,pep->nconv);
1146:   BVNormalize(pjd->X,pep->eigi);
1147:   MatDestroy(&U);
1148: #if !defined(PETSC_USE_COMPLEX)
1149:   PetscFree2(Z,wr);
1150: #else
1151:   PetscFree3(Z,wr,w);
1152: #endif
1153:   return(0);
1154: }

1156: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv,PetscInt sz)
1157: {
1159:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1160:   PetscInt       j,i,*P,ldds,rk=0,nvv=*nv;
1161:   Vec            v,x,w;
1162:   PetscScalar    *R,*r,*pX,target[2];
1163:   Mat            X;
1164:   PetscBLASInt   sz_,rk_,nv_,info;
1165:   PetscMPIInt    np;

1168:   /* update AX and XpX */
1169:   for (i=sz;i>0;i--) {
1170:     BVGetColumn(pjd->X,pjd->nlock-i,&x);
1171:     for (j=0;j<pep->nmat;j++) {
1172:       BVGetColumn(pjd->AX[j],pjd->nlock-i,&v);
1173:       MatMult(pep->A[j],x,v);
1174:       BVRestoreColumn(pjd->AX[j],pjd->nlock-i,&v);
1175:       BVSetActiveColumns(pjd->AX[j],0,pjd->nlock-i+1);
1176:     }
1177:     BVRestoreColumn(pjd->X,pjd->nlock-i,&x);
1178:     BVDotColumn(pjd->X,(pjd->nlock-i),pjd->XpX+(pjd->nlock-i)*(pjd->ld));
1179:     pjd->XpX[(pjd->nlock-i)*(1+pjd->ld)] = 1.0;
1180:     for (j=0;j<pjd->nlock-i;j++) pjd->XpX[j*(pjd->ld)+pjd->nlock-i] = PetscConj(pjd->XpX[(pjd->nlock-i)*(pjd->ld)+j]);
1181:   }

1183:   /* minimality index */
1184:   pjd->midx = PetscMin(pjd->mmidx,pjd->nlock);

1186:   /* evaluate the polynomial basis in T */
1187:   PetscArrayzero(pjd->Tj,pjd->ld*pjd->ld*pep->nmat);
1188:   for (j=0;j<pep->nmat;j++) {
1189:     PEPEvaluateBasisMat(pep,pjd->nlock,pjd->T,pep->ncv,j,(j>1)?pjd->Tj+(j-2)*pjd->ld*pjd->ld:NULL,pjd->ld,j?pjd->Tj+(j-1)*pjd->ld*pjd->ld:NULL,pjd->ld,pjd->Tj+j*pjd->ld*pjd->ld,pjd->ld);
1190:   }

1192:   /* Extend search space */
1193:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1194:   PetscCalloc3(nvv,&P,nvv*nvv,&R,nvv*sz,&r);
1195:   DSGetLeadingDimension(pep->ds,&ldds);
1196:   DSGetArray(pep->ds,DS_MAT_X,&pX);
1197:   PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
1198:   for (j=0;j<sz;j++) {
1199:     for (i=0;i<rk;i++) r[i*sz+j] = PetscConj(R[nvv*i+j]*pep->eigr[P[i]]); /* first row scaled with permuted diagonal */
1200:   }
1201:   PetscBLASIntCast(rk,&rk_);
1202:   PetscBLASIntCast(sz,&sz_);
1203:   PetscBLASIntCast(nvv,&nv_);
1204:   PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
1205:   SlepcCheckLapackInfo("trtri",info);
1206:   for (i=0;i<sz;i++) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r+i,&sz_));
1207:   for (i=0;i<sz*rk;i++) r[i] = PetscConj(r[i])/PetscSqrtReal(np); /* revert */
1208:   BVSetActiveColumns(pjd->V,0,nvv);
1209:   rk -= sz;
1210:   for (j=0;j<rk;j++) {
1211:     PetscArraycpy(R+j*nvv,pX+(j+sz)*ldds,nvv);
1212:   }
1213:   DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1214:   MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk,R,&X);
1215:   BVMultInPlace(pjd->V,X,0,rk);
1216:   MatDestroy(&X);
1217:   BVSetActiveColumns(pjd->V,0,rk);
1218:   for (j=0;j<rk;j++) {
1219:     BVGetColumn(pjd->V,j,&v);
1220:     PEPJDCopyToExtendedVec(pep,NULL,r+sz*(j+sz),sz,pjd->nlock-sz,v,PETSC_FALSE);
1221:     BVRestoreColumn(pjd->V,j,&v);
1222:   }
1223:   BVOrthogonalize(pjd->V,NULL);

1225:   if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1226:     for (j=0;j<rk;j++) {
1227:       /* W = P(target)*V */
1228:       BVGetColumn(pjd->W,j,&w);
1229:       BVGetColumn(pjd->V,j,&v);
1230:       target[0] = pep->target; target[1] = 0.0;
1231:       PEPJDComputeResidual(pep,PETSC_FALSE,1,&v,target,&w,pep->work);
1232:       BVRestoreColumn(pjd->V,j,&v);
1233:       BVRestoreColumn(pjd->W,j,&w);
1234:     }
1235:     BVSetActiveColumns(pjd->W,0,rk);
1236:     BVOrthogonalize(pjd->W,NULL);
1237:   }
1238:   *nv = rk;
1239:   PetscFree3(P,R,r);
1240:   return(0);
1241: }

1243: PetscErrorCode PEPJDSystemSetUp(PEP pep,PetscInt sz,PetscScalar *theta,Vec *u,Vec *p,Vec *ww)
1244: {
1246:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1247:   PEP_JD_PCSHELL *pcctx;
1248: #if !defined(PETSC_USE_COMPLEX)
1249:   PetscScalar    s[2];
1250: #endif

1253:   PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1254:   PEPJDMatSetUp(pep,sz,theta);
1255:   pcctx->u[0] = u[0]; pcctx->u[1] = u[1];
1256:   /* Compute r'. p is a work space vector */
1257:   PEPJDComputeResidual(pep,PETSC_TRUE,sz,u,theta,p,ww);
1258:   PEPJDExtendedPCApply(pjd->pcshell,p[0],pcctx->Bp[0]);
1259:   VecDot(pcctx->Bp[0],u[0],pcctx->gamma);
1260: #if !defined(PETSC_USE_COMPLEX)
1261:   if (sz==2) {
1262:     PEPJDExtendedPCApply(pjd->pcshell,p[1],pcctx->Bp[1]);
1263:     VecDot(pcctx->Bp[0],u[1],pcctx->gamma+1);
1264:     VecMDot(pcctx->Bp[1],2,u,s);
1265:     pcctx->gamma[0] += s[1];
1266:     pcctx->gamma[1] = -pcctx->gamma[1]+s[0];
1267:   }
1268: #endif
1269:   if (sz==1) {
1270:     VecZeroEntries(pcctx->Bp[1]);
1271:     pcctx->gamma[1] = 0.0;
1272:   }
1273:   return(0);
1274: }

1276: PetscErrorCode PEPSolve_JD(PEP pep)
1277: {
1278:   PetscErrorCode  ierr;
1279:   PEP_JD          *pjd = (PEP_JD*)pep->data;
1280:   PetscInt        k,nv,nvc,ld,minv,dim,bupdated=0,sz=1,kspsf=1,idx,off,maxits,nloc;
1281:   PetscMPIInt     np,count;
1282:   PetscScalar     theta[2]={0.0,0.0},ritz[2]={0.0,0.0},*pX,*eig,*eigi,*array;
1283:   PetscReal       norm,*res,tol=0.0,rtol,abstol, dtol;
1284:   PetscBool       lindep,ini=PETSC_TRUE;
1285:   Vec             tc,t[2]={NULL,NULL},u[2]={NULL,NULL},p[2]={NULL,NULL};
1286:   Vec             rc,rr[2],r[2]={NULL,NULL},*ww=pep->work,v[2];
1287:   Mat             G,X,Y;
1288:   KSP             ksp;
1289:   PEP_JD_PCSHELL  *pcctx;
1290:   PEP_JD_MATSHELL *matctx;
1291: #if !defined(PETSC_USE_COMPLEX)
1292:   PetscReal       norm1;
1293: #endif

1296:   PetscCitationsRegister(citation,&cited);
1297:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1298:   BVGetSizes(pep->V,&nloc,NULL,NULL);
1299:   DSGetLeadingDimension(pep->ds,&ld);
1300:   PetscCalloc3(pep->ncv+pep->nev,&eig,pep->ncv+pep->nev,&eigi,pep->ncv+pep->nev,&res);
1301:   pjd->nlock = 0;
1302:   STGetKSP(pep->st,&ksp);
1303:   KSPGetTolerances(ksp,&rtol,&abstol,&dtol,&maxits);
1304: #if !defined (PETSC_USE_COMPLEX)
1305:   kspsf = 2;
1306: #endif
1307:   PEPJDProcessInitialSpace(pep,ww);
1308:   nv = (pep->nini)?pep->nini:1;

1310:   /* Replace preconditioner with one containing projectors */
1311:   PEPJDCreateShellPC(pep,ww);
1312:   PCShellGetContext(pjd->pcshell,(void**)&pcctx);

1314:   /* Create auxiliar vectors */
1315:   BVCreateVec(pjd->V,&u[0]);
1316:   VecDuplicate(u[0],&p[0]);
1317:   VecDuplicate(u[0],&r[0]);
1318: #if !defined (PETSC_USE_COMPLEX)
1319:   VecDuplicate(u[0],&u[1]);
1320:   VecDuplicate(u[0],&p[1]);
1321:   VecDuplicate(u[0],&r[1]);
1322: #endif

1324:   /* Restart loop */
1325:   while (pep->reason == PEP_CONVERGED_ITERATING) {
1326:     pep->its++;
1327:     DSSetDimensions(pep->ds,nv,0,0,0);
1328:     BVSetActiveColumns(pjd->V,bupdated,nv);
1329:     PEPJDUpdateTV(pep,bupdated,nv,ww);
1330:     if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVSetActiveColumns(pjd->W,bupdated,nv); }
1331:     for (k=0;k<pep->nmat;k++) {
1332:       BVSetActiveColumns(pjd->TV[k],bupdated,nv);
1333:       DSGetMat(pep->ds,DSMatExtra[k],&G);
1334:       BVMatProject(pjd->TV[k],NULL,pjd->W,G);
1335:       DSRestoreMat(pep->ds,DSMatExtra[k],&G);
1336:     }
1337:     BVSetActiveColumns(pjd->V,0,nv);
1338:     BVSetActiveColumns(pjd->W,0,nv);

1340:     /* Solve projected problem */
1341:     DSSetState(pep->ds,DS_STATE_RAW);
1342:     DSSolve(pep->ds,pep->eigr,pep->eigi);
1343:     DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1344:     DSSynchronize(pep->ds,pep->eigr,pep->eigi);
1345:     idx = 0;
1346:     do {
1347:       ritz[0] = pep->eigr[idx];
1348: #if !defined(PETSC_USE_COMPLEX)
1349:       ritz[1] = pep->eigi[idx];
1350:       sz = (ritz[1]==0.0)?1:2;
1351: #endif
1352:       /* Compute Ritz vector u=V*X(:,1) */
1353:       DSGetArray(pep->ds,DS_MAT_X,&pX);
1354:       BVSetActiveColumns(pjd->V,0,nv);
1355:       BVMultVec(pjd->V,1.0,0.0,u[0],pX+idx*ld);
1356: #if !defined(PETSC_USE_COMPLEX)
1357:       if (sz==2) {
1358:         BVMultVec(pjd->V,1.0,0.0,u[1],pX+(idx+1)*ld);
1359:       }
1360: #endif
1361:       DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1362:       PEPJDComputeResidual(pep,PETSC_FALSE,sz,u,ritz,r,ww);
1363:       /* Check convergence */
1364:       VecNorm(r[0],NORM_2,&norm);
1365: #if !defined(PETSC_USE_COMPLEX)
1366:       if (sz==2) {
1367:         VecNorm(r[1],NORM_2,&norm1);
1368:         norm = SlepcAbs(norm,norm1);
1369:       }
1370: #endif
1371:       (*pep->converged)(pep,ritz[0],ritz[1],norm,&pep->errest[pep->nconv],pep->convergedctx);
1372:       if (sz==2) pep->errest[pep->nconv+1] = pep->errest[pep->nconv];
1373:       if (ini) {
1374:         tol = PetscMin(.1,pep->errest[pep->nconv]); ini = PETSC_FALSE;
1375:       } else tol = PetscMin(pep->errest[pep->nconv],tol/2);
1376:       (*pep->stopping)(pep,pep->its,pep->max_it,(pep->errest[pep->nconv]<pep->tol)?pep->nconv+sz:pep->nconv,pep->nev,&pep->reason,pep->stoppingctx);
1377:       if (pep->errest[pep->nconv]<pep->tol) {
1378:         /* Ritz pair converged */
1379:         ini = PETSC_TRUE;
1380:         minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1381:         if (pjd->ld>1) {
1382:           BVGetColumn(pjd->X,pep->nconv,&v[0]);
1383:           PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*pep->nconv,pjd->ld-1,0,u[0],PETSC_TRUE);
1384:           BVRestoreColumn(pjd->X,pep->nconv,&v[0]);
1385:           BVSetActiveColumns(pjd->X,0,pep->nconv+1);
1386:           BVNormColumn(pjd->X,pep->nconv,NORM_2,&norm);
1387:           BVScaleColumn(pjd->X,pep->nconv,1.0/norm);
1388:           for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*pep->nconv+k] *= PetscSqrtReal(np)/norm;
1389:           pjd->T[(pep->ncv+1)*pep->nconv] = ritz[0];
1390:           eig[pep->nconv] = ritz[0];
1391:           idx++;
1392: #if !defined(PETSC_USE_COMPLEX)
1393:           if (sz==2) {
1394:             BVGetColumn(pjd->X,pep->nconv+1,&v[0]);
1395:             PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*(pep->nconv+1),pjd->ld-1,0,u[1],PETSC_TRUE);
1396:             BVRestoreColumn(pjd->X,pep->nconv+1,&v[0]);
1397:             BVSetActiveColumns(pjd->X,0,pep->nconv+2);
1398:             BVNormColumn(pjd->X,pep->nconv+1,NORM_2,&norm1);
1399:             BVScaleColumn(pjd->X,pep->nconv+1,1.0/norm1);
1400:             for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*(pep->nconv+1)+k] *= PetscSqrtReal(np)/norm1;
1401:             pjd->T[(pep->ncv+1)*(pep->nconv+1)] = ritz[0];
1402:             pjd->T[(pep->ncv+1)*pep->nconv+1] = -ritz[1]*norm1/norm;
1403:             pjd->T[(pep->ncv+1)*(pep->nconv+1)-1] = ritz[1]*norm/norm1;
1404:             eig[pep->nconv+1] = ritz[0];
1405:             eigi[pep->nconv] = ritz[1]; eigi[pep->nconv+1] = -ritz[1];
1406:             idx++;
1407:           }
1408: #endif
1409:         } else {
1410:           BVInsertVec(pep->V,pep->nconv,u[0]);
1411:         }
1412:         pep->nconv += sz;
1413:       }
1414:     } while (pep->errest[pep->nconv]<pep->tol && pep->nconv<nv);

1416:     if (pep->reason==PEP_CONVERGED_ITERATING) {
1417:       nvc = nv;
1418:       if (idx) {
1419:         pjd->nlock +=idx;
1420:         PEPJDLockConverged(pep,&nv,idx);
1421:       }
1422:       if (nv+sz>=pep->ncv-1) {
1423:         /* Basis full, force restart */
1424:         minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1425:         DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
1426:         DSGetArray(pep->ds,DS_MAT_X,&pX);
1427:         PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1428:         DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1429:         DSGetArray(pep->ds,DS_MAT_Y,&pX);
1430:         PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1431:         DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
1432:         DSGetMat(pep->ds,DS_MAT_X,&X);
1433:         BVMultInPlace(pjd->V,X,0,minv);
1434:         if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1435:          DSGetMat(pep->ds,DS_MAT_Y,&Y);
1436:          BVMultInPlace(pjd->W,Y,pep->nconv,minv);
1437:          DSRestoreMat(pep->ds,DS_MAT_Y,&Y);
1438:         }
1439:         MatDestroy(&X);
1440:         nv = minv;
1441:         bupdated = 0;
1442:       } else {
1443:         if (!idx && pep->errest[pep->nconv]<pjd->fix) {theta[0] = ritz[0]; theta[1] = ritz[1];}
1444:         else {theta[0] = pep->target; theta[1] = 0.0;}
1445:         /* Update system mat */
1446:         PEPJDSystemSetUp(pep,sz,theta,u,p,ww);
1447:         /* Solve correction equation to expand basis */
1448:         BVGetColumn(pjd->V,nv,&t[0]);
1449:         rr[0] = r[0];
1450:         if (sz==2) {
1451:           BVGetColumn(pjd->V,nv+1,&t[1]);
1452:           rr[1] = r[1];
1453:         } else {
1454:           t[1] = NULL;
1455:           rr[1] = NULL;
1456:         }
1457:         VecCreateCompWithVecs(t,kspsf,pjd->vtempl,&tc);
1458:         VecCreateCompWithVecs(rr,kspsf,pjd->vtempl,&rc);
1459:         VecCompSetSubVecs(pjd->vtempl,sz,NULL);
1460:         tol  = PetscMax(rtol,tol/2);
1461:         KSPSetTolerances(ksp,tol,abstol,dtol,maxits);
1462:         KSPSolve(ksp,rc,tc);
1463:         VecDestroy(&tc);
1464:         VecDestroy(&rc);
1465:         VecGetArray(t[0],&array);
1466:         PetscMPIIntCast(pep->nconv,&count);
1467:         MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1468:         VecRestoreArray(t[0],&array);
1469:         BVRestoreColumn(pjd->V,nv,&t[0]);
1470:         BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
1471:         if (lindep || norm==0.0) {
1472:           if (sz==1) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1473:           else off = 1;
1474:         } else {
1475:           off = 0;
1476:           BVScaleColumn(pjd->V,nv,1.0/norm);
1477:         }
1478: #if !defined(PETSC_USE_COMPLEX)
1479:         if (sz==2) {
1480:           VecGetArray(t[1],&array);
1481:           MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1482:           VecRestoreArray(t[1],&array);
1483:           BVRestoreColumn(pjd->V,nv+1,&t[1]);
1484:           if (off) {
1485:             BVCopyColumn(pjd->V,nv+1,nv);
1486:           }
1487:           BVOrthogonalizeColumn(pjd->V,nv+1-off,NULL,&norm,&lindep);
1488:           if (lindep || norm==0.0) {
1489:             if (off) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1490:             else off = 1;
1491:           } else {
1492:             BVScaleColumn(pjd->V,nv+1-off,1.0/norm);
1493:           }
1494:         }
1495: #endif
1496:         if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1497:           BVInsertVec(pjd->W,nv,r[0]);
1498:           if (sz==2 && !off) {
1499:             BVInsertVec(pjd->W,nv+1,r[1]);
1500:           }
1501:           BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
1502:           if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1503:           BVScaleColumn(pjd->W,nv,1.0/norm);
1504:           if (sz==2 && !off) {
1505:             BVOrthogonalizeColumn(pjd->W,nv+1,NULL,&norm,&lindep);
1506:             if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1507:             BVScaleColumn(pjd->W,nv+1,1.0/norm);
1508:           }
1509:         }
1510:         bupdated = idx?0:nv;
1511:         nv += sz-off;
1512:       }
1513:       for (k=0;k<nvc;k++) {
1514:         eig[pep->nconv-idx+k] = pep->eigr[k];
1515: #if !defined(PETSC_USE_COMPLEX)
1516:         eigi[pep->nconv-idx+k] = pep->eigi[k];
1517: #endif
1518:       }
1519:       PEPMonitor(pep,pep->its,pep->nconv,eig,eigi,pep->errest,pep->nconv+1);
1520:     }
1521:   }
1522:   if (pjd->ld>1) {
1523:     for (k=0;k<pep->nconv;k++) {
1524:       pep->eigr[k] = eig[k];
1525:       pep->eigi[k] = eigi[k];
1526:     }
1527:     if (pep->nconv>0) { PEPJDEigenvectors(pep); }
1528:     PetscFree2(pcctx->M,pcctx->ps);
1529:   }
1530:   VecDestroy(&u[0]);
1531:   VecDestroy(&r[0]);
1532:   VecDestroy(&p[0]);
1533: #if !defined (PETSC_USE_COMPLEX)
1534:   VecDestroy(&u[1]);
1535:   VecDestroy(&r[1]);
1536:   VecDestroy(&p[1]);
1537: #endif
1538:   KSPSetTolerances(ksp,rtol,abstol,dtol,maxits);
1539:   KSPSetPC(ksp,pcctx->pc);
1540:   VecDestroy(&pcctx->Bp[0]);
1541:   VecDestroy(&pcctx->Bp[1]);
1542:   MatShellGetContext(pjd->Pshell,(void**)&matctx);
1543:   MatDestroy(&matctx->Pr);
1544:   MatDestroy(&matctx->Pi);
1545:   MatDestroy(&pjd->Pshell);
1546:   MatDestroy(&pcctx->PPr);
1547:   PCDestroy(&pcctx->pc);
1548:   PetscFree(pcctx);
1549:   PetscFree(matctx);
1550:   PCDestroy(&pjd->pcshell);
1551:   PetscFree3(eig,eigi,res);
1552:   VecDestroy(&pjd->vtempl);
1553:   return(0);
1554: }

1556: PetscErrorCode PEPJDSetRestart_JD(PEP pep,PetscReal keep)
1557: {
1558:   PEP_JD *pjd = (PEP_JD*)pep->data;

1561:   if (keep==PETSC_DEFAULT) pjd->keep = 0.5;
1562:   else {
1563:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
1564:     pjd->keep = keep;
1565:   }
1566:   return(0);
1567: }

1569: /*@
1570:    PEPJDSetRestart - Sets the restart parameter for the Jacobi-Davidson
1571:    method, in particular the proportion of basis vectors that must be kept
1572:    after restart.

1574:    Logically Collective on pep

1576:    Input Parameters:
1577: +  pep  - the eigenproblem solver context
1578: -  keep - the number of vectors to be kept at restart

1580:    Options Database Key:
1581: .  -pep_jd_restart - Sets the restart parameter

1583:    Notes:
1584:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

1586:    Level: advanced

1588: .seealso: PEPJDGetRestart()
1589: @*/
1590: PetscErrorCode PEPJDSetRestart(PEP pep,PetscReal keep)
1591: {

1597:   PetscTryMethod(pep,"PEPJDSetRestart_C",(PEP,PetscReal),(pep,keep));
1598:   return(0);
1599: }

1601: PetscErrorCode PEPJDGetRestart_JD(PEP pep,PetscReal *keep)
1602: {
1603:   PEP_JD *pjd = (PEP_JD*)pep->data;

1606:   *keep = pjd->keep;
1607:   return(0);
1608: }

1610: /*@
1611:    PEPJDGetRestart - Gets the restart parameter used in the Jacobi-Davidson method.

1613:    Not Collective

1615:    Input Parameter:
1616: .  pep - the eigenproblem solver context

1618:    Output Parameter:
1619: .  keep - the restart parameter

1621:    Level: advanced

1623: .seealso: PEPJDSetRestart()
1624: @*/
1625: PetscErrorCode PEPJDGetRestart(PEP pep,PetscReal *keep)
1626: {

1632:   PetscUseMethod(pep,"PEPJDGetRestart_C",(PEP,PetscReal*),(pep,keep));
1633:   return(0);
1634: }

1636: PetscErrorCode PEPJDSetFix_JD(PEP pep,PetscReal fix)
1637: {
1638:   PEP_JD *pjd = (PEP_JD*)pep->data;

1641:   if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) pjd->fix = 0.01;
1642:   else {
1643:     if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
1644:     pjd->fix = fix;
1645:   }
1646:   return(0);
1647: }

1649: /*@
1650:    PEPJDSetFix - Sets the threshold for changing the target in the correction
1651:    equation.

1653:    Logically Collective on pep

1655:    Input Parameters:
1656: +  pep - the eigenproblem solver context
1657: -  fix - threshold for changing the target

1659:    Options Database Key:
1660: .  -pep_jd_fix - the fix value

1662:    Note:
1663:    The target in the correction equation is fixed at the first iterations.
1664:    When the norm of the residual vector is lower than the fix value,
1665:    the target is set to the corresponding eigenvalue.

1667:    Level: advanced

1669: .seealso: PEPJDGetFix()
1670: @*/
1671: PetscErrorCode PEPJDSetFix(PEP pep,PetscReal fix)
1672: {

1678:   PetscTryMethod(pep,"PEPJDSetFix_C",(PEP,PetscReal),(pep,fix));
1679:   return(0);
1680: }

1682: PetscErrorCode PEPJDGetFix_JD(PEP pep,PetscReal *fix)
1683: {
1684:   PEP_JD *pjd = (PEP_JD*)pep->data;

1687:   *fix = pjd->fix;
1688:   return(0);
1689: }

1691: /*@
1692:    PEPJDGetFix - Returns the threshold for changing the target in the correction
1693:    equation.

1695:    Not Collective

1697:    Input Parameter:
1698: .  pep - the eigenproblem solver context

1700:    Output Parameter:
1701: .  fix - threshold for changing the target

1703:    Note:
1704:    The target in the correction equation is fixed at the first iterations.
1705:    When the norm of the residual vector is lower than the fix value,
1706:    the target is set to the corresponding eigenvalue.

1708:    Level: advanced

1710: .seealso: PEPJDSetFix()
1711: @*/
1712: PetscErrorCode PEPJDGetFix(PEP pep,PetscReal *fix)
1713: {

1719:   PetscUseMethod(pep,"PEPJDGetFix_C",(PEP,PetscReal*),(pep,fix));
1720:   return(0);
1721: }

1723: PetscErrorCode PEPJDSetReusePreconditioner_JD(PEP pep,PetscBool reusepc)
1724: {
1725:   PEP_JD *pjd = (PEP_JD*)pep->data;

1728:   pjd->reusepc = reusepc;
1729:   return(0);
1730: }

1732: /*@
1733:    PEPJDSetReusePreconditioner - Sets a flag indicating whether the preconditioner
1734:    must be reused or not.

1736:    Logically Collective on pep

1738:    Input Parameters:
1739: +  pep     - the eigenproblem solver context
1740: -  reusepc - the reuse flag

1742:    Options Database Key:
1743: .  -pep_jd_reuse_preconditioner - the reuse flag

1745:    Note:
1746:    The default value is False. If set to True, the preconditioner is built
1747:    only at the beginning, using the target value. Otherwise, it may be rebuilt
1748:    (depending on the fix parameter) at each iteration from the Ritz value.

1750:    Level: advanced

1752: .seealso: PEPJDGetReusePreconditioner(), PEPJDSetFix()
1753: @*/
1754: PetscErrorCode PEPJDSetReusePreconditioner(PEP pep,PetscBool reusepc)
1755: {

1761:   PetscTryMethod(pep,"PEPJDSetReusePreconditioner_C",(PEP,PetscBool),(pep,reusepc));
1762:   return(0);
1763: }

1765: PetscErrorCode PEPJDGetReusePreconditioner_JD(PEP pep,PetscBool *reusepc)
1766: {
1767:   PEP_JD *pjd = (PEP_JD*)pep->data;

1770:   *reusepc = pjd->reusepc;
1771:   return(0);
1772: }

1774: /*@
1775:    PEPJDGetReusePreconditioner - Returns the flag for reusing the preconditioner.

1777:    Not Collective

1779:    Input Parameter:
1780: .  pep - the eigenproblem solver context

1782:    Output Parameter:
1783: .  reusepc - the reuse flag

1785:    Level: advanced

1787: .seealso: PEPJDSetReusePreconditioner()
1788: @*/
1789: PetscErrorCode PEPJDGetReusePreconditioner(PEP pep,PetscBool *reusepc)
1790: {

1796:   PetscUseMethod(pep,"PEPJDGetReusePreconditioner_C",(PEP,PetscBool*),(pep,reusepc));
1797:   return(0);
1798: }

1800: PetscErrorCode PEPJDSetMinimalityIndex_JD(PEP pep,PetscInt mmidx)
1801: {
1802:   PEP_JD *pjd = (PEP_JD*)pep->data;

1805:   if (mmidx == PETSC_DEFAULT || mmidx == PETSC_DECIDE) pjd->mmidx = 1;
1806:   else {
1807:     if (mmidx < 1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mmidx value");
1808:     pjd->mmidx = mmidx;
1809:     pep->state = PEP_STATE_INITIAL;
1810:   }
1811:   return(0);
1812: }

1814: /*@
1815:    PEPJDSetMinimalityIndex - Sets the maximum allowed value for the minimality index.

1817:    Logically Collective on pep

1819:    Input Parameters:
1820: +  pep   - the eigenproblem solver context
1821: -  mmidx - maximum minimality index

1823:    Options Database Key:
1824: .  -pep_jd_minimality_index - the minimality index value

1826:    Note:
1827:    The default value is equal to the degree of the polynomial. A smaller value
1828:    can be used if the wanted eigenvectors are known to be linearly independent.

1830:    Level: advanced

1832: .seealso: PEPJDGetMinimalityIndex()
1833: @*/
1834: PetscErrorCode PEPJDSetMinimalityIndex(PEP pep,PetscInt mmidx)
1835: {

1841:   PetscTryMethod(pep,"PEPJDSetMinimalityIndex_C",(PEP,PetscInt),(pep,mmidx));
1842:   return(0);
1843: }

1845: PetscErrorCode PEPJDGetMinimalityIndex_JD(PEP pep,PetscInt *mmidx)
1846: {
1847:   PEP_JD *pjd = (PEP_JD*)pep->data;

1850:   *mmidx = pjd->mmidx;
1851:   return(0);
1852: }

1854: /*@
1855:    PEPJDGetMinimalityIndex - Returns the maximum allowed value of the minimality
1856:    index.

1858:    Not Collective

1860:    Input Parameter:
1861: .  pep - the eigenproblem solver context

1863:    Output Parameter:
1864: .  mmidx - minimality index

1866:    Level: advanced

1868: .seealso: PEPJDSetMinimalityIndex()
1869: @*/
1870: PetscErrorCode PEPJDGetMinimalityIndex(PEP pep,PetscInt *mmidx)
1871: {

1877:   PetscUseMethod(pep,"PEPJDGetMinimalityIndex_C",(PEP,PetscInt*),(pep,mmidx));
1878:   return(0);
1879: }

1881: PetscErrorCode PEPJDSetProjection_JD(PEP pep,PEPJDProjection proj)
1882: {
1883:   PEP_JD *pjd = (PEP_JD*)pep->data;

1886:   switch (proj) {
1887:     case PEP_JD_PROJECTION_HARMONIC:
1888:     case PEP_JD_PROJECTION_ORTHOGONAL:
1889:       if (pjd->proj != proj) {
1890:         pep->state = PEP_STATE_INITIAL;
1891:         pjd->proj = proj;
1892:       }
1893:       break;
1894:     default:
1895:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'proj' value");
1896:   }
1897:   return(0);
1898: }

1900: /*@
1901:    PEPJDSetProjection - Sets the type of projection to be used in the Jacobi-Davidson solver.

1903:    Logically Collective on pep

1905:    Input Parameters:
1906: +  pep  - the eigenproblem solver context
1907: -  proj - the type of projection

1909:    Options Database Key:
1910: .  -pep_jd_projection - the projection type, either orthogonal or harmonic

1912:    Level: advanced

1914: .seealso: PEPJDGetProjection()
1915: @*/
1916: PetscErrorCode PEPJDSetProjection(PEP pep,PEPJDProjection proj)
1917: {

1923:   PetscTryMethod(pep,"PEPJDSetProjection_C",(PEP,PEPJDProjection),(pep,proj));
1924:   return(0);
1925: }

1927: PetscErrorCode PEPJDGetProjection_JD(PEP pep,PEPJDProjection *proj)
1928: {
1929:   PEP_JD *pjd = (PEP_JD*)pep->data;

1932:   *proj = pjd->proj;
1933:   return(0);
1934: }

1936: /*@
1937:    PEPJDGetProjection - Returns the type of projection used by the Jacobi-Davidson solver.

1939:    Not Collective

1941:    Input Parameter:
1942: .  pep - the eigenproblem solver context

1944:    Output Parameter:
1945: .  proj - the type of projection

1947:    Level: advanced

1949: .seealso: PEPJDSetProjection()
1950: @*/
1951: PetscErrorCode PEPJDGetProjection(PEP pep,PEPJDProjection *proj)
1952: {

1958:   PetscUseMethod(pep,"PEPJDGetProjection_C",(PEP,PEPJDProjection*),(pep,proj));
1959:   return(0);
1960: }

1962: PetscErrorCode PEPSetFromOptions_JD(PetscOptionItems *PetscOptionsObject,PEP pep)
1963: {
1964:   PetscErrorCode  ierr;
1965:   PetscBool       flg,b1;
1966:   PetscReal       r1;
1967:   PetscInt        i1;
1968:   PEPJDProjection proj;

1971:   PetscOptionsHead(PetscOptionsObject,"PEP JD Options");

1973:     PetscOptionsReal("-pep_jd_restart","Proportion of vectors kept after restart","PEPJDSetRestart",0.5,&r1,&flg);
1974:     if (flg) { PEPJDSetRestart(pep,r1); }

1976:     PetscOptionsReal("-pep_jd_fix","Tolerance for changing the target in the correction equation","PEPJDSetFix",0.01,&r1,&flg);
1977:     if (flg) { PEPJDSetFix(pep,r1); }

1979:     PetscOptionsBool("-pep_jd_reuse_preconditioner","Whether to reuse the preconditioner","PEPJDSetReusePreconditoiner",PETSC_FALSE,&b1,&flg);
1980:     if (flg) { PEPJDSetReusePreconditioner(pep,b1); }

1982:     PetscOptionsInt("-pep_jd_minimality_index","Maximum allowed minimality index","PEPJDSetMinimalityIndex",1,&i1,&flg);
1983:     if (flg) { PEPJDSetMinimalityIndex(pep,i1); }

1985:     PetscOptionsEnum("-pep_jd_projection","Type of projection","PEPJDSetProjection",PEPJDProjectionTypes,(PetscEnum)PEP_JD_PROJECTION_HARMONIC,(PetscEnum*)&proj,&flg);
1986:     if (flg) { PEPJDSetProjection(pep,proj); }

1988:   PetscOptionsTail();
1989:   return(0);
1990: }

1992: PetscErrorCode PEPView_JD(PEP pep,PetscViewer viewer)
1993: {
1995:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1996:   PetscBool      isascii;

1999:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
2000:   if (isascii) {
2001:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*pjd->keep));
2002:     PetscViewerASCIIPrintf(viewer,"  threshold for changing the target in the correction equation (fix): %g\n",(double)pjd->fix);
2003:     PetscViewerASCIIPrintf(viewer,"  projection type: %s\n",PEPJDProjectionTypes[pjd->proj]);
2004:     PetscViewerASCIIPrintf(viewer,"  maximum allowed minimality index: %d\n",pjd->mmidx);
2005:     if (pjd->reusepc) { PetscViewerASCIIPrintf(viewer,"  reusing the preconditioner\n"); }
2006:   }
2007:   return(0);
2008: }

2010: PetscErrorCode PEPSetDefaultST_JD(PEP pep)
2011: {
2013:   KSP            ksp;

2016:   if (!((PetscObject)pep->st)->type_name) {
2017:     STSetType(pep->st,STPRECOND);
2018:     STPrecondSetKSPHasMat(pep->st,PETSC_TRUE);
2019:   }
2020:   STSetTransform(pep->st,PETSC_FALSE);
2021:   STGetKSP(pep->st,&ksp);
2022:   if (!((PetscObject)ksp)->type_name) {
2023:     KSPSetType(ksp,KSPBCGSL);
2024:     KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
2025:   }
2026:   return(0);
2027: }

2029: PetscErrorCode PEPReset_JD(PEP pep)
2030: {
2032:   PEP_JD         *pjd = (PEP_JD*)pep->data;
2033:   PetscInt       i;

2036:   for (i=0;i<pep->nmat;i++) {
2037:     BVDestroy(pjd->TV+i);
2038:   }
2039:   if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVDestroy(&pjd->W); }
2040:   if (pjd->ld>1) {
2041:     BVDestroy(&pjd->V);
2042:     for (i=0;i<pep->nmat;i++) {
2043:       BVDestroy(pjd->AX+i);
2044:     }
2045:     BVDestroy(&pjd->N[0]);
2046:     BVDestroy(&pjd->N[1]);
2047:     PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
2048:   }
2049:   PetscFree2(pjd->TV,pjd->AX);
2050:   return(0);
2051: }

2053: PetscErrorCode PEPDestroy_JD(PEP pep)
2054: {

2058:   PetscFree(pep->data);
2059:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
2060:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
2061:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",NULL);
2062:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",NULL);
2063:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",NULL);
2064:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",NULL);
2065:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",NULL);
2066:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",NULL);
2067:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",NULL);
2068:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",NULL);
2069:   return(0);
2070: }

2072: SLEPC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
2073: {
2074:   PEP_JD         *pjd;

2078:   PetscNewLog(pep,&pjd);
2079:   pep->data = (void*)pjd;

2081:   pep->lineariz = PETSC_FALSE;
2082:   pjd->fix      = 0.01;
2083:   pjd->mmidx    = 0;

2085:   pep->ops->solve          = PEPSolve_JD;
2086:   pep->ops->setup          = PEPSetUp_JD;
2087:   pep->ops->setfromoptions = PEPSetFromOptions_JD;
2088:   pep->ops->destroy        = PEPDestroy_JD;
2089:   pep->ops->reset          = PEPReset_JD;
2090:   pep->ops->view           = PEPView_JD;
2091:   pep->ops->setdefaultst   = PEPSetDefaultST_JD;

2093:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
2094:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
2095:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",PEPJDSetFix_JD);
2096:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",PEPJDGetFix_JD);
2097:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",PEPJDSetReusePreconditioner_JD);
2098:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",PEPJDGetReusePreconditioner_JD);
2099:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",PEPJDSetMinimalityIndex_JD);
2100:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",PEPJDGetMinimalityIndex_JD);
2101:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",PEPJDSetProjection_JD);
2102:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",PEPJDGetProjection_JD);
2103:   return(0);
2104: }