Actual source code: pjd.c

slepc-3.13.4 2020-09-02
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. (in press), 2019.

 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 = \"IP\",\n"
 45:   "   number = \"-\",\n"
 46:   "   pages = \"1--24\",\n"
 47:   "   year = \"2019,\"\n"
 48:   "   doi = \"https://doi.org/10.1007/s10543-019-00778-z\"\n"
 49:   "}\n";

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

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

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

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

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

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

136:   pep->lineariz = PETSC_FALSE;
137:   PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
138:   if (pep->max_it==PETSC_DEFAULT) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
139:   if (!pep->which) pep->which = PEP_TARGET_MAGNITUDE;
140:   if (pep->which!=PEP_TARGET_MAGNITUDE && pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Wrong value of pep->which");

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

145:   STGetTransform(pep->st,&flg);
146:   if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

1121: static PetscErrorCode PEPJDEigenvectors(PEP pep)
1122: {
1124:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1125:   PetscBLASInt   ld,nconv,info,nc;
1126:   PetscScalar    *Z,*w;
1127:   PetscReal      *wr,norm;
1128:   PetscInt       i;
1129:   Mat            U;
1130: #if !defined(PETSC_USE_COMPLEX)
1131:   Vec            v,v1;
1132: #endif

1135:   PetscMalloc3(pep->nconv*pep->nconv,&Z,3*pep->ncv,&wr,2*pep->ncv,&w);
1136:   PetscBLASIntCast(pep->ncv,&ld);
1137:   PetscBLASIntCast(pep->nconv,&nconv);
1138: #if !defined(PETSC_USE_COMPLEX)
1139:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
1140: #else
1141:   PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
1142: #endif
1143:   SlepcCheckLapackInfo("trevc",info);
1144:   MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
1145:   BVSetActiveColumns(pjd->X,0,pep->nconv);
1146:   BVMultInPlace(pjd->X,U,0,pep->nconv);
1147:   for (i=0;i<pep->nconv;i++) {
1148: #if !defined(PETSC_USE_COMPLEX)
1149:     if (pep->eigi[i]!=0.0) {   /* first eigenvalue of a complex conjugate pair */
1150:       BVGetColumn(pjd->X,i,&v);
1151:       BVGetColumn(pjd->X,i+1,&v1);
1152:       VecNormalizeComplex(v,v1,PETSC_TRUE,NULL);
1153:       BVRestoreColumn(pjd->X,i,&v);
1154:       BVRestoreColumn(pjd->X,i+1,&v1);
1155:       i++;
1156:     } else   /* real eigenvalue */
1157: #endif
1158:     {
1159:       BVNormColumn(pjd->X,i,NORM_2,&norm);
1160:       BVScaleColumn(pjd->X,i,1.0/norm);
1161:     }
1162:   }
1163:   MatDestroy(&U);
1164:   PetscFree3(Z,wr,w);
1165:   return(0);
1166: }

1168: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv,PetscInt sz)
1169: {
1171:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1172:   PetscInt       j,i,*P,ldds,rk=0,nvv=*nv;
1173:   Vec            v,x,w;
1174:   PetscScalar    *R,*r,*pX,target[2];
1175:   Mat            X;
1176:   PetscBLASInt   sz_,rk_,nv_,info;
1177:   PetscMPIInt    np;

1180:   /* update AX and XpX */
1181:   for (i=sz;i>0;i--) {
1182:     BVGetColumn(pjd->X,pjd->nlock-i,&x);
1183:     for (j=0;j<pep->nmat;j++) {
1184:       BVGetColumn(pjd->AX[j],pjd->nlock-i,&v);
1185:       MatMult(pep->A[j],x,v);
1186:       BVRestoreColumn(pjd->AX[j],pjd->nlock-i,&v);
1187:       BVSetActiveColumns(pjd->AX[j],0,pjd->nlock-i+1);
1188:     }
1189:     BVRestoreColumn(pjd->X,pjd->nlock-i,&x);
1190:     BVDotColumn(pjd->X,(pjd->nlock-i),pjd->XpX+(pjd->nlock-i)*(pjd->ld));
1191:     pjd->XpX[(pjd->nlock-i)*(1+pjd->ld)] = 1.0;
1192:     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]);
1193:   }

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

