Actual source code: pjd.c
slepc-3.11.0 2019-03-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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.
18: Based on code contributed by the authors of [2] below.
20: References:
22: [1] G.L.G. Sleijpen et al., "Jacobi-Davidson type methods for
23: generalized eigenproblems and polynomial eigenproblems", BIT
24: 36(3):595-633, 1996.
26: [2] Feng-Nan Hwang, Zih-Hao Wei, Tsung-Ming Huang, Weichung Wang,
27: "A Parallel Additive Schwarz Preconditioned Jacobi-Davidson
28: Algorithm for Polynomial Eigenvalue Problems in Quantum Dot
29: Simulation", J. Comput. Phys. 229(8):2932-2947, 2010.
30: */
32: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
33: #include <slepcblaslapack.h>
35: typedef struct {
36: PetscReal keep; /* restart parameter */
37: PetscReal fix; /* fix parameter */
38: PetscBool reusepc; /* flag indicating whether pc is rebuilt or not */
39: BV V; /* work basis vectors to store the search space */
40: BV W; /* work basis vectors to store the test space */
41: 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) */
42: BV *AX; /* work basis vectors to store A_i*X for locked eigenvectors */
43: BV N[2]; /* auxiliary work BVs */
44: BV X; /* locked eigenvectors */
45: PetscScalar *T; /* matrix of the invariant pair */
46: PetscScalar *Tj; /* matrix containing the powers of the invariant pair matrix */
47: PetscScalar *XpX; /* X^H*X */
48: PetscInt ld; /* leading dimension for Tj and XpX */
49: PC pcshell; /* preconditioner including basic precond+projector */
50: Mat Pshell; /* auxiliary shell matrix */
51: PetscInt nlock; /* number of locked vectors in the invariant pair */
52: Vec vtempl; /* reference nested vector */
53: PetscInt midx; /* minimality index */
54: PetscInt mmidx; /* maximum allowed minimality index */
55: PEPJDProjection proj; /* projection type (orthogonal, harmonic) */
56: } PEP_JD;
58: typedef struct {
59: PEP pep;
60: PC pc; /* basic preconditioner */
61: Vec Bp[2]; /* preconditioned residual of derivative polynomial, B\p */
62: Vec u[2]; /* Ritz vector */
63: PetscScalar gamma[2]; /* precomputed scalar u'*B\p */
64: PetscScalar theta;
65: PetscScalar *M;
66: PetscScalar *ps;
67: PetscInt ld;
68: Vec *work;
69: Mat PPr;
70: BV X;
71: PetscInt n;
72: } PEP_JD_PCSHELL;
74: typedef struct {
75: Mat Pr,Pi; /* matrix polynomial evaluated at theta */
76: PEP pep;
77: Vec *work;
78: PetscScalar theta[2];
79: } PEP_JD_MATSHELL;
81: /*
82: Duplicate and resize auxiliary basis
83: */
84: static PetscErrorCode PEPJDDuplicateBasis(PEP pep,BV *basis)
85: {
86: PetscErrorCode ierr;
87: PEP_JD *pjd = (PEP_JD*)pep->data;
88: PetscInt nloc,m;
89: BVType type;
90: BVOrthogType otype;
91: BVOrthogRefineType oref;
92: PetscReal oeta;
93: BVOrthogBlockType oblock;
96: if (pjd->ld>1) {
97: BVCreate(PetscObjectComm((PetscObject)pep),basis);
98: BVGetSizes(pep->V,&nloc,NULL,&m);
99: nloc += pjd->ld-1;
100: BVSetSizes(*basis,nloc,PETSC_DECIDE,m);
101: BVGetType(pep->V,&type);
102: BVSetType(*basis,type);
103: BVGetOrthogonalization(pep->V,&otype,&oref,&oeta,&oblock);
104: BVSetOrthogonalization(*basis,otype,oref,oeta,oblock);
105: PetscObjectStateIncrease((PetscObject)*basis);
106: } else {
107: BVDuplicate(pep->V,basis);
108: }
109: return(0);
110: }
112: PetscErrorCode PEPSetUp_JD(PEP pep)
113: {
115: PEP_JD *pjd = (PEP_JD*)pep->data;
116: PetscBool isprecond,flg;
117: PetscInt i;
120: pep->lineariz = PETSC_FALSE;
121: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
122: if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
123: if (!pep->which) pep->which = PEP_TARGET_MAGNITUDE;
125: PetscObjectTypeCompare((PetscObject)pep->st,STPRECOND,&isprecond);
126: if (!isprecond) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"JD only works with PRECOND spectral transformation");
128: STGetTransform(pep->st,&flg);
129: if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");
131: if (!pjd->mmidx) pjd->mmidx = pep->nmat-1;
132: pjd->mmidx = PetscMin(pjd->mmidx,pep->nmat-1);
133: if (!pjd->keep) pjd->keep = 0.5;
134: PEPBasisCoefficients(pep,pep->pbc);
135: PEPAllocateSolution(pep,0);
136: PEPSetWorkVecs(pep,5);
137: pjd->ld = pep->nev;
138: #if !defined (PETSC_USE_COMPLEX)
139: pjd->ld++;
140: #endif
141: PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
142: for (i=0;i<pep->nmat;i++) {
143: PEPJDDuplicateBasis(pep,pjd->TV+i);
144: }
145: if (pjd->ld>1) {
146: PEPJDDuplicateBasis(pep,&pjd->V);
147: BVSetFromOptions(pjd->V);
148: for (i=0;i<pep->nmat;i++) {
149: BVDuplicateResize(pep->V,pjd->ld-1,pjd->AX+i);
150: }
151: BVDuplicateResize(pep->V,pjd->ld-1,pjd->N);
152: BVDuplicateResize(pep->V,pjd->ld-1,pjd->N+1);
153: pjd->X = pep->V;
154: PetscCalloc3((pjd->ld)*(pjd->ld),&pjd->XpX,pep->ncv*pep->ncv,&pjd->T,pjd->ld*pjd->ld*pep->nmat,&pjd->Tj);
155: } else pjd->V = pep->V;
156: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { PEPJDDuplicateBasis(pep,&pjd->W); }
157: else pjd->W = pjd->V;
158: DSSetType(pep->ds,DSPEP);
159: DSPEPSetDegree(pep->ds,pep->nmat-1);
160: if (pep->basis!=PEP_BASIS_MONOMIAL) {
161: DSPEPSetCoefficients(pep->ds,pep->pbc);
162: }
163: DSAllocate(pep->ds,pep->ncv);
164: return(0);
165: }
167: /*
168: Updates columns (low to (high-1)) of TV[i]
169: */
170: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
171: {
173: PEP_JD *pjd = (PEP_JD*)pep->data;
174: PetscInt pp,col,i,nloc,nconv;
175: Vec v1,v2,t1,t2;
176: PetscScalar *array1,*array2,*x2,*xx,*N,*Np,*y2=NULL,zero=0.0,sone=1.0,*pT,fact,*psc;
177: PetscReal *cg,*ca,*cb;
178: PetscMPIInt rk,np;
179: PetscBLASInt n_,ld_,one=1;
180: Mat T;
181: BV pbv;
184: ca = pep->pbc; cb = ca+pep->nmat; cg = cb + pep->nmat;
185: nconv = pjd->nlock;
186: PetscMalloc5(nconv,&x2,nconv,&xx,nconv*nconv,&pT,nconv*nconv,&N,nconv*nconv,&Np);
187: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
188: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
189: BVGetSizes(pep->V,&nloc,NULL,NULL);
190: t1 = w[0];
191: t2 = w[1];
192: PetscBLASIntCast(pjd->nlock,&n_);
193: PetscBLASIntCast(pjd->ld,&ld_);
194: if (nconv){
195: for (i=0;i<nconv;i++) {
196: PetscMemcpy(pT+i*nconv,pjd->T+i*pep->ncv,nconv*sizeof(PetscScalar));
197: }
198: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,pT,&T);
199: }
200: for (col=low;col<high;col++) {
201: BVGetColumn(pjd->V,col,&v1);
202: VecGetArray(v1,&array1);
203: if (nconv>0) {
204: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
205: }
206: VecPlaceArray(t1,array1);
207: if (nconv) {
208: BVSetActiveColumns(pjd->N[0],0,nconv);
209: BVSetActiveColumns(pjd->N[1],0,nconv);
210: BVDotVec(pjd->X,t1,xx);
211: }
212: for (pp=pep->nmat-1;pp>=0;pp--) {
213: BVGetColumn(pjd->TV[pp],col,&v2);
214: VecGetArray(v2,&array2);
215: VecPlaceArray(t2,array2);
216: MatMult(pep->A[pp],t1,t2);
217: if (nconv) {
218: if (pp<pep->nmat-3) {
219: BVMult(pjd->N[0],1.0,-cg[pp+2],pjd->AX[pp+1],NULL);
220: MatShift(T,-cb[pp+1]);
221: BVMult(pjd->N[0],1.0/ca[pp],1.0/ca[pp],pjd->N[1],T);
222: pbv = pjd->N[0]; pjd->N[0] = pjd->N[1]; pjd->N[1] = pbv;
223: BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
224: MatShift(T,cb[pp+1]);
225: } else if (pp==pep->nmat-3) {
226: BVCopy(pjd->AX[pp+2],pjd->N[0]);
227: BVScale(pjd->N[0],1/ca[pp+1]);
228: BVCopy(pjd->AX[pp+1],pjd->N[1]);
229: MatShift(T,-cb[pp+1]);
230: BVMult(pjd->N[1],1.0/ca[pp],1.0/ca[pp],pjd->N[0],T);
231: BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
232: MatShift(T,cb[pp+1]);
233: } else if (pp==pep->nmat-2) {
234: BVMultVec(pjd->AX[pp+1],1.0/ca[pp],1.0,t2,x2);
235: }
236: if (pp<pjd->midx) {
237: y2 = array2+nloc;
238: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&n_,&sone,pjd->Tj+pjd->ld*pjd->ld*pp,&ld_,xx,&one,&zero,y2,&one));
239: if (pp<pjd->midx-2) {
240: fact = -cg[pp+2];
241: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&fact,Np,&n_));
242: fact = 1/ca[pp];
243: MatShift(T,-cb[pp+1]);
244: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&fact,N,&n_,pT,&n_,&fact,Np,&n_));
245: MatShift(T,cb[pp+1]);
246: psc = Np; Np = N; N = psc;
247: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
248: } else if (pp==pjd->midx-2) {
249: fact = 1/ca[pp];
250: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&fact,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&zero,N,&n_));
251: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
252: } else if (pp==pjd->midx-1) {
253: PetscMemzero(Np,nconv*nconv*sizeof(PetscScalar));
254: }
255: }
256: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
257: }
258: VecResetArray(t2);
259: VecRestoreArray(v2,&array2);
260: BVRestoreColumn(pjd->TV[pp],col,&v2);
261: }
262: VecResetArray(t1);
263: VecRestoreArray(v1,&array1);
264: BVRestoreColumn(pjd->V,col,&v1);
265: }
266: if (nconv) {MatDestroy(&T);}
267: PetscFree5(x2,xx,pT,N,Np);
268: return(0);
269: }
271: /*
272: RRQR of X. Xin*P=Xou*R. Rank of R is rk
273: */
274: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
275: {
276: #if defined(SLEPC_MISSING_LAPACK_GEQP3) || defined(PETSC_MISSING_LAPACK_ORGQR)
278: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQP3/QRGQR - Lapack routines are unavailable");
279: #else
281: PetscInt i,j,n,r;
282: PetscBLASInt row_,col_,ldx_,*p,lwork,info,n_;
283: PetscScalar *tau,*work;
284: PetscReal tol,*rwork;
287: PetscBLASIntCast(row,&row_);
288: PetscBLASIntCast(col,&col_);
289: PetscBLASIntCast(ldx,&ldx_);
290: n = PetscMin(row,col);
291: PetscBLASIntCast(n,&n_);
292: lwork = 3*col_+1;
293: PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
294: for (i=1;i<col;i++) p[i] = 0;
295: p[0] = 1;
297: /* rank revealing QR */
298: #if defined(PETSC_USE_COMPLEX)
299: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
300: #else
301: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
302: #endif
303: SlepcCheckLapackInfo("geqp3",info);
304: if (P) for (i=0;i<col;i++) P[i] = p[i]-1;
306: /* rank computation */
307: tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
308: r = 1;
309: for (i=1;i<n;i++) {
310: if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
311: else break;
312: }
313: if (rk) *rk=r;
315: /* copy upper triangular matrix if requested */
316: if (R) {
317: for (i=0;i<r;i++) {
318: PetscMemzero(R+i*ldr,r*sizeof(PetscScalar));
319: for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
320: }
321: }
322: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
323: SlepcCheckLapackInfo("orgqr",info);
324: PetscFree4(p,tau,work,rwork);
325: return(0);
326: #endif
327: }
329: /*
330: Application of extended preconditioner
331: */
332: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
333: {
334: PetscInt i,j,nloc,n,ld=0;
335: PetscMPIInt np;
336: Vec tx,ty;
337: PEP_JD_PCSHELL *ctx;
338: PetscErrorCode ierr;
339: const PetscScalar *array1;
340: PetscScalar *x2=NULL,*t=NULL,*ps=NULL,*array2,zero=0.0,sone=1.0;
341: PetscBLASInt one=1.0,ld_,n_,ncv_;
342: PEP_JD *pjd=NULL;
345: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
346: PCShellGetContext(pc,(void**)&ctx);
347: n = ctx->n;
348: if (n) {
349: pjd = (PEP_JD*)ctx->pep->data;
350: ps = ctx->ps;
351: ld = pjd->ld;
352: PetscMalloc2(n,&x2,n,&t);
353: VecGetLocalSize(ctx->work[0],&nloc);
354: VecGetArrayRead(x,&array1);
355: for (i=0;i<n;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
356: VecRestoreArrayRead(x,&array1);
357: }
359: /* y = B\x apply PC */
360: tx = ctx->work[0];
361: ty = ctx->work[1];
362: VecGetArrayRead(x,&array1);
363: VecPlaceArray(tx,array1);
364: VecGetArray(y,&array2);
365: VecPlaceArray(ty,array2);
366: PCApply(ctx->pc,tx,ty);
367: if (n) {
368: PetscBLASIntCast(ld,&ld_);
369: PetscBLASIntCast(n,&n_);
370: for (i=0;i<n;i++) {
371: t[i] = 0.0;
372: for (j=0;j<n;j++) t[i] += ctx->M[i+j*ld]*x2[j];
373: }
374: if (pjd->midx==1) {
375: PetscBLASIntCast(ctx->pep->ncv,&ncv_);
376: for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] -= ctx->theta;
377: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,pjd->T,&ncv_,t,&one,&zero,x2,&one));
378: for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] += ctx->theta;
379: for (i=0;i<n;i++) array2[nloc+i] = x2[i];
380: for (i=0;i<n;i++) x2[i] = -t[i];
381: } else {
382: for (i=0;i<n;i++) array2[nloc+i] = t[i];
383: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,ps,&ld_,t,&one,&zero,x2,&one));
384: }
385: for (i=0;i<n;i++) array2[nloc+i] /= PetscSqrtReal(np);
386: BVSetActiveColumns(pjd->X,0,n);
387: BVMultVec(pjd->X,-1.0,1.0,ty,x2);
388: PetscFree2(x2,t);
389: }
390: VecResetArray(tx);
391: VecResetArray(ty);
392: VecRestoreArrayRead(x,&array1);
393: VecRestoreArray(y,&array2);
394: return(0);
395: }
397: /*
398: Application of shell preconditioner:
399: y = B\x - eta*B\p, with eta = (u'*B\x)/(u'*B\p)
400: */
401: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
402: {
404: PetscScalar rr,rx,xr,xx,eta;
405: PEP_JD_PCSHELL *ctx;
406: PetscInt sz;
407: const Vec *xs,*ys;
410: PCShellGetContext(pc,(void**)&ctx);
411: VecCompGetSubVecs(x,&sz,&xs);
412: VecCompGetSubVecs(y,NULL,&ys);
413: /* y = B\x apply extended PC */
414: PEPJDExtendedPCApply(pc,xs[0],ys[0]);
415: if (sz==2) {
416: PEPJDExtendedPCApply(pc,xs[1],ys[1]);
417: }
419: /* Compute eta = u'*y / u'*Bp */
420: VecDot(ys[0],ctx->u[0],&rr);
421: eta = -rr*ctx->gamma[0];
423: if (sz==2) {
424: VecDot(ys[0],ctx->u[1],&xr);
425: VecDot(ys[1],ctx->u[0],&rx);
426: VecDot(ys[1],ctx->u[1],&xx);
427: eta += -ctx->gamma[0]*xx-ctx->gamma[1]*(-xr+rx);
428: }
429: eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];
431: /* y = y - eta*Bp */
432: VecAXPY(ys[0],eta,ctx->Bp[0]);
433: if (sz==2) {
434: VecAXPY(ys[1],eta,ctx->Bp[1]);
435: eta = -ctx->gamma[1]*(rr+xx)+ctx->gamma[0]*(-xr+rx);
436: eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];
437: VecAXPY(ys[0],eta,ctx->Bp[1]);
438: VecAXPY(ys[1],-eta,ctx->Bp[0]);
439: }
440: return(0);
441: }
443: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
444: {
446: PetscMPIInt np,rk,count;
447: PetscScalar *array1,*array2;
448: PetscInt nloc;
451: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
452: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
453: BVGetSizes(pep->V,&nloc,NULL,NULL);
454: if (v) {
455: VecGetArray(v,&array1);
456: VecGetArray(vex,&array2);
457: if (back) {
458: PetscMemcpy(array1,array2,nloc*sizeof(PetscScalar));
459: } else {
460: PetscMemcpy(array2,array1,nloc*sizeof(PetscScalar));
461: }
462: VecRestoreArray(v,&array1);
463: VecRestoreArray(vex,&array2);
464: }
465: if (a) {
466: VecGetArray(vex,&array2);
467: if (back) {
468: PetscMemcpy(a,array2+nloc+off,na*sizeof(PetscScalar));
469: PetscMPIIntCast(na,&count);
470: MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
471: } else {
472: PetscMemcpy(array2+nloc+off,a,na*sizeof(PetscScalar));
473: PetscMPIIntCast(na,&count);
474: MPI_Bcast(array2+nloc+off,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
475: }
476: VecRestoreArray(vex,&array2);
477: }
478: return(0);
479: }
481: /* Computes Phi^hat(lambda) times a vector or its derivative (depends on beval)
482: if no vector is provided returns a matrix
483: */
484: static PetscErrorCode PEPJDEvaluateHatBasis(PEP pep,PetscInt n,PetscScalar *H,PetscInt ldh,PetscScalar *beval,PetscScalar *t,PetscInt idx,PetscScalar *qpp,PetscScalar *qp,PetscScalar *q)
485: {
487: PetscInt j,i;
488: PetscBLASInt n_,ldh_,one=1;
489: PetscReal *a,*b,*g;
490: PetscScalar sone=1.0,zero=0.0;
493: a = pep->pbc; b=a+pep->nmat; g=b+pep->nmat;
494: PetscBLASIntCast(n,&n_);
495: PetscBLASIntCast(ldh,&ldh_);
496: if (idx<1) {
497: PetscMemzero(q,(t?n:n*n)*sizeof(PetscScalar));
498: } else if (idx==1) {
499: if (t) {for (j=0;j<n;j++) q[j] = t[j]*beval[idx-1]/a[0];}
500: else {
501: PetscMemzero(q,n*n*sizeof(PetscScalar));
502: for (j=0;j<n;j++) q[(j+1)*n] = beval[idx-1]/a[0];
503: }
504: } else {
505: if (t) {
506: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,H,&ldh_,qp,&one,&zero,q,&one));
507: for (j=0;j<n;j++) {
508: q[j] += beval[idx-1]*t[j]-b[idx-1]*qp[j]-g[idx-1]*qpp[j];
509: q[j] /= a[idx-1];
510: }
511: } else {
512: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,H,&ldh_,qp,&n_,&zero,q,&n_));
513: for (j=0;j<n;j++) {
514: q[j+n*j] += beval[idx-1];
515: for (i=0;i<n;i++) {
516: q[i+n*j] += -b[idx-1]*qp[j*n+i]-g[idx-1]*qpp[j*n+i];
517: q[i+n*j] /= a[idx-1];
518: }
519: }
520: }
521: }
522: return(0);
523: }
525: static PetscErrorCode PEPJDComputeResidual(PEP pep,PetscBool derivative,PetscInt sz,Vec *u,PetscScalar *theta,Vec *p,Vec *work)
526: {
527: PEP_JD *pjd = (PEP_JD*)pep->data;
529: PetscMPIInt rk,np,count;
530: Vec tu,tui=NULL,tp,tpi=NULL,w;
531: PetscScalar *dval,*dvali,*array1,*array2,*arrayi1,*arrayi2;
532: PetscScalar *x2=NULL,*x2i=NULL,*y2,*y2i,*qj=NULL,*qji=NULL,*qq,*tt=NULL,*xx=NULL,*xxi=NULL,sone=1.0;
533: PetscInt i,j,nconv,nloc;
534: PetscBLASInt n,ld,one=1;
537: nconv = pjd->nlock;
538: if (!nconv) {
539: PetscMalloc1(2*sz*pep->nmat,&dval);
540: } else {
541: PetscMalloc5(2*pep->nmat,&dval,2*nconv,&xx,nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*pep->nmat,&qj);
542: if (sz==2) x2i = x2+nconv;
543: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
544: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
545: BVGetSizes(pep->V,&nloc,NULL,NULL);
546: VecGetArray(u[0],&array1);
547: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
548: VecRestoreArray(u[0],&array1);
549: if (sz==2) {
550: VecGetArray(u[1],&arrayi1);
551: for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]* PetscSqrtReal(np);
552: VecRestoreArray(u[1],&arrayi1);
553: }
554: }
555: dvali = dval+pep->nmat;
556: tu = work[0];
557: tp = work[1];
558: w = work[2];
559: VecGetArray(u[0],&array1);
560: VecPlaceArray(tu,array1);
561: VecGetArray(p[0],&array2);
562: VecPlaceArray(tp,array2);
563: VecSet(tp,0.0);
564: if (sz==2) {
565: tui = work[3];
566: tpi = work[4];
567: VecGetArray(u[1],&arrayi1);
568: VecPlaceArray(tui,arrayi1);
569: VecGetArray(p[1],&arrayi2);
570: VecPlaceArray(tpi,arrayi2);
571: VecSet(tpi,0.0);
572: }
573: if (derivative) {
574: PEPEvaluateBasisDerivative(pep,theta[0],theta[1],dval,dvali);
575: } else {
576: PEPEvaluateBasis(pep,theta[0],theta[1],dval,dvali);
577: }
578: for (i=derivative?1:0;i<pep->nmat;i++) {
579: MatMult(pep->A[i],tu,w);
580: VecAXPY(tp,dval[i],w);
581: if (sz==2) {
582: VecAXPY(tpi,dvali[i],w);
583: MatMult(pep->A[i],tui,w);
584: VecAXPY(tpi,dval[i],w);
585: VecAXPY(tp,-dvali[i],w);
586: }
587: }
588: if (nconv) {
589: for (i=0;i<pep->nmat;i++) {
590: 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);
591: }
592: if (sz==2) {
593: qji = qj+nconv*pep->nmat;
594: qq = qji+nconv*pep->nmat;
595: for (i=0;i<pep->nmat;i++) {
596: 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);
597: }
598: for (i=0;i<nconv*pep->nmat;i++) qj[i] -= qji[i];
599: for (i=0;i<pep->nmat;i++) {
600: 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);
601: 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);
602: }
603: for (i=0;i<nconv*pep->nmat;i++) qji[i] += qq[i];
604: for (i=derivative?2:1;i<pep->nmat;i++) {
605: BVMultVec(pjd->AX[i],1.0,1.0,tpi,qji+i*nconv);
606: }
607: }
608: for (i=derivative?2:1;i<pep->nmat;i++) {
609: BVMultVec(pjd->AX[i],1.0,1.0,tp,qj+i*nconv);
610: }
612: /* extended vector part */
613: BVSetActiveColumns(pjd->X,0,nconv);
614: BVDotVec(pjd->X,tu,xx);
615: xxi = xx+nconv;
616: if (sz==2) {
617: BVDotVec(pjd->X,tui,xxi);
618: } else {
619: PetscMemzero(xxi,nconv*sizeof(PetscScalar));
620: }
621: if (rk==np-1) {
622: PetscBLASIntCast(nconv,&n);
623: PetscBLASIntCast(pjd->ld,&ld);
624: y2 = array2+nloc;
625: PetscMemzero(y2,nconv*sizeof(PetscScalar));
626: for (j=derivative?1:0;j<pjd->midx;j++) {
627: for (i=0;i<nconv;i++) tt[i] = dval[j]*xx[i]-dvali[j]*xxi[i];
628: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
629: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
630: }
631: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
632: if (sz==2) {
633: y2i = arrayi2+nloc;
634: PetscMemzero(y2i,nconv*sizeof(PetscScalar));
635: for (j=derivative?1:0;j<pjd->midx;j++) {
636: for (i=0;i<nconv;i++) tt[i] = dval[j]*xxi[i]+dvali[j]*xx[i];
637: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
638: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
639: }
640: for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
641: }
642: }
643: PetscMPIIntCast(nconv,&count);
644: MPI_Bcast(array2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
645: if (sz==2) {MPI_Bcast(arrayi2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));}
646: }
647: if (nconv) {
648: PetscFree5(dval,xx,tt,x2,qj);
649: } else {
650: PetscFree(dval);
651: }
652: VecResetArray(tu);
653: VecRestoreArray(u[0],&array1);
654: VecResetArray(tp);
655: VecRestoreArray(p[0],&array2);
656: if (sz==2) {
657: VecResetArray(tui);
658: VecRestoreArray(u[1],&arrayi1);
659: VecResetArray(tpi);
660: VecRestoreArray(p[1],&arrayi2);
661: }
662: return(0);
663: }
665: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
666: {
667: PEP_JD *pjd = (PEP_JD*)pep->data;
669: PetscScalar *tt,target[2];
670: Vec vg,wg;
671: PetscInt i;
672: PetscReal norm;
675: PetscMalloc1(pjd->ld-1,&tt);
676: if (pep->nini==0) {
677: BVSetRandomColumn(pjd->V,0);
678: for (i=0;i<pjd->ld-1;i++) tt[i] = 0.0;
679: BVGetColumn(pjd->V,0,&vg);
680: PEPJDCopyToExtendedVec(pep,NULL,tt,pjd->ld-1,0,vg,PETSC_FALSE);
681: BVRestoreColumn(pjd->V,0,&vg);
682: BVNormColumn(pjd->V,0,NORM_2,&norm);
683: BVScaleColumn(pjd->V,0,1.0/norm);
684: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
685: BVGetColumn(pjd->V,0,&vg);
686: BVGetColumn(pjd->W,0,&wg);
687: VecSet(wg,0.0);
688: target[0] = pep->target; target[1] = 0.0;
689: PEPJDComputeResidual(pep,PETSC_TRUE,1,&vg,target,&wg,w);
690: BVRestoreColumn(pjd->W,0,&wg);
691: BVRestoreColumn(pjd->V,0,&vg);
692: BVNormColumn(pjd->W,0,NORM_2,&norm);
693: BVScaleColumn(pjd->W,0,1.0/norm);
694: }
695: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for initial vectors not implemented yet");
696: PetscFree(tt);
697: return(0);
698: }
700: static PetscErrorCode PEPJDShellMatMult(Mat P,Vec x,Vec y)
701: {
702: PetscErrorCode ierr;
703: PEP_JD_MATSHELL *matctx;
704: PEP_JD *pjd;
705: PetscInt i,j,nconv,nloc,nmat,ldt,ncv,sz;
706: Vec tx,ty,txi=NULL,tyi=NULL;
707: const Vec *xs,*ys;
708: PetscScalar *array2,*arrayi1,*arrayi2,*array1;
709: PetscScalar *x2=NULL,*x2i=NULL,*y2,*y2i,*tt=NULL,*xx=NULL,*xxi,theta[2],sone=1.0,*qj,*qji=NULL,*qq,*val,*vali=NULL;
710: PetscBLASInt n,ld,one=1;
711: PetscMPIInt np;
714: MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
715: MatShellGetContext(P,(void**)&matctx);
716: pjd = (PEP_JD*)(matctx->pep->data);
717: nconv = pjd->nlock;
718: nmat = matctx->pep->nmat;
719: ncv = matctx->pep->ncv;
720: ldt = pjd->ld;
721: VecCompGetSubVecs(x,&sz,&xs);
722: VecCompGetSubVecs(y,NULL,&ys);
723: theta[0] = matctx->theta[0];
724: theta[1] = (sz==2)?matctx->theta[1]:0.0;
725: if (nconv>0) {
726: PetscMalloc5(nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*nmat,&qj,2*nconv,&xx,2*nmat,&val);
727: if (sz==2) x2i = x2+nconv;
728: BVGetSizes(matctx->pep->V,&nloc,NULL,NULL);
729: VecGetArray(xs[0],&array1);
730: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
731: VecRestoreArray(xs[0],&array1);
732: if (sz==2) {
733: VecGetArray(xs[1],&arrayi1);
734: for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]* PetscSqrtReal(np);
735: VecRestoreArray(xs[1],&arrayi1);
736: }
737: vali = val+nmat;
738: }
739: tx = matctx->work[0];
740: ty = matctx->work[1];
741: VecGetArray(xs[0],&array1);
742: VecPlaceArray(tx,array1);
743: VecGetArray(ys[0],&array2);
744: VecPlaceArray(ty,array2);
745: MatMult(matctx->Pr,tx,ty);
746: if (sz==2) {
747: txi = matctx->work[2];
748: tyi = matctx->work[3];
749: VecGetArray(xs[1],&arrayi1);
750: VecPlaceArray(txi,arrayi1);
751: VecGetArray(ys[1],&arrayi2);
752: VecPlaceArray(tyi,arrayi2);
753: MatMult(matctx->Pr,txi,tyi);
754: if (theta[1]!=0.0) {
755: MatMult(matctx->Pi,txi,matctx->work[4]);
756: VecAXPY(ty,-1.0,matctx->work[4]);
757: MatMult(matctx->Pi,tx,matctx->work[4]);
758: VecAXPY(tyi,1.0,matctx->work[4]);
759: }
760: }
761: if (nconv) {
762: PEPEvaluateBasis(matctx->pep,theta[0],theta[1],val,vali);
763: for (i=0;i<nmat;i++) {
764: 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);
765: }
766: if (sz==2) {
767: qji = qj+nconv*nmat;
768: qq = qji+nconv*nmat;
769: for (i=0;i<nmat;i++) {
770: 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);
771: }
772: for (i=0;i<nconv*nmat;i++) qj[i] -= qji[i];
773: for (i=0;i<nmat;i++) {
774: 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);
775: 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);
776: }
777: for (i=0;i<nconv*nmat;i++) qji[i] += qq[i];
778: for (i=1;i<matctx->pep->nmat;i++) {
779: BVMultVec(pjd->AX[i],1.0,1.0,tyi,qji+i*nconv);
780: }
781: }
782: for (i=1;i<nmat;i++) {
783: BVMultVec(pjd->AX[i],1.0,1.0,ty,qj+i*nconv);
784: }
786: /* extended vector part */
787: BVSetActiveColumns(pjd->X,0,nconv);
788: BVDotVec(pjd->X,tx,xx);
789: xxi = xx+nconv;
790: if (sz==2) {
791: BVDotVec(pjd->X,txi,xxi);
792: } else {
793: PetscMemzero(xxi,nconv*sizeof(PetscScalar));
794: }
795: PetscBLASIntCast(pjd->nlock,&n);
796: PetscBLASIntCast(ldt,&ld);
797: y2 = array2+nloc;
798: PetscMemzero(y2,nconv*sizeof(PetscScalar));
799: for (j=0;j<pjd->midx;j++) {
800: for (i=0;i<nconv;i++) tt[i] = val[j]*xx[i]-vali[j]*xxi[i];
801: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
802: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
803: }
804: if (sz==2) {
805: y2i = arrayi2+nloc;
806: PetscMemzero(y2i,nconv*sizeof(PetscScalar));
807: for (j=0;j<pjd->midx;j++) {
808: for (i=0;i<nconv;i++) tt[i] = val[j]*xxi[i]+vali[j]*xx[i];
809: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
810: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
811: }
812: }
813: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
814: PetscFree5(tt,x2,qj,xx,val);
815: if (sz==2) {
816: for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
817: }
818: }
819: VecResetArray(tx);
820: VecRestoreArray(xs[0],&array1);
821: VecResetArray(ty);
822: VecRestoreArray(ys[0],&array2);
823: if (sz==2) {
824: VecResetArray(txi);
825: VecRestoreArray(xs[1],&arrayi1);
826: VecResetArray(tyi);
827: VecRestoreArray(ys[1],&arrayi2);
828: }
829: return(0);
830: }
832: static PetscErrorCode PEPJDSellMatCreateVecs(Mat A,Vec *right,Vec *left)
833: {
834: PetscErrorCode ierr;
835: PEP_JD_MATSHELL *matctx;
836: PEP_JD *pjd;
837: PetscInt kspsf=1,i;
838: Vec v[2];
841: MatShellGetContext(A,(void**)&matctx);
842: pjd = (PEP_JD*)(matctx->pep->data);
843: #if !defined (PETSC_USE_COMPLEX)
844: kspsf = 2;
845: #endif
846: for (i=0;i<kspsf;i++){
847: BVCreateVec(pjd->V,v+i);
848: }
849: if (right) {
850: VecCreateCompWithVecs(v,kspsf,pjd->vtempl,right);
851: }
852: if (left) {
853: VecCreateCompWithVecs(v,kspsf,pjd->vtempl,left);
854: }
855: for (i=0;i<kspsf;i++) {
856: VecDestroy(&v[i]);
857: }
858: return(0);
859: }
861: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
862: {
863: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GETRI) || defined(PETSC_MISSING_LAPACK_GETRF)
865: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GETRI/GETRF - Lapack routines are unavailable");
866: #else
868: PEP_JD *pjd = (PEP_JD*)pep->data;
869: PEP_JD_PCSHELL *pcctx;
870: PetscInt i,j,k,n=pjd->nlock,ld=pjd->ld,deg=pep->nmat-1;
871: PetscScalar *M,*ps,*work,*U,*V,*S,*Sp,*Spp,snone=-1.0,sone=1.0,zero=0.0,*val;
872: PetscReal tol,maxeig=0.0,*sg,*rwork;
873: PetscBLASInt n_,info,ld_,*p,lw_,rk=0;
876: if (n) {
877: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
878: pcctx->theta = theta;
879: pcctx->n = n;
880: M = pcctx->M;
881: PetscBLASIntCast(n,&n_);
882: PetscBLASIntCast(ld,&ld_);
883: if (pjd->midx==1) {
884: PetscMemcpy(M,pjd->XpX,ld*ld*sizeof(PetscScalar));
885: PetscCalloc2(10*n,&work,n,&p);
886: } else {
887: ps = pcctx->ps;
888: PetscCalloc7(2*n*n,&U,3*n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p,deg+1,&val);
889: V = U+n*n;
890: /* pseudo-inverse */
891: for (j=0;j<n;j++) {
892: for (i=0;i<n;i++) S[n*j+i] = -pjd->T[pep->ncv*j+i];
893: S[n*j+j] += theta;
894: }
895: lw_ = 10*n_;
896: #if !defined (PETSC_USE_COMPLEX)
897: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
898: #else
899: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
900: #endif
901: SlepcCheckLapackInfo("gesvd",info);
902: for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
903: tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
904: for (j=0;j<n;j++) {
905: if (sg[j]>tol) {
906: for (i=0;i<n;i++) U[j*n+i] /= sg[j];
907: rk++;
908: } else break;
909: }
910: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));
912: /* compute M */
913: PEPEvaluateBasis(pep,theta,0.0,val,NULL);
914: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&snone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
915: PetscMemzero(S,2*n*n*sizeof(PetscScalar));
916: Sp = S+n*n;
917: for (j=0;j<n;j++) S[j*(n+1)] = 1.0;
918: for (k=1;k<pjd->midx;k++) {
919: 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];
920: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
921: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
922: Spp = Sp; Sp = S;
923: PEPJDEvaluateHatBasis(pep,n,pjd->T,pep->ncv,val,NULL,k+1,Spp,Sp,S);
924: }
925: }
926: /* inverse */
927: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
928: SlepcCheckLapackInfo("getrf",info);
929: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
930: SlepcCheckLapackInfo("getri",info);
931: if (pjd->midx==1) {
932: PetscFree2(work,p);
933: } else {
934: PetscFree7(U,S,sg,work,rwork,p,val);
935: }
936: }
937: return(0);
938: #endif
939: }
941: static PetscErrorCode PEPJDMatSetUp(PEP pep,PetscInt sz,PetscScalar *theta)
942: {
943: PetscErrorCode ierr;
944: PEP_JD *pjd = (PEP_JD*)pep->data;
945: PEP_JD_MATSHELL *matctx;
946: PEP_JD_PCSHELL *pcctx;
947: MatStructure str;
948: PetscScalar *vals,*valsi;
949: PetscBool skipmat=PETSC_FALSE;
950: PetscInt i;
951: Mat Pr=NULL;
954: if (sz==2 && theta[1]==0.0) sz = 1;
955: MatShellGetContext(pjd->Pshell,(void**)&matctx);
956: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
957: if (matctx->Pr && matctx->theta[0]==theta[0] && ((!matctx->Pi && sz==1) || (sz==2 && matctx->theta[1]==theta[1]))) {
958: if (pcctx->n == pjd->nlock) return(0);
959: skipmat = PETSC_TRUE;
960: }
961: if (!skipmat) {
962: PetscMalloc2(pep->nmat,&vals,pep->nmat,&valsi);
963: STGetMatStructure(pep->st,&str);
964: PEPEvaluateBasis(pep,theta[0],theta[1],vals,valsi);
965: if (!matctx->Pr) {
966: MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pr);
967: } else {
968: MatCopy(pep->A[0],matctx->Pr,str);
969: }
970: for (i=1;i<pep->nmat;i++) {
971: MatAXPY(matctx->Pr,vals[i],pep->A[i],str);
972: }
973: if (!pjd->reusepc) {
974: if (pcctx->PPr && sz==2) {
975: MatCopy(matctx->Pr,pcctx->PPr,str);
976: Pr = pcctx->PPr;
977: } else Pr = matctx->Pr;
978: }
979: matctx->theta[0] = theta[0];
980: if (sz==2) {
981: if (!matctx->Pi ) {
982: MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pi);
983: } else {
984: MatCopy(pep->A[1],matctx->Pi,str);
985: }
986: MatScale(matctx->Pi,valsi[1]);
987: for (i=2;i<pep->nmat;i++) {
988: MatAXPY(matctx->Pi,valsi[i],pep->A[i],str);
989: }
990: matctx->theta[1] = theta[1];
991: }
992: PetscFree2(vals,valsi);
993: }
994: if (!pjd->reusepc) {
995: if (!skipmat) {
996: PCSetOperators(pcctx->pc,Pr,Pr);
997: PCSetUp(pcctx->pc);
998: }
999: PEPJDUpdateExtendedPC(pep,theta[0]);
1000: }
1001: return(0);
1002: }
1004: static PetscErrorCode PEPJDCreateShellPC(PEP pep,Vec *ww)
1005: {
1006: PEP_JD *pjd = (PEP_JD*)pep->data;
1007: PEP_JD_PCSHELL *pcctx;
1008: PEP_JD_MATSHELL *matctx;
1009: KSP ksp;
1010: PetscInt nloc,mloc,kspsf=1;
1011: Vec v[2];
1012: PetscScalar target[2];
1013: PetscErrorCode ierr;
1014: Mat Pr;
1017: /* Create the reference vector */
1018: BVGetColumn(pjd->V,0,&v[0]);
1019: v[1] = v[0];
1020: #if !defined (PETSC_USE_COMPLEX)
1021: kspsf = 2;
1022: #endif
1023: VecCreateCompWithVecs(v,kspsf,NULL,&pjd->vtempl);
1024: BVRestoreColumn(pjd->V,0,&v[0]);
1025: PetscLogObjectParent((PetscObject)pep,(PetscObject)pjd->vtempl);
1027: /* Replace preconditioner with one containing projectors */
1028: PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
1029: PCSetType(pjd->pcshell,PCSHELL);
1030: PCShellSetName(pjd->pcshell,"PCPEPJD");
1031: PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
1032: PetscNew(&pcctx);
1033: PCShellSetContext(pjd->pcshell,pcctx);
1034: STGetKSP(pep->st,&ksp);
1035: BVCreateVec(pjd->V,&pcctx->Bp[0]);
1036: VecDuplicate(pcctx->Bp[0],&pcctx->Bp[1]);
1037: KSPGetPC(ksp,&pcctx->pc);
1038: PetscObjectReference((PetscObject)pcctx->pc);
1039: MatGetLocalSize(pep->A[0],&mloc,&nloc);
1040: if (pjd->ld>1) {
1041: nloc += pjd->ld-1; mloc += pjd->ld-1;
1042: }
1043: PetscNew(&matctx);
1044: MatCreateShell(PetscObjectComm((PetscObject)pep),kspsf*nloc,kspsf*mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
1045: MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)(void))PEPJDShellMatMult);
1046: MatShellSetOperation(pjd->Pshell,MATOP_CREATE_VECS,(void(*)(void))PEPJDSellMatCreateVecs);
1047: matctx->pep = pep;
1048: target[0] = pep->target; target[1] = 0.0;
1049: PEPJDMatSetUp(pep,1,target);
1050: Pr = matctx->Pr;
1051: pcctx->PPr = NULL;
1052: #if !defined(PETSC_USE_COMPLEX)
1053: if (!pjd->reusepc) {
1054: MatDuplicate(matctx->Pr,MAT_COPY_VALUES,&pcctx->PPr);
1055: Pr = pcctx->PPr;
1056: }
1057: #endif
1058: PCSetOperators(pcctx->pc,Pr,Pr);
1059: PCSetErrorIfFailure(pcctx->pc,PETSC_TRUE);
1060: KSPSetPC(ksp,pjd->pcshell);
1061: if (pjd->reusepc) {
1062: PCSetReusePreconditioner(pcctx->pc,PETSC_TRUE);
1063: KSPSetReusePreconditioner(ksp,PETSC_TRUE);
1064: }
1065: KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
1066: KSPSetUp(ksp);
1067: if (pjd->ld>1) {
1068: PetscMalloc2(pjd->ld*pjd->ld,&pcctx->M,pjd->ld*pjd->ld,&pcctx->ps);
1069: pcctx->pep = pep;
1070: }
1071: matctx->work = ww;
1072: pcctx->work = ww;
1073: return(0);
1074: }
1076: static PetscErrorCode PEPJDEigenvectors(PEP pep)
1077: {
1078: #if defined(SLEPC_MISSING_LAPACK_TREVC)
1080: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
1081: #else
1083: PEP_JD *pjd = (PEP_JD*)pep->data;
1084: PetscBLASInt ld,nconv,info,nc;
1085: PetscScalar *Z,*w;
1086: PetscReal *wr,norm;
1087: PetscInt i;
1088: Mat U;
1089: #if !defined(PETSC_USE_COMPLEX)
1090: Vec v,v1;
1091: #endif
1094: PetscMalloc3(pep->nconv*pep->nconv,&Z,3*pep->ncv,&wr,2*pep->ncv,&w);
1095: PetscBLASIntCast(pep->ncv,&ld);
1096: PetscBLASIntCast(pep->nconv,&nconv);
1097: #if !defined(PETSC_USE_COMPLEX)
1098: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
1099: #else
1100: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
1101: #endif
1102: SlepcCheckLapackInfo("trevc",info);
1103: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
1104: BVSetActiveColumns(pjd->X,0,pep->nconv);
1105: BVMultInPlace(pjd->X,U,0,pep->nconv);
1106: for (i=0;i<pep->nconv;i++) {
1107: #if !defined(PETSC_USE_COMPLEX)
1108: if (pep->eigi[i]!=0.0) { /* first eigenvalue of a complex conjugate pair */
1109: BVGetColumn(pjd->X,i,&v);
1110: BVGetColumn(pjd->X,i+1,&v1);
1111: VecNormalizeComplex(v,v1,PETSC_TRUE,NULL);
1112: BVRestoreColumn(pjd->X,i,&v);
1113: BVRestoreColumn(pjd->X,i+1,&v1);
1114: i++;
1115: } else /* real eigenvalue */
1116: #endif
1117: {
1118: BVNormColumn(pjd->X,i,NORM_2,&norm);
1119: BVScaleColumn(pjd->X,i,1.0/norm);
1120: }
1121: }
1122: MatDestroy(&U);
1123: PetscFree3(Z,wr,w);
1124: return(0);
1125: #endif
1126: }
1128: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv,PetscInt sz)
1129: {
1131: PEP_JD *pjd = (PEP_JD*)pep->data;
1132: PetscInt j,i,*P,ldds,rk=0,nvv=*nv;
1133: Vec v,x,w;
1134: PetscScalar *R,*r,*pX,target[2];
1135: Mat X;
1136: PetscBLASInt sz_,rk_,nv_,info;
1137: PetscMPIInt np;
1140: /* update AX and XpX */
1141: for (i=sz;i>0;i--) {
1142: BVGetColumn(pjd->X,pjd->nlock-i,&x);
1143: for (j=0;j<pep->nmat;j++) {
1144: BVGetColumn(pjd->AX[j],pjd->nlock-i,&v);
1145: MatMult(pep->A[j],x,v);
1146: BVRestoreColumn(pjd->AX[j],pjd->nlock-i,&v);
1147: BVSetActiveColumns(pjd->AX[j],0,pjd->nlock-i+1);
1148: }
1149: BVRestoreColumn(pjd->X,pjd->nlock-i,&x);
1150: BVDotColumn(pjd->X,(pjd->nlock-i),pjd->XpX+(pjd->nlock-i)*(pjd->ld));
1151: pjd->XpX[(pjd->nlock-i)*(1+pjd->ld)] = 1.0;
1152: 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]);
1153: }
1155: /* minimality index */
1156: pjd->midx = PetscMin(pjd->mmidx,pjd->nlock);
1158: /* evaluate the polynomial basis in T */
1159: PetscMemzero(pjd->Tj,pjd->ld*pjd->ld*pep->nmat*sizeof(PetscScalar));
1160: for (j=0;j<pep->nmat;j++) {
1161: 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);
1162: }
1164: /* Extend search space */
1165: MPI_Comm_size(PETSC_COMM_WORLD,&np);
1166: PetscCalloc3(nvv,&P,nvv*nvv,&R,nvv*sz,&r);
1167: DSGetLeadingDimension(pep->ds,&ldds);
1168: DSGetArray(pep->ds,DS_MAT_X,&pX);
1169: PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
1170: for (j=0;j<sz;j++) {
1171: 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 */
1172: }
1173: PetscBLASIntCast(rk,&rk_);
1174: PetscBLASIntCast(sz,&sz_);
1175: PetscBLASIntCast(nvv,&nv_);
1176: PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
1177: SlepcCheckLapackInfo("trtri",info);
1178: for (i=0;i<sz;i++) {
1179: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r+i,&sz_));
1180: }
1181: for (i=0;i<sz*rk;i++) r[i] = PetscConj(r[i])/PetscSqrtReal(np); /* revert */
1182: BVSetActiveColumns(pjd->V,0,nvv);
1183: rk -= sz;
1184: for (j=0;j<rk;j++) {
1185: PetscMemcpy(R+j*nvv,pX+(j+sz)*ldds,nvv*sizeof(PetscScalar));
1186: }
1187: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1188: MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk,R,&X);
1189: BVMultInPlace(pjd->V,X,0,rk);
1190: MatDestroy(&X);
1191: BVSetActiveColumns(pjd->V,0,rk);
1192: for (j=0;j<rk;j++) {
1193: BVGetColumn(pjd->V,j,&v);
1194: PEPJDCopyToExtendedVec(pep,NULL,r+sz*(j+sz),sz,pjd->nlock-sz,v,PETSC_FALSE);
1195: BVRestoreColumn(pjd->V,j,&v);
1196: }
1197: BVOrthogonalize(pjd->V,NULL);
1199: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1200: for (j=0;j<rk;j++) {
1201: /* W = P(target)*V */
1202: BVGetColumn(pjd->W,j,&w);
1203: BVGetColumn(pjd->V,j,&v);
1204: target[0] = pep->target; target[1] = 0.0;
1205: PEPJDComputeResidual(pep,PETSC_FALSE,1,&v,target,&w,pep->work);
1206: BVRestoreColumn(pjd->V,j,&v);
1207: BVRestoreColumn(pjd->W,j,&w);
1208: }
1209: BVSetActiveColumns(pjd->W,0,rk);
1210: BVOrthogonalize(pjd->W,NULL);
1211: }
1212: *nv = rk;
1213: PetscFree3(P,R,r);
1214: return(0);
1215: }
1217: PetscErrorCode PEPJDSystemSetUp(PEP pep,PetscInt sz,PetscScalar *theta,Vec *u,Vec *p,Vec *ww)
1218: {
1219: PetscErrorCode ierr;
1220: PetscScalar s[2];
1221: PEP_JD *pjd = (PEP_JD*)pep->data;
1222: PEP_JD_PCSHELL *pcctx;
1225: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1226: PEPJDMatSetUp(pep,sz,theta);
1227: pcctx->u[0] = u[0]; pcctx->u[1] = u[1];
1228: /* Compute r'. p is a work space vector */
1229: PEPJDComputeResidual(pep,PETSC_TRUE,sz,u,theta,p,ww);
1230: PEPJDExtendedPCApply(pjd->pcshell,p[0],pcctx->Bp[0]);
1231: VecDot(pcctx->Bp[0],u[0],pcctx->gamma);
1232: if (sz==2) {
1233: PEPJDExtendedPCApply(pjd->pcshell,p[1],pcctx->Bp[1]);
1234: VecDot(pcctx->Bp[0],u[1],pcctx->gamma+1);
1235: VecMDot(pcctx->Bp[1],2,u,s);
1236: pcctx->gamma[0] += s[1];
1237: pcctx->gamma[1] = -pcctx->gamma[1]+s[0];
1238: } else {
1239: VecZeroEntries(pcctx->Bp[1]);
1240: pcctx->gamma[1] = 0.0;
1241: }
1242: return(0);
1243: }
1245: PetscErrorCode PEPSolve_JD(PEP pep)
1246: {
1247: PetscErrorCode ierr;
1248: PEP_JD *pjd = (PEP_JD*)pep->data;
1249: PetscInt k,nv,nvc,ld,minv,dim,bupdated=0,sz=1,kspsf=1,idx,off,maxits,nloc;
1250: PetscMPIInt np,count;
1251: PetscScalar theta[2]={0.0,0.0},ritz[2]={0.0,0.0},*pX,*eig,*eigi,*array;
1252: PetscReal norm,norm1,*res,tol=0.0,rtol,abstol, dtol;
1253: PetscBool lindep,ini=PETSC_TRUE;
1254: Vec tc,t[2]={NULL,NULL},u[2]={NULL,NULL},p[2]={NULL,NULL};
1255: Vec rc,rr[2],r[2]={NULL,NULL},*ww=pep->work,v[2];
1256: Mat G,X,Y;
1257: KSP ksp;
1258: PEP_JD_PCSHELL *pcctx;
1259: PEP_JD_MATSHELL *matctx;
1262: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1263: BVGetSizes(pep->V,&nloc,NULL,NULL);
1264: DSGetLeadingDimension(pep->ds,&ld);
1265: PetscCalloc3(pep->ncv+pep->nev,&eig,pep->ncv+pep->nev,&eigi,pep->ncv+pep->nev,&res);
1266: pjd->nlock = 0;
1267: STGetKSP(pep->st,&ksp);
1268: KSPGetTolerances(ksp,&rtol,&abstol,&dtol,&maxits);
1269: #if !defined (PETSC_USE_COMPLEX)
1270: kspsf = 2;
1271: #endif
1272: PEPJDProcessInitialSpace(pep,ww);
1273: nv = (pep->nini)?pep->nini:1;
1275: /* Replace preconditioner with one containing projectors */
1276: PEPJDCreateShellPC(pep,ww);
1277: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1279: /* Create auxiliar vectors */
1280: BVCreateVec(pjd->V,&u[0]);
1281: VecDuplicate(u[0],&p[0]);
1282: VecDuplicate(u[0],&r[0]);
1283: #if !defined (PETSC_USE_COMPLEX)
1284: VecDuplicate(u[0],&u[1]);
1285: VecDuplicate(u[0],&p[1]);
1286: VecDuplicate(u[0],&r[1]);
1287: #endif
1289: /* Restart loop */
1290: while (pep->reason == PEP_CONVERGED_ITERATING) {
1291: pep->its++;
1292: DSSetDimensions(pep->ds,nv,0,0,0);
1293: BVSetActiveColumns(pjd->V,bupdated,nv);
1294: PEPJDUpdateTV(pep,bupdated,nv,ww);
1295: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVSetActiveColumns(pjd->W,bupdated,nv); }
1296: for (k=0;k<pep->nmat;k++) {
1297: BVSetActiveColumns(pjd->TV[k],bupdated,nv);
1298: DSGetMat(pep->ds,DSMatExtra[k],&G);
1299: BVMatProject(pjd->TV[k],NULL,pjd->W,G);
1300: DSRestoreMat(pep->ds,DSMatExtra[k],&G);
1301: }
1302: BVSetActiveColumns(pjd->V,0,nv);
1303: BVSetActiveColumns(pjd->W,0,nv);
1305: /* Solve projected problem */
1306: DSSetState(pep->ds,DS_STATE_RAW);
1307: DSSolve(pep->ds,pep->eigr,pep->eigi);
1308: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1309: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
1310: idx = 0;
1311: do {
1312: ritz[0] = pep->eigr[idx];
1313: #if !defined(PETSC_USE_COMPLEX)
1314: ritz[1] = pep->eigi[idx];
1315: #endif
1316: sz = (ritz[1]==0.0)?1:2;
1317: /* Compute Ritz vector u=V*X(:,1) */
1318: DSGetArray(pep->ds,DS_MAT_X,&pX);
1319: BVSetActiveColumns(pjd->V,0,nv);
1320: BVMultVec(pjd->V,1.0,0.0,u[0],pX+idx*ld);
1321: if (sz==2) {
1322: BVMultVec(pjd->V,1.0,0.0,u[1],pX+(idx+1)*ld);
1323: }
1324: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1325: PEPJDComputeResidual(pep,PETSC_FALSE,sz,u,ritz,r,ww);
1326: /* Check convergence */
1327: VecNorm(r[0],NORM_2,&norm);
1328: if (sz==2) {
1329: VecNorm(r[1],NORM_2,&norm1);
1330: norm = SlepcAbs(norm,norm1);
1331: }
1332: (*pep->converged)(pep,ritz[0],ritz[1],norm,&pep->errest[pep->nconv],pep->convergedctx);
1333: if (sz==2) pep->errest[pep->nconv+1] = pep->errest[pep->nconv];
1334: if (ini) {
1335: tol = PetscMin(.1,pep->errest[pep->nconv]); ini = PETSC_FALSE;
1336: } else tol = PetscMin(pep->errest[pep->nconv],tol/2);
1337: (*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);
1338: if (pep->errest[pep->nconv]<pep->tol) {
1339: /* Ritz pair converged */
1340: ini = PETSC_TRUE;
1341: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1342: if (pjd->ld>1) {
1343: BVGetColumn(pjd->X,pep->nconv,&v[0]);
1344: PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*pep->nconv,pjd->ld-1,0,u[0],PETSC_TRUE);
1345: BVRestoreColumn(pjd->X,pep->nconv,&v[0]);
1346: BVSetActiveColumns(pjd->X,0,pep->nconv+1);
1347: BVNormColumn(pjd->X,pep->nconv,NORM_2,&norm);
1348: BVScaleColumn(pjd->X,pep->nconv,1.0/norm);
1349: for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*pep->nconv+k] *= PetscSqrtReal(np)/norm;
1350: pjd->T[(pep->ncv+1)*pep->nconv] = ritz[0];
1351: eig[pep->nconv] = ritz[0];
1352: idx++;
1353: if (sz==2) {
1354: BVGetColumn(pjd->X,pep->nconv+1,&v[0]);
1355: PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*(pep->nconv+1),pjd->ld-1,0,u[1],PETSC_TRUE);
1356: BVRestoreColumn(pjd->X,pep->nconv+1,&v[0]);
1357: BVSetActiveColumns(pjd->X,0,pep->nconv+2);
1358: BVNormColumn(pjd->X,pep->nconv+1,NORM_2,&norm1);
1359: BVScaleColumn(pjd->X,pep->nconv+1,1.0/norm1);
1360: for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*(pep->nconv+1)+k] *= PetscSqrtReal(np)/norm1;
1361: pjd->T[(pep->ncv+1)*(pep->nconv+1)] = ritz[0];
1362: pjd->T[(pep->ncv+1)*pep->nconv+1] = -ritz[1]*norm1/norm;
1363: pjd->T[(pep->ncv+1)*(pep->nconv+1)-1] = ritz[1]*norm/norm1;
1364: eig[pep->nconv+1] = ritz[0];
1365: eigi[pep->nconv] = ritz[1]; eigi[pep->nconv+1] = -ritz[1];
1366: idx++;
1367: }
1368: } else {
1369: BVInsertVec(pep->V,pep->nconv,u[0]);
1370: }
1371: pep->nconv += sz;
1372: }
1373: } while (pep->errest[pep->nconv]<pep->tol && pep->nconv<nv);
1375: if (pep->reason==PEP_CONVERGED_ITERATING) {
1376: nvc = nv;
1377: if (idx) {
1378: pjd->nlock +=idx;
1379: PEPJDLockConverged(pep,&nv,idx);
1380: }
1381: if (nv+sz>=pep->ncv-1) {
1382: /* Basis full, force restart */
1383: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1384: DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
1385: DSGetArray(pep->ds,DS_MAT_X,&pX);
1386: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1387: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1388: DSGetArray(pep->ds,DS_MAT_Y,&pX);
1389: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1390: DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
1391: DSGetMat(pep->ds,DS_MAT_X,&X);
1392: BVMultInPlace(pjd->V,X,0,minv);
1393: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1394: DSGetMat(pep->ds,DS_MAT_Y,&Y);
1395: BVMultInPlace(pjd->W,Y,pep->nconv,minv);
1396: DSRestoreMat(pep->ds,DS_MAT_Y,&Y);
1397: }
1398: MatDestroy(&X);
1399: nv = minv;
1400: bupdated = 0;
1401: } else {
1402: if (!idx && pep->errest[pep->nconv]<pjd->fix) {theta[0] = ritz[0]; theta[1] = ritz[1];}
1403: else {theta[0] = pep->target; theta[1] = 0.0;}
1404: /* Update system mat */
1405: PEPJDSystemSetUp(pep,sz,theta,u,p,ww);
1406: /* Solve correction equation to expand basis */
1407: BVGetColumn(pjd->V,nv,&t[0]);
1408: rr[0] = r[0];
1409: if (sz==2) {
1410: BVGetColumn(pjd->V,nv+1,&t[1]);
1411: rr[1] = r[1];
1412: } else {
1413: t[1] = NULL;
1414: rr[1] = NULL;
1415: }
1416: VecCreateCompWithVecs(t,kspsf,pjd->vtempl,&tc);
1417: VecCreateCompWithVecs(rr,kspsf,pjd->vtempl,&rc);
1418: VecCompSetSubVecs(pjd->vtempl,sz,NULL);
1419: tol = PetscMax(rtol,tol/2);
1420: KSPSetTolerances(ksp,tol,abstol,dtol,maxits);
1421: KSPSolve(ksp,rc,tc);
1422: VecDestroy(&tc);
1423: VecDestroy(&rc);
1424: VecGetArray(t[0],&array);
1425: PetscMPIIntCast(pep->nconv,&count);
1426: MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1427: VecRestoreArray(t[0],&array);
1428: BVRestoreColumn(pjd->V,nv,&t[0]);
1429: BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
1430: if (lindep || norm==0.0) {
1431: if (sz==1) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1432: else off = 1;
1433: } else {
1434: off = 0;
1435: BVScaleColumn(pjd->V,nv,1.0/norm);
1436: }
1437: if (sz==2) {
1438: VecGetArray(t[1],&array);
1439: MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1440: VecRestoreArray(t[1],&array);
1441: BVRestoreColumn(pjd->V,nv+1,&t[1]);
1442: if (off) {
1443: BVCopyColumn(pjd->V,nv+1,nv);
1444: }
1445: BVOrthogonalizeColumn(pjd->V,nv+1-off,NULL,&norm,&lindep);
1446: if (lindep || norm==0.0) {
1447: if (off) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1448: else off = 1;
1449: } else {
1450: BVScaleColumn(pjd->V,nv+1-off,1.0/norm);
1451: }
1452: }
1453: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1454: BVInsertVec(pjd->W,nv,r[0]);
1455: if (sz==2 && !off) {
1456: BVInsertVec(pjd->W,nv+1,r[1]);
1457: }
1458: BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
1459: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1460: BVScaleColumn(pjd->W,nv,1.0/norm);
1461: if (sz==2 && !off) {
1462: BVOrthogonalizeColumn(pjd->W,nv+1,NULL,&norm,&lindep);
1463: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1464: BVScaleColumn(pjd->W,nv+1,1.0/norm);
1465: }
1466: }
1467: bupdated = idx?0:nv;
1468: nv += sz-off;
1469: }
1470: for (k=0;k<nvc;k++) {
1471: eig[pep->nconv-idx+k] = pep->eigr[k];
1472: #if !defined(PETSC_USE_COMPLEX)
1473: eigi[pep->nconv-idx+k] = pep->eigi[k];
1474: #endif
1475: }
1476: PEPMonitor(pep,pep->its,pep->nconv,eig,eigi,pep->errest,pep->nconv+1);
1477: }
1478: }
1479: if (pjd->ld>1) {
1480: for (k=0;k<pep->nconv;k++) {
1481: pep->eigr[k] = eig[k];
1482: pep->eigi[k] = eigi[k];
1483: }
1484: if (pep->nconv>0) { PEPJDEigenvectors(pep); }
1485: PetscFree2(pcctx->M,pcctx->ps);
1486: }
1487: VecDestroy(&u[0]);
1488: VecDestroy(&r[0]);
1489: VecDestroy(&p[0]);
1490: #if !defined (PETSC_USE_COMPLEX)
1491: VecDestroy(&u[1]);
1492: VecDestroy(&r[1]);
1493: VecDestroy(&p[1]);
1494: #endif
1495: KSPSetTolerances(ksp,rtol,abstol,dtol,maxits);
1496: KSPSetPC(ksp,pcctx->pc);
1497: VecDestroy(&pcctx->Bp[0]);
1498: VecDestroy(&pcctx->Bp[1]);
1499: MatShellGetContext(pjd->Pshell,(void**)&matctx);
1500: MatDestroy(&matctx->Pr);
1501: MatDestroy(&matctx->Pi);
1502: MatDestroy(&pjd->Pshell);
1503: MatDestroy(&pcctx->PPr);
1504: PCDestroy(&pcctx->pc);
1505: PetscFree(pcctx);
1506: PetscFree(matctx);
1507: PCDestroy(&pjd->pcshell);
1508: PetscFree3(eig,eigi,res);
1509: VecDestroy(&pjd->vtempl);
1510: return(0);
1511: }
1513: PetscErrorCode PEPJDSetRestart_JD(PEP pep,PetscReal keep)
1514: {
1515: PEP_JD *pjd = (PEP_JD*)pep->data;
1518: if (keep==PETSC_DEFAULT) pjd->keep = 0.5;
1519: else {
1520: 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]");
1521: pjd->keep = keep;
1522: }
1523: return(0);
1524: }
1526: /*@
1527: PEPJDSetRestart - Sets the restart parameter for the Jacobi-Davidson
1528: method, in particular the proportion of basis vectors that must be kept
1529: after restart.
1531: Logically Collective on PEP
1533: Input Parameters:
1534: + pep - the eigenproblem solver context
1535: - keep - the number of vectors to be kept at restart
1537: Options Database Key:
1538: . -pep_jd_restart - Sets the restart parameter
1540: Notes:
1541: Allowed values are in the range [0.1,0.9]. The default is 0.5.
1543: Level: advanced
1545: .seealso: PEPJDGetRestart()
1546: @*/
1547: PetscErrorCode PEPJDSetRestart(PEP pep,PetscReal keep)
1548: {
1554: PetscTryMethod(pep,"PEPJDSetRestart_C",(PEP,PetscReal),(pep,keep));
1555: return(0);
1556: }
1558: PetscErrorCode PEPJDGetRestart_JD(PEP pep,PetscReal *keep)
1559: {
1560: PEP_JD *pjd = (PEP_JD*)pep->data;
1563: *keep = pjd->keep;
1564: return(0);
1565: }
1567: /*@
1568: PEPJDGetRestart - Gets the restart parameter used in the Jacobi-Davidson method.
1570: Not Collective
1572: Input Parameter:
1573: . pep - the eigenproblem solver context
1575: Output Parameter:
1576: . keep - the restart parameter
1578: Level: advanced
1580: .seealso: PEPJDSetRestart()
1581: @*/
1582: PetscErrorCode PEPJDGetRestart(PEP pep,PetscReal *keep)
1583: {
1589: PetscUseMethod(pep,"PEPJDGetRestart_C",(PEP,PetscReal*),(pep,keep));
1590: return(0);
1591: }
1593: PetscErrorCode PEPJDSetFix_JD(PEP pep,PetscReal fix)
1594: {
1595: PEP_JD *pjd = (PEP_JD*)pep->data;
1598: if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) pjd->fix = 0.01;
1599: else {
1600: if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
1601: pjd->fix = fix;
1602: }
1603: return(0);
1604: }
1606: /*@
1607: PEPJDSetFix - Sets the threshold for changing the target in the correction
1608: equation.
1610: Logically Collective on PEP
1612: Input Parameters:
1613: + pep - the eigenproblem solver context
1614: - fix - threshold for changing the target
1616: Options Database Key:
1617: . -pep_jd_fix - the fix value
1619: Note:
1620: The target in the correction equation is fixed at the first iterations.
1621: When the norm of the residual vector is lower than the fix value,
1622: the target is set to the corresponding eigenvalue.
1624: Level: advanced
1626: .seealso: PEPJDGetFix()
1627: @*/
1628: PetscErrorCode PEPJDSetFix(PEP pep,PetscReal fix)
1629: {
1635: PetscTryMethod(pep,"PEPJDSetFix_C",(PEP,PetscReal),(pep,fix));
1636: return(0);
1637: }
1639: PetscErrorCode PEPJDGetFix_JD(PEP pep,PetscReal *fix)
1640: {
1641: PEP_JD *pjd = (PEP_JD*)pep->data;
1644: *fix = pjd->fix;
1645: return(0);
1646: }
1648: /*@
1649: PEPJDGetFix - Returns the threshold for changing the target in the correction
1650: equation.
1652: Not Collective
1654: Input Parameter:
1655: . pep - the eigenproblem solver context
1657: Output Parameter:
1658: . fix - threshold for changing the target
1660: Note:
1661: The target in the correction equation is fixed at the first iterations.
1662: When the norm of the residual vector is lower than the fix value,
1663: the target is set to the corresponding eigenvalue.
1665: Level: advanced
1667: .seealso: PEPJDSetFix()
1668: @*/
1669: PetscErrorCode PEPJDGetFix(PEP pep,PetscReal *fix)
1670: {
1676: PetscUseMethod(pep,"PEPJDGetFix_C",(PEP,PetscReal*),(pep,fix));
1677: return(0);
1678: }
1680: PetscErrorCode PEPJDSetReusePreconditioner_JD(PEP pep,PetscBool reusepc)
1681: {
1682: PEP_JD *pjd = (PEP_JD*)pep->data;
1685: pjd->reusepc = reusepc;
1686: return(0);
1687: }
1689: /*@
1690: PEPJDSetReusePreconditioner - Sets a flag indicating whether the preconditioner
1691: must be reused or not.
1693: Logically Collective on PEP
1695: Input Parameters:
1696: + pep - the eigenproblem solver context
1697: - reusepc - the reuse flag
1699: Options Database Key:
1700: . -pep_jd_reuse_preconditioner - the reuse flag
1702: Note:
1703: The default value is False. If set to True, the preconditioner is built
1704: only at the beginning, using the target value. Otherwise, it may be rebuilt
1705: (depending on the fix parameter) at each iteration from the Ritz value.
1707: Level: advanced
1709: .seealso: PEPJDGetReusePreconditioner(), PEPJDSetFix()
1710: @*/
1711: PetscErrorCode PEPJDSetReusePreconditioner(PEP pep,PetscBool reusepc)
1712: {
1718: PetscTryMethod(pep,"PEPJDSetReusePreconditioner_C",(PEP,PetscBool),(pep,reusepc));
1719: return(0);
1720: }
1722: PetscErrorCode PEPJDGetReusePreconditioner_JD(PEP pep,PetscBool *reusepc)
1723: {
1724: PEP_JD *pjd = (PEP_JD*)pep->data;
1727: *reusepc = pjd->reusepc;
1728: return(0);
1729: }
1731: /*@
1732: PEPJDGetReusePreconditioner - Returns the flag for reusing the preconditioner.
1734: Not Collective
1736: Input Parameter:
1737: . pep - the eigenproblem solver context
1739: Output Parameter:
1740: . reusepc - the reuse flag
1742: Level: advanced
1744: .seealso: PEPJDSetReusePreconditioner()
1745: @*/
1746: PetscErrorCode PEPJDGetReusePreconditioner(PEP pep,PetscBool *reusepc)
1747: {
1753: PetscUseMethod(pep,"PEPJDGetReusePreconditioner_C",(PEP,PetscBool*),(pep,reusepc));
1754: return(0);
1755: }
1757: PetscErrorCode PEPJDSetMinimalityIndex_JD(PEP pep,PetscInt mmidx)
1758: {
1759: PEP_JD *pjd = (PEP_JD*)pep->data;
1762: if (mmidx == PETSC_DEFAULT || mmidx == PETSC_DECIDE) pjd->mmidx = 1;
1763: else {
1764: if (mmidx < 1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mmidx value");
1765: pjd->mmidx = mmidx;
1766: pep->state = PEP_STATE_INITIAL;
1767: }
1768: return(0);
1769: }
1771: /*@
1772: PEPJDSetMinimalityIndex - Sets the maximum allowed value for the minimality index.
1774: Logically Collective on PEP
1776: Input Parameters:
1777: + pep - the eigenproblem solver context
1778: - mmidx - maximum minimality index
1780: Options Database Key:
1781: . -pep_jd_minimality_index - the minimality index value
1783: Note:
1784: The default value is equal to the degree of the polynomial. A smaller value
1785: can be used if the wanted eigenvectors are known to be linearly independent.
1787: Level: advanced
1789: .seealso: PEPJDGetMinimalityIndex()
1790: @*/
1791: PetscErrorCode PEPJDSetMinimalityIndex(PEP pep,PetscInt mmidx)
1792: {
1798: PetscTryMethod(pep,"PEPJDSetMinimalityIndex_C",(PEP,PetscInt),(pep,mmidx));
1799: return(0);
1800: }
1802: PetscErrorCode PEPJDGetMinimalityIndex_JD(PEP pep,PetscInt *mmidx)
1803: {
1804: PEP_JD *pjd = (PEP_JD*)pep->data;
1807: *mmidx = pjd->mmidx;
1808: return(0);
1809: }
1811: /*@
1812: PEPJDGetMinimalityIndex - Returns the maximum allowed value of the minimality
1813: index.
1815: Not Collective
1817: Input Parameter:
1818: . pep - the eigenproblem solver context
1820: Output Parameter:
1821: . mmidx - minimality index
1823: Level: advanced
1825: .seealso: PEPJDSetMinimalityIndex()
1826: @*/
1827: PetscErrorCode PEPJDGetMinimalityIndex(PEP pep,PetscInt *mmidx)
1828: {
1834: PetscUseMethod(pep,"PEPJDGetMinimalityIndex_C",(PEP,PetscInt*),(pep,mmidx));
1835: return(0);
1836: }
1838: PetscErrorCode PEPJDSetProjection_JD(PEP pep,PEPJDProjection proj)
1839: {
1840: PEP_JD *pjd = (PEP_JD*)pep->data;
1843: switch (proj) {
1844: case PEP_JD_PROJECTION_HARMONIC:
1845: case PEP_JD_PROJECTION_ORTHOGONAL:
1846: if (pjd->proj != proj) {
1847: pep->state = PEP_STATE_INITIAL;
1848: pjd->proj = proj;
1849: }
1850: break;
1851: default:
1852: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'proj' value");
1853: }
1854: return(0);
1855: }
1857: /*@
1858: PEPJDSetProjection - Sets the type of projection to be used in the Jacobi-Davidson solver.
1860: Logically Collective on PEP
1862: Input Parameters:
1863: + pep - the eigenproblem solver context
1864: - proj - the type of projection
1866: Options Database Key:
1867: . -pep_jd_projection - the projection type, either orthogonal or harmonic
1869: Level: advanced
1871: .seealso: PEPJDGetProjection()
1872: @*/
1873: PetscErrorCode PEPJDSetProjection(PEP pep,PEPJDProjection proj)
1874: {
1880: PetscTryMethod(pep,"PEPJDSetProjection_C",(PEP,PEPJDProjection),(pep,proj));
1881: return(0);
1882: }
1884: PetscErrorCode PEPJDGetProjection_JD(PEP pep,PEPJDProjection *proj)
1885: {
1886: PEP_JD *pjd = (PEP_JD*)pep->data;
1889: *proj = pjd->proj;
1890: return(0);
1891: }
1893: /*@
1894: PEPJDGetProjection - Returns the type of projection used by the Jacobi-Davidson solver.
1896: Not Collective
1898: Input Parameter:
1899: . pep - the eigenproblem solver context
1901: Output Parameter:
1902: . proj - the type of projection
1904: Level: advanced
1906: .seealso: PEPJDSetProjection()
1907: @*/
1908: PetscErrorCode PEPJDGetProjection(PEP pep,PEPJDProjection *proj)
1909: {
1915: PetscUseMethod(pep,"PEPJDGetProjection_C",(PEP,PEPJDProjection*),(pep,proj));
1916: return(0);
1917: }
1919: PetscErrorCode PEPSetFromOptions_JD(PetscOptionItems *PetscOptionsObject,PEP pep)
1920: {
1921: PetscErrorCode ierr;
1922: PetscBool flg,b1;
1923: PetscReal r1;
1924: PetscInt i1;
1925: PEPJDProjection proj;
1928: PetscOptionsHead(PetscOptionsObject,"PEP JD Options");
1930: PetscOptionsReal("-pep_jd_restart","Proportion of vectors kept after restart","PEPJDSetRestart",0.5,&r1,&flg);
1931: if (flg) { PEPJDSetRestart(pep,r1); }
1933: PetscOptionsReal("-pep_jd_fix","Tolerance for changing the target in the correction equation","PEPJDSetFix",0.01,&r1,&flg);
1934: if (flg) { PEPJDSetFix(pep,r1); }
1936: PetscOptionsBool("-pep_jd_reuse_preconditioner","Whether to reuse the preconditioner","PEPJDSetReusePreconditoiner",PETSC_FALSE,&b1,&flg);
1937: if (flg) { PEPJDSetReusePreconditioner(pep,b1); }
1939: PetscOptionsInt("-pep_jd_minimality_index","Maximum allowed minimality index","PEPJDSetMinimalityIndex",1,&i1,&flg);
1940: if (flg) { PEPJDSetMinimalityIndex(pep,i1); }
1942: PetscOptionsEnum("-pep_jd_projection","Type of projection","PEPJDSetProjection",PEPJDProjectionTypes,(PetscEnum)PEP_JD_PROJECTION_HARMONIC,(PetscEnum*)&proj,&flg);
1943: if (flg) { PEPJDSetProjection(pep,proj); }
1945: PetscOptionsTail();
1946: return(0);
1947: }
1949: PetscErrorCode PEPView_JD(PEP pep,PetscViewer viewer)
1950: {
1952: PEP_JD *pjd = (PEP_JD*)pep->data;
1953: PetscBool isascii;
1956: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1957: if (isascii) {
1958: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*pjd->keep));
1959: PetscViewerASCIIPrintf(viewer," threshold for changing the target in the correction equation (fix): %g\n",(double)pjd->fix);
1960: PetscViewerASCIIPrintf(viewer," projection type: %s\n",PEPJDProjectionTypes[pjd->proj]);
1961: PetscViewerASCIIPrintf(viewer," maximum allowed minimality index: %d\n",pjd->mmidx);
1962: if (pjd->reusepc) { PetscViewerASCIIPrintf(viewer," reusing the preconditioner\n"); }
1963: }
1964: return(0);
1965: }
1967: PetscErrorCode PEPSetDefaultST_JD(PEP pep)
1968: {
1970: KSP ksp;
1973: if (!((PetscObject)pep->st)->type_name) {
1974: STSetType(pep->st,STPRECOND);
1975: STPrecondSetKSPHasMat(pep->st,PETSC_TRUE);
1976: }
1977: STSetTransform(pep->st,PETSC_FALSE);
1978: STGetKSP(pep->st,&ksp);
1979: if (!((PetscObject)ksp)->type_name) {
1980: KSPSetType(ksp,KSPBCGSL);
1981: KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
1982: }
1983: return(0);
1984: }
1986: PetscErrorCode PEPReset_JD(PEP pep)
1987: {
1989: PEP_JD *pjd = (PEP_JD*)pep->data;
1990: PetscInt i;
1993: for (i=0;i<pep->nmat;i++) {
1994: BVDestroy(pjd->TV+i);
1995: }
1996: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVDestroy(&pjd->W); }
1997: if (pjd->ld>1) {
1998: BVDestroy(&pjd->V);
1999: for (i=0;i<pep->nmat;i++) {
2000: BVDestroy(pjd->AX+i);
2001: }
2002: BVDestroy(&pjd->N[0]);
2003: BVDestroy(&pjd->N[1]);
2004: PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
2005: }
2006: PetscFree2(pjd->TV,pjd->AX);
2007: return(0);
2008: }
2010: PetscErrorCode PEPDestroy_JD(PEP pep)
2011: {
2015: PetscFree(pep->data);
2016: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
2017: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
2018: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",NULL);
2019: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",NULL);
2020: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",NULL);
2021: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",NULL);
2022: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",NULL);
2023: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",NULL);
2024: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",NULL);
2025: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",NULL);
2026: return(0);
2027: }
2029: SLEPC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
2030: {
2031: PEP_JD *pjd;
2035: PetscNewLog(pep,&pjd);
2036: pep->data = (void*)pjd;
2038: pjd->fix = 0.01;
2039: pjd->mmidx = 0;
2041: pep->ops->solve = PEPSolve_JD;
2042: pep->ops->setup = PEPSetUp_JD;
2043: pep->ops->setfromoptions = PEPSetFromOptions_JD;
2044: pep->ops->destroy = PEPDestroy_JD;
2045: pep->ops->reset = PEPReset_JD;
2046: pep->ops->view = PEPView_JD;
2047: pep->ops->setdefaultst = PEPSetDefaultST_JD;
2049: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
2050: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
2051: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",PEPJDSetFix_JD);
2052: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",PEPJDGetFix_JD);
2053: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",PEPJDSetReusePreconditioner_JD);
2054: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",PEPJDGetReusePreconditioner_JD);
2055: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",PEPJDSetMinimalityIndex_JD);
2056: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",PEPJDGetMinimalityIndex_JD);
2057: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",PEPJDSetProjection_JD);
2058: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",PEPJDGetProjection_JD);
2059: return(0);
2060: }