1198:   /* evaluate the polynomial basis in T */
1199:   PetscArrayzero(pjd->Tj,pjd->ld*pjd->ld*pep->nmat);
1200:   for (j=0;j<pep->nmat;j++) {
1201:     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);
1202:   }

1204:   /* Extend search space */
1205:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1206:   PetscCalloc3(nvv,&P,nvv*nvv,&R,nvv*sz,&r);
1207:   DSGetLeadingDimension(pep->ds,&ldds);
1208:   DSGetArray(pep->ds,DS_MAT_X,&pX);
1209:   PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
1210:   for (j=0;j<sz;j++) {
1211:     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 */
1212:   }
1213:   PetscBLASIntCast(rk,&rk_);
1214:   PetscBLASIntCast(sz,&sz_);
1215:   PetscBLASIntCast(nvv,&nv_);
1216:   PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
1217:   SlepcCheckLapackInfo("trtri",info);
1218:   for (i=0;i<sz;i++) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r+i,&sz_));
1219:   for (i=0;i<sz*rk;i++) r[i] = PetscConj(r[i])/PetscSqrtReal(np); /* revert */
1220:   BVSetActiveColumns(pjd->V,0,nvv);
1221:   rk -= sz;
1222:   for (j=0;j<rk;j++) {
1223:     PetscArraycpy(R+j*nvv,pX+(j+sz)*ldds,nvv);
1224:   }
1225:   DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1226:   MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk,R,&X);
1227:   BVMultInPlace(pjd->V,X,0,rk);
1228:   MatDestroy(&X);
1229:   BVSetActiveColumns(pjd->V,0,rk);
1230:   for (j=0;j<rk;j++) {
1231:     BVGetColumn(pjd->V,j,&v);
1232:     PEPJDCopyToExtendedVec(pep,NULL,r+sz*(j+sz),sz,pjd->nlock-sz,v,PETSC_FALSE);
1233:     BVRestoreColumn(pjd->V,j,&v);
1234:   }
1235:   BVOrthogonalize(pjd->V,NULL);

1237:   if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1238:     for (j=0;j<rk;j++) {
1239:       /* W = P(target)*V */
1240:       BVGetColumn(pjd->W,j,&w);
1241:       BVGetColumn(pjd->V,j,&v);
1242:       target[0] = pep->target; target[1] = 0.0;
1243:       PEPJDComputeResidual(pep,PETSC_FALSE,1,&v,target,&w,pep->work);
1244:       BVRestoreColumn(pjd->V,j,&v);
1245:       BVRestoreColumn(pjd->W,j,&w);
1246:     }
1247:     BVSetActiveColumns(pjd->W,0,rk);
1248:     BVOrthogonalize(pjd->W,NULL);
1249:   }
1250:   *nv = rk;
1251:   PetscFree3(P,R,r);
1252:   return(0);
1253: }

1255: PetscErrorCode PEPJDSystemSetUp(PEP pep,PetscInt sz,PetscScalar *theta,Vec *u,Vec *p,Vec *ww)
1256: {
1258:   PEP_JD         *pjd = (PEP_JD*)pep->data;
1259:   PEP_JD_PCSHELL *pcctx;
1260: #if !defined(PETSC_USE_COMPLEX)
1261:   PetscScalar    s[2];
1262: #endif

1265:   PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1266:   PEPJDMatSetUp(pep,sz,theta);
1267:   pcctx->u[0] = u[0]; pcctx->u[1] = u[1];
1268:   /* Compute r'. p is a work space vector */
1269:   PEPJDComputeResidual(pep,PETSC_TRUE,sz,u,theta,p,ww);
1270:   PEPJDExtendedPCApply(pjd->pcshell,p[0],pcctx->Bp[0]);
1271:   VecDot(pcctx->Bp[0],u[0],pcctx->gamma);
1272: #if !defined(PETSC_USE_COMPLEX)
1273:   if (sz==2) {
1274:     PEPJDExtendedPCApply(pjd->pcshell,p[1],pcctx->Bp[1]);
1275:     VecDot(pcctx->Bp[0],u[1],pcctx->gamma+1);
1276:     VecMDot(pcctx->Bp[1],2,u,s);
1277:     pcctx->gamma[0] += s[1];
1278:     pcctx->gamma[1] = -pcctx->gamma[1]+s[0];
1279:   }
1280: #endif
1281:   if (sz==1) {
1282:     VecZeroEntries(pcctx->Bp[1]);
1283:     pcctx->gamma[1] = 0.0;
1284:   }
1285:   return(0);
1286: }

1288: PetscErrorCode PEPSolve_JD(PEP pep)
1289: {
1290:   PetscErrorCode  ierr;
1291:   PEP_JD          *pjd = (PEP_JD*)pep->data;
1292:   PetscInt        k,nv,nvc,ld,minv,dim,bupdated=0,sz=1,kspsf=1,idx,off,maxits,nloc;
1293:   PetscMPIInt     np,count;
1294:   PetscScalar     theta[2]={0.0,0.0},ritz[2]={0.0,0.0},*pX,*eig,*eigi,*array;
1295:   PetscReal       norm,*res,tol=0.0,rtol,abstol, dtol;
1296:   PetscBool       lindep,ini=PETSC_TRUE;
1297:   Vec             tc,t[2]={NULL,NULL},u[2]={NULL,NULL},p[2]={NULL,NULL};
1298:   Vec             rc,rr[2],r[2]={NULL,NULL},*ww=pep->work,v[2];
1299:   Mat             G,X,Y;
1300:   KSP             ksp;
1301:   PEP_JD_PCSHELL  *pcctx;
1302:   PEP_JD_MATSHELL *matctx;
1303: #if !defined(PETSC_USE_COMPLEX)
1304:   PetscReal       norm1;
1305: #endif

1308:   PetscCitationsRegister(citation,&cited);
1309:   MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1310:   BVGetSizes(pep->V,&nloc,NULL,NULL);
1311:   DSGetLeadingDimension(pep->ds,&ld);
1312:   PetscCalloc3(pep->ncv+pep->nev,&eig,pep->ncv+pep->nev,&eigi,pep->ncv+pep->nev,&res);
1313:   pjd->nlock = 0;
1314:   STGetKSP(pep->st,&ksp);
1315:   KSPGetTolerances(ksp,&rtol,&abstol,&dtol,&maxits);
1316: #if !defined (PETSC_USE_COMPLEX)
1317:   kspsf = 2;
1318: #endif
1319:   PEPJDProcessInitialSpace(pep,ww);
1320:   nv = (pep->nini)?pep->nini:1;

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

1326:   /* Create auxiliar vectors */
1327:   BVCreateVec(pjd->V,&u[0]);
1328:   VecDuplicate(u[0],&p[0]);
1329:   VecDuplicate(u[0],&r[0]);
1330: #if !defined (PETSC_USE_COMPLEX)
1331:   VecDuplicate(u[0],&u[1]);
1332:   VecDuplicate(u[0],&p[1]);
1333:   VecDuplicate(u[0],&r[1]);
1334: #endif

1336:   /* Restart loop */
1337:   while (pep->reason == PEP_CONVERGED_ITERATING) {
1338:     pep->its++;
1339:     DSSetDimensions(pep->ds,nv,0,0,0);
1340:     BVSetActiveColumns(pjd->V,bupdated,nv);
1341:     PEPJDUpdateTV(pep,bupdated,nv,ww);
1342:     if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVSetActiveColumns(pjd->W,bupdated,nv); }
1343:     for (k=0;k<pep->nmat;k++) {
1344:       BVSetActiveColumns(pjd->TV[k],bupdated,nv);
1345:       DSGetMat(pep->ds,DSMatExtra[k],&G);
1346:       BVMatProject(pjd->TV[k],NULL,pjd->W,G);
1347:       DSRestoreMat(pep->ds,DSMatExtra[k],&G);
1348:     }
1349:     BVSetActiveColumns(pjd->V,0,nv);
1350:     BVSetActiveColumns(pjd->W,0,nv);

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

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

1568: PetscErrorCode PEPJDSetRestart_JD(PEP pep,PetscReal keep)
1569: {
1570:   PEP_JD *pjd = (PEP_JD*)pep->data;

1573:   if (keep==PETSC_DEFAULT) pjd->keep = 0.5;
1574:   else {
1575:     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]");
1576:     pjd->keep = keep;
1577:   }
1578:   return(0);
1579: }

1581: /*@
1582:    PEPJDSetRestart - Sets the restart parameter for the Jacobi-Davidson
1583:    method, in particular the proportion of basis vectors that must be kept
1584:    after restart.

1586:    Logically Collective on pep

1588:    Input Parameters:
1589: +  pep  - the eigenproblem solver context
1590: -  keep - the number of vectors to be kept at restart

1592:    Options Database Key:
1593: .  -pep_jd_restart - Sets the restart parameter

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

1598:    Level: advanced

1600: .seealso: PEPJDGetRestart()
1601: @*/
1602: PetscErrorCode PEPJDSetRestart(PEP pep,PetscReal keep)
1603: {

1609:   PetscTryMethod(pep,"PEPJDSetRestart_C",(PEP,PetscReal),(pep,keep));
1610:   return(0);
1611: }

1613: PetscErrorCode PEPJDGetRestart_JD(PEP pep,PetscReal *keep)
1614: {
1615:   PEP_JD *pjd = (PEP_JD*)pep->data;

1618:   *keep = pjd->keep;
1619:   return(0);
1620: }

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

1625:    Not Collective

1627:    Input Parameter:
1628: .  pep - the eigenproblem solver context

1630:    Output Parameter:
1631: .  keep - the restart parameter

1633:    Level: advanced

1635: .seealso: PEPJDSetRestart()
1636: @*/
1637: PetscErrorCode PEPJDGetRestart(PEP pep,PetscReal *keep)
1638: {

1644:   PetscUseMethod(pep,"PEPJDGetRestart_C",(PEP,PetscReal*),(pep,keep));
1645:   return(0);
1646: }

1648: PetscErrorCode PEPJDSetFix_JD(PEP pep,PetscReal fix)
1649: {
1650:   PEP_JD *pjd = (PEP_JD*)pep->data;

1653:   if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) pjd->fix = 0.01;
1654:   else {
1655:     if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
1656:     pjd->fix = fix;
1657:   }
1658:   return(0);
1659: }

1661: /*@
1662:    PEPJDSetFix - Sets the threshold for changing the target in the correction
1663:    equation.

1665:    Logically Collective on pep

1667:    Input Parameters:
1668: +  pep - the eigenproblem solver context
1669: -  fix - threshold for changing the target

1671:    Options Database Key:
1672: .  -pep_jd_fix - the fix value

1674:    Note:
1675:    The target in the correction equation is fixed at the first iterations.
1676:    When the norm of the residual vector is lower than the fix value,
1677:    the target is set to the corresponding eigenvalue.

1679:    Level: advanced

1681: .seealso: PEPJDGetFix()
1682: @*/
1683: PetscErrorCode PEPJDSetFix(PEP pep,PetscReal fix)
1684: {

1690:   PetscTryMethod(pep,"PEPJDSetFix_C",(PEP,PetscReal),(pep,fix));
1691:   return(0);
1692: }

1694: PetscErrorCode PEPJDGetFix_JD(PEP pep,PetscReal *fix)
1695: {
1696:   PEP_JD *pjd = (PEP_JD*)pep->data;

1699:   *fix = pjd->fix;
1700:   return(0);
1701: }

1703: /*@
1704:    PEPJDGetFix - Returns the threshold for changing the target in the correction
1705:    equation.

1707:    Not Collective

1709:    Input Parameter:
1710: .  pep - the eigenproblem solver context

1712:    Output Parameter:
1713: .  fix - threshold for changing the target

1715:    Note:
1716:    The target in the correction equation is fixed at the first iterations.
1717:    When the norm of the residual vector is lower than the fix value,
1718:    the target is set to the corresponding eigenvalue.

1720:    Level: advanced

1722: .seealso: PEPJDSetFix()
1723: @*/
1724: PetscErrorCode PEPJDGetFix(PEP pep,PetscReal *fix)
1725: {

1731:   PetscUseMethod(pep,"PEPJDGetFix_C",(PEP,PetscReal*),(pep,fix));
1732:   return(0);
1733: }

1735: PetscErrorCode PEPJDSetReusePreconditioner_JD(PEP pep,PetscBool reusepc)
1736: {
1737:   PEP_JD *pjd = (PEP_JD*)pep->data;

1740:   pjd->reusepc = reusepc;
1741:   return(0);
1742: }

1744: /*@
1745:    PEPJDSetReusePreconditioner - Sets a flag indicating whether the preconditioner
1746:    must be reused or not.

1748:    Logically Collective on pep

1750:    Input Parameters:
1751: +  pep     - the eigenproblem solver context
1752: -  reusepc - the reuse flag

1754:    Options Database Key:
1755: .  -pep_jd_reuse_preconditioner - the reuse flag

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

1762:    Level: advanced

1764: .seealso: PEPJDGetReusePreconditioner(), PEPJDSetFix()
1765: @*/
1766: PetscErrorCode PEPJDSetReusePreconditioner(PEP pep,PetscBool reusepc)
1767: {

1773:   PetscTryMethod(pep,"PEPJDSetReusePreconditioner_C",(PEP,PetscBool),(pep,reusepc));
1774:   return(0);
1775: }

1777: PetscErrorCode PEPJDGetReusePreconditioner_JD(PEP pep,PetscBool *reusepc)
1778: {
1779:   PEP_JD *pjd = (PEP_JD*)pep->data;

1782:   *reusepc = pjd->reusepc;
1783:   return(0);
1784: }

1786: /*@
1787:    PEPJDGetReusePreconditioner - Returns the flag for reusing the preconditioner.

1789:    Not Collective

1791:    Input Parameter:
1792: .  pep - the eigenproblem solver context

1794:    Output Parameter:
1795: .  reusepc - the reuse flag

1797:    Level: advanced

1799: .seealso: PEPJDSetReusePreconditioner()
1800: @*/
1801: PetscErrorCode PEPJDGetReusePreconditioner(PEP pep,PetscBool *reusepc)
1802: {

1808:   PetscUseMethod(pep,"PEPJDGetReusePreconditioner_C",(PEP,PetscBool*),(pep,reusepc));
1809:   return(0);
1810: }

1812: PetscErrorCode PEPJDSetMinimalityIndex_JD(PEP pep,PetscInt mmidx)
1813: {
1814:   PEP_JD *pjd = (PEP_JD*)pep->data;

1817:   if (mmidx == PETSC_DEFAULT || mmidx == PETSC_DECIDE) pjd->mmidx = 1;
1818:   else {
1819:     if (mmidx < 1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mmidx value");
1820:     pjd->mmidx = mmidx;
1821:     pep->state = PEP_STATE_INITIAL;
1822:   }
1823:   return(0);
1824: }

1826: /*@
1827:    PEPJDSetMinimalityIndex - Sets the maximum allowed value for the minimality index.

1829:    Logically Collective on pep

1831:    Input Parameters:
1832: +  pep   - the eigenproblem solver context
1833: -  mmidx - maximum minimality index

1835:    Options Database Key:
1836: .  -pep_jd_minimality_index - the minimality index value

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

1842:    Level: advanced

1844: .seealso: PEPJDGetMinimalityIndex()
1845: @*/
1846: PetscErrorCode PEPJDSetMinimalityIndex(PEP pep,PetscInt mmidx)
1847: {

1853:   PetscTryMethod(pep,"PEPJDSetMinimalityIndex_C",(PEP,PetscInt),(pep,mmidx));
1854:   return(0);
1855: }

1857: PetscErrorCode PEPJDGetMinimalityIndex_JD(PEP pep,PetscInt *mmidx)
1858: {
1859:   PEP_JD *pjd = (PEP_JD*)pep->data;

1862:   *mmidx = pjd->mmidx;
1863:   return(0);
1864: }

1866: /*@
1867:    PEPJDGetMinimalityIndex - Returns the maximum allowed value of the minimality
1868:    index.

1870:    Not Collective

1872:    Input Parameter:
1873: .  pep - the eigenproblem solver context

1875:    Output Parameter:
1876: .  mmidx - minimality index

1878:    Level: advanced

1880: .seealso: PEPJDSetMinimalityIndex()
1881: @*/
1882: PetscErrorCode PEPJDGetMinimalityIndex(PEP pep,PetscInt *mmidx)
1883: {

1889:   PetscUseMethod(pep,"PEPJDGetMinimalityIndex_C",(PEP,PetscInt*),(pep,mmidx));
1890:   return(0);
1891: }

1893: PetscErrorCode PEPJDSetProjection_JD(PEP pep,PEPJDProjection proj)
1894: {
1895:   PEP_JD *pjd = (PEP_JD*)pep->data;

1898:   switch (proj) {
1899:     case PEP_JD_PROJECTION_HARMONIC:
1900:     case PEP_JD_PROJECTION_ORTHOGONAL:
1901:       if (pjd->proj != proj) {
1902:         pep->state = PEP_STATE_INITIAL;
1903:         pjd->proj = proj;
1904:       }
1905:       break;
1906:     default:
1907:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'proj' value");
1908:   }
1909:   return(0);
1910: }

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

1915:    Logically Collective on pep

1917:    Input Parameters:
1918: +  pep  - the eigenproblem solver context
1919: -  proj - the type of projection

1921:    Options Database Key:
1922: .  -pep_jd_projection - the projection type, either orthogonal or harmonic

1924:    Level: advanced

1926: .seealso: PEPJDGetProjection()
1927: @*/
1928: PetscErrorCode PEPJDSetProjection(PEP pep,PEPJDProjection proj)
1929: {

1935:   PetscTryMethod(pep,"PEPJDSetProjection_C",(PEP,PEPJDProjection),(pep,proj));
1936:   return(0);
1937: }

1939: PetscErrorCode PEPJDGetProjection_JD(PEP pep,PEPJDProjection *proj)
1940: {
1941:   PEP_JD *pjd = (PEP_JD*)pep->data;

1944:   *proj = pjd->proj;
1945:   return(0);
1946: }

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

1951:    Not Collective

1953:    Input Parameter:
1954: .  pep - the eigenproblem solver context

1956:    Output Parameter:
1957: .  proj - the type of projection

1959:    Level: advanced

1961: .seealso: PEPJDSetProjection()
1962: @*/
1963: PetscErrorCode PEPJDGetProjection(PEP pep,PEPJDProjection *proj)
1964: {

1970:   PetscUseMethod(pep,"PEPJDGetProjection_C",(PEP,PEPJDProjection*),(pep,proj));
1971:   return(0);
1972: }

1974: PetscErrorCode PEPSetFromOptions_JD(PetscOptionItems *PetscOptionsObject,PEP pep)
1975: {
1976:   PetscErrorCode  ierr;
1977:   PetscBool       flg,b1;
1978:   PetscReal       r1;
1979:   PetscInt        i1;
1980:   PEPJDProjection proj;

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

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

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

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

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

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

2000:   PetscOptionsTail();
2001:   return(0);
2002: }

2004: PetscErrorCode PEPView_JD(PEP pep,PetscViewer viewer)
2005: {
2007:   PEP_JD         *pjd = (PEP_JD*)pep->data;
2008:   PetscBool      isascii;

2011:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
2012:   if (isascii) {
2013:     PetscViewerASCIIPrintf(viewer,"  %d%% of basis vectors kept after restart\n",(int)(100*pjd->keep));
2014:     PetscViewerASCIIPrintf(viewer,"  threshold for changing the target in the correction equation (fix): %g\n",(double)pjd->fix);
2015:     PetscViewerASCIIPrintf(viewer,"  projection type: %s\n",PEPJDProjectionTypes[pjd->proj]);
2016:     PetscViewerASCIIPrintf(viewer,"  maximum allowed minimality index: %d\n",pjd->mmidx);
2017:     if (pjd->reusepc) { PetscViewerASCIIPrintf(viewer,"  reusing the preconditioner\n"); }
2018:   }
2019:   return(0);
2020: }

2022: PetscErrorCode PEPSetDefaultST_JD(PEP pep)
2023: {
2025:   KSP            ksp;

2028:   if (!((PetscObject)pep->st)->type_name) {
2029:     STSetType(pep->st,STPRECOND);
2030:     STPrecondSetKSPHasMat(pep->st,PETSC_TRUE);
2031:   }
2032:   STSetTransform(pep->st,PETSC_FALSE);
2033:   STGetKSP(pep->st,&ksp);
2034:   if (!((PetscObject)ksp)->type_name) {
2035:     KSPSetType(ksp,KSPBCGSL);
2036:     KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
2037:   }
2038:   return(0);
2039: }

2041: PetscErrorCode PEPReset_JD(PEP pep)
2042: {
2044:   PEP_JD         *pjd = (PEP_JD*)pep->data;
2045:   PetscInt       i;

2048:   for (i=0;i<pep->nmat;i++) {
2049:     BVDestroy(pjd->TV+i);
2050:   }
2051:   if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVDestroy(&pjd->W); }
2052:   if (pjd->ld>1) {
2053:     BVDestroy(&pjd->V);
2054:     for (i=0;i<pep->nmat;i++) {
2055:       BVDestroy(pjd->AX+i);
2056:     }
2057:     BVDestroy(&pjd->N[0]);
2058:     BVDestroy(&pjd->N[1]);
2059:     PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
2060:   }
2061:   PetscFree2(pjd->TV,pjd->AX);
2062:   return(0);
2063: }

2065: PetscErrorCode PEPDestroy_JD(PEP pep)
2066: {

2070:   PetscFree(pep->data);
2071:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
2072:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
2073:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",NULL);
2074:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",NULL);
2075:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",NULL);
2076:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",NULL);
2077:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",NULL);
2078:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",NULL);
2079:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",NULL);
2080:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",NULL);
2081:   return(0);
2082: }

2084: SLEPC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
2085: {
2086:   PEP_JD         *pjd;

2090:   PetscNewLog(pep,&pjd);
2091:   pep->data = (void*)pjd;

2093:   pjd->fix   = 0.01;
2094:   pjd->mmidx = 0;

2096:   pep->ops->solve          = PEPSolve_JD;
2097:   pep->ops->setup          = PEPSetUp_JD;
2098:   pep->ops->setfromoptions = PEPSetFromOptions_JD;
2099:   pep->ops->destroy        = PEPDestroy_JD;
2100:   pep->ops->reset          = PEPReset_JD;
2101:   pep->ops->view           = PEPView_JD;
2102:   pep->ops->setdefaultst   = PEPSetDefaultST_JD;

2104:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
2105:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
2106:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",PEPJDSetFix_JD);
2107:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",PEPJDGetFix_JD);
2108:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",PEPJDSetReusePreconditioner_JD);
2109:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",PEPJDGetReusePreconditioner_JD);
2110:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",PEPJDSetMinimalityIndex_JD);
2111:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",PEPJDGetMinimalityIndex_JD);
2112:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",PEPJDSetProjection_JD);
2113:   PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",PEPJDGetProjection_JD);
2114:   return(0);
2115: }