Actual source code: pjd.c
slepc-3.12.0 2019-09-30
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;
124: 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");
126: PetscObjectTypeCompare((PetscObject)pep->st,STPRECOND,&isprecond);
127: if (!isprecond) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"JD only works with PRECOND spectral transformation");
129: STGetTransform(pep->st,&flg);
130: if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");
132: if (!pjd->mmidx) pjd->mmidx = pep->nmat-1;
133: pjd->mmidx = PetscMin(pjd->mmidx,pep->nmat-1);
134: if (!pjd->keep) pjd->keep = 0.5;
135: PEPBasisCoefficients(pep,pep->pbc);
136: PEPAllocateSolution(pep,0);
137: PEPSetWorkVecs(pep,5);
138: pjd->ld = pep->nev;
139: #if !defined (PETSC_USE_COMPLEX)
140: pjd->ld++;
141: #endif
142: PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
143: for (i=0;i<pep->nmat;i++) {
144: PEPJDDuplicateBasis(pep,pjd->TV+i);
145: }
146: if (pjd->ld>1) {
147: PEPJDDuplicateBasis(pep,&pjd->V);
148: BVSetFromOptions(pjd->V);
149: for (i=0;i<pep->nmat;i++) {
150: BVDuplicateResize(pep->V,pjd->ld-1,pjd->AX+i);
151: }
152: BVDuplicateResize(pep->V,pjd->ld-1,pjd->N);
153: BVDuplicateResize(pep->V,pjd->ld-1,pjd->N+1);
154: pjd->X = pep->V;
155: PetscCalloc3((pjd->ld)*(pjd->ld),&pjd->XpX,pep->ncv*pep->ncv,&pjd->T,pjd->ld*pjd->ld*pep->nmat,&pjd->Tj);
156: } else pjd->V = pep->V;
157: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { PEPJDDuplicateBasis(pep,&pjd->W); }
158: else pjd->W = pjd->V;
159: DSSetType(pep->ds,DSPEP);
160: DSPEPSetDegree(pep->ds,pep->nmat-1);
161: if (pep->basis!=PEP_BASIS_MONOMIAL) {
162: DSPEPSetCoefficients(pep->ds,pep->pbc);
163: }
164: DSAllocate(pep->ds,pep->ncv);
165: return(0);
166: }
168: /*
169: Updates columns (low to (high-1)) of TV[i]
170: */
171: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
172: {
174: PEP_JD *pjd = (PEP_JD*)pep->data;
175: PetscInt pp,col,i,nloc,nconv;
176: Vec v1,v2,t1,t2;
177: PetscScalar *array1,*array2,*x2,*xx,*N,*Np,*y2=NULL,zero=0.0,sone=1.0,*pT,fact,*psc;
178: PetscReal *cg,*ca,*cb;
179: PetscMPIInt rk,np;
180: PetscBLASInt n_,ld_,one=1;
181: Mat T;
182: BV pbv;
185: ca = pep->pbc; cb = ca+pep->nmat; cg = cb + pep->nmat;
186: nconv = pjd->nlock;
187: PetscMalloc5(nconv,&x2,nconv,&xx,nconv*nconv,&pT,nconv*nconv,&N,nconv*nconv,&Np);
188: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
189: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
190: BVGetSizes(pep->V,&nloc,NULL,NULL);
191: t1 = w[0];
192: t2 = w[1];
193: PetscBLASIntCast(pjd->nlock,&n_);
194: PetscBLASIntCast(pjd->ld,&ld_);
195: if (nconv){
196: for (i=0;i<nconv;i++) {
197: PetscArraycpy(pT+i*nconv,pjd->T+i*pep->ncv,nconv);
198: }
199: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,pT,&T);
200: }
201: for (col=low;col<high;col++) {
202: BVGetColumn(pjd->V,col,&v1);
203: VecGetArray(v1,&array1);
204: if (nconv>0) {
205: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
206: }
207: VecPlaceArray(t1,array1);
208: if (nconv) {
209: BVSetActiveColumns(pjd->N[0],0,nconv);
210: BVSetActiveColumns(pjd->N[1],0,nconv);
211: BVDotVec(pjd->X,t1,xx);
212: }
213: for (pp=pep->nmat-1;pp>=0;pp--) {
214: BVGetColumn(pjd->TV[pp],col,&v2);
215: VecGetArray(v2,&array2);
216: VecPlaceArray(t2,array2);
217: MatMult(pep->A[pp],t1,t2);
218: if (nconv) {
219: if (pp<pep->nmat-3) {
220: BVMult(pjd->N[0],1.0,-cg[pp+2],pjd->AX[pp+1],NULL);
221: MatShift(T,-cb[pp+1]);
222: BVMult(pjd->N[0],1.0/ca[pp],1.0/ca[pp],pjd->N[1],T);
223: pbv = pjd->N[0]; pjd->N[0] = pjd->N[1]; pjd->N[1] = pbv;
224: BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
225: MatShift(T,cb[pp+1]);
226: } else if (pp==pep->nmat-3) {
227: BVCopy(pjd->AX[pp+2],pjd->N[0]);
228: BVScale(pjd->N[0],1/ca[pp+1]);
229: BVCopy(pjd->AX[pp+1],pjd->N[1]);
230: MatShift(T,-cb[pp+1]);
231: BVMult(pjd->N[1],1.0/ca[pp],1.0/ca[pp],pjd->N[0],T);
232: BVMultVec(pjd->N[1],1.0,1.0,t2,x2);
233: MatShift(T,cb[pp+1]);
234: } else if (pp==pep->nmat-2) {
235: BVMultVec(pjd->AX[pp+1],1.0/ca[pp],1.0,t2,x2);
236: }
237: if (pp<pjd->midx) {
238: y2 = array2+nloc;
239: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n_,&n_,&sone,pjd->Tj+pjd->ld*pjd->ld*pp,&ld_,xx,&one,&zero,y2,&one));
240: if (pp<pjd->midx-2) {
241: fact = -cg[pp+2];
242: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&fact,Np,&n_));
243: fact = 1/ca[pp];
244: MatShift(T,-cb[pp+1]);
245: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&fact,N,&n_,pT,&n_,&fact,Np,&n_));
246: MatShift(T,cb[pp+1]);
247: psc = Np; Np = N; N = psc;
248: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
249: } else if (pp==pjd->midx-2) {
250: fact = 1/ca[pp];
251: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&fact,pjd->Tj+(pp+1)*pjd->ld*pjd->ld,&ld_,pjd->XpX,&ld_,&zero,N,&n_));
252: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,N,&n_,x2,&one,&sone,y2,&one));
253: } else if (pp==pjd->midx-1) {
254: PetscArrayzero(Np,nconv*nconv);
255: }
256: }
257: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
258: }
259: VecResetArray(t2);
260: VecRestoreArray(v2,&array2);
261: BVRestoreColumn(pjd->TV[pp],col,&v2);
262: }
263: VecResetArray(t1);
264: VecRestoreArray(v1,&array1);
265: BVRestoreColumn(pjd->V,col,&v1);
266: }
267: if (nconv) {MatDestroy(&T);}
268: PetscFree5(x2,xx,pT,N,Np);
269: return(0);
270: }
272: /*
273: RRQR of X. Xin*P=Xou*R. Rank of R is rk
274: */
275: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
276: {
277: #if defined(SLEPC_MISSING_LAPACK_GEQP3) || defined(PETSC_MISSING_LAPACK_ORGQR)
279: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQP3/QRGQR - Lapack routines are unavailable");
280: #else
282: PetscInt i,j,n,r;
283: PetscBLASInt row_,col_,ldx_,*p,lwork,info,n_;
284: PetscScalar *tau,*work;
285: PetscReal tol,*rwork;
288: PetscBLASIntCast(row,&row_);
289: PetscBLASIntCast(col,&col_);
290: PetscBLASIntCast(ldx,&ldx_);
291: n = PetscMin(row,col);
292: PetscBLASIntCast(n,&n_);
293: lwork = 3*col_+1;
294: PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
295: for (i=1;i<col;i++) p[i] = 0;
296: p[0] = 1;
298: /* rank revealing QR */
299: #if defined(PETSC_USE_COMPLEX)
300: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
301: #else
302: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
303: #endif
304: SlepcCheckLapackInfo("geqp3",info);
305: if (P) for (i=0;i<col;i++) P[i] = p[i]-1;
307: /* rank computation */
308: tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
309: r = 1;
310: for (i=1;i<n;i++) {
311: if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
312: else break;
313: }
314: if (rk) *rk=r;
316: /* copy upper triangular matrix if requested */
317: if (R) {
318: for (i=0;i<r;i++) {
319: PetscArrayzero(R+i*ldr,r);
320: for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
321: }
322: }
323: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
324: SlepcCheckLapackInfo("orgqr",info);
325: PetscFree4(p,tau,work,rwork);
326: return(0);
327: #endif
328: }
330: /*
331: Application of extended preconditioner
332: */
333: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
334: {
335: PetscInt i,j,nloc,n,ld=0;
336: PetscMPIInt np;
337: Vec tx,ty;
338: PEP_JD_PCSHELL *ctx;
339: PetscErrorCode ierr;
340: const PetscScalar *array1;
341: PetscScalar *x2=NULL,*t=NULL,*ps=NULL,*array2,zero=0.0,sone=1.0;
342: PetscBLASInt one=1.0,ld_,n_,ncv_;
343: PEP_JD *pjd=NULL;
346: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
347: PCShellGetContext(pc,(void**)&ctx);
348: n = ctx->n;
349: if (n) {
350: pjd = (PEP_JD*)ctx->pep->data;
351: ps = ctx->ps;
352: ld = pjd->ld;
353: PetscMalloc2(n,&x2,n,&t);
354: VecGetLocalSize(ctx->work[0],&nloc);
355: VecGetArrayRead(x,&array1);
356: for (i=0;i<n;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
357: VecRestoreArrayRead(x,&array1);
358: }
360: /* y = B\x apply PC */
361: tx = ctx->work[0];
362: ty = ctx->work[1];
363: VecGetArrayRead(x,&array1);
364: VecPlaceArray(tx,array1);
365: VecGetArray(y,&array2);
366: VecPlaceArray(ty,array2);
367: PCApply(ctx->pc,tx,ty);
368: if (n) {
369: PetscBLASIntCast(ld,&ld_);
370: PetscBLASIntCast(n,&n_);
371: for (i=0;i<n;i++) {
372: t[i] = 0.0;
373: for (j=0;j<n;j++) t[i] += ctx->M[i+j*ld]*x2[j];
374: }
375: if (pjd->midx==1) {
376: PetscBLASIntCast(ctx->pep->ncv,&ncv_);
377: for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] -= ctx->theta;
378: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,pjd->T,&ncv_,t,&one,&zero,x2,&one));
379: for (i=0;i<n;i++) pjd->T[i*(1+ctx->pep->ncv)] += ctx->theta;
380: for (i=0;i<n;i++) array2[nloc+i] = x2[i];
381: for (i=0;i<n;i++) x2[i] = -t[i];
382: } else {
383: for (i=0;i<n;i++) array2[nloc+i] = t[i];
384: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,ps,&ld_,t,&one,&zero,x2,&one));
385: }
386: for (i=0;i<n;i++) array2[nloc+i] /= PetscSqrtReal(np);
387: BVSetActiveColumns(pjd->X,0,n);
388: BVMultVec(pjd->X,-1.0,1.0,ty,x2);
389: PetscFree2(x2,t);
390: }
391: VecResetArray(tx);
392: VecResetArray(ty);
393: VecRestoreArrayRead(x,&array1);
394: VecRestoreArray(y,&array2);
395: return(0);
396: }
398: /*
399: Application of shell preconditioner:
400: y = B\x - eta*B\p, with eta = (u'*B\x)/(u'*B\p)
401: */
402: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
403: {
405: PetscScalar rr,eta;
406: PEP_JD_PCSHELL *ctx;
407: PetscInt sz;
408: const Vec *xs,*ys;
409: #if !defined(PETSC_USE_COMPLEX)
410: PetscScalar rx,xr,xx;
411: #endif
414: PCShellGetContext(pc,(void**)&ctx);
415: VecCompGetSubVecs(x,&sz,&xs);
416: VecCompGetSubVecs(y,NULL,&ys);
417: /* y = B\x apply extended PC */
418: PEPJDExtendedPCApply(pc,xs[0],ys[0]);
419: #if !defined(PETSC_USE_COMPLEX)
420: if (sz==2) {
421: PEPJDExtendedPCApply(pc,xs[1],ys[1]);
422: }
423: #endif
425: /* Compute eta = u'*y / u'*Bp */
426: VecDot(ys[0],ctx->u[0],&rr);
427: eta = -rr*ctx->gamma[0];
428: #if !defined(PETSC_USE_COMPLEX)
429: if (sz==2) {
430: VecDot(ys[0],ctx->u[1],&xr);
431: VecDot(ys[1],ctx->u[0],&rx);
432: VecDot(ys[1],ctx->u[1],&xx);
433: eta += -ctx->gamma[0]*xx-ctx->gamma[1]*(-xr+rx);
434: }
435: #endif
436: eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];
438: /* y = y - eta*Bp */
439: VecAXPY(ys[0],eta,ctx->Bp[0]);
440: #if !defined(PETSC_USE_COMPLEX)
441: if (sz==2) {
442: VecAXPY(ys[1],eta,ctx->Bp[1]);
443: eta = -ctx->gamma[1]*(rr+xx)+ctx->gamma[0]*(-xr+rx);
444: eta /= ctx->gamma[0]*ctx->gamma[0]+ctx->gamma[1]*ctx->gamma[1];
445: VecAXPY(ys[0],eta,ctx->Bp[1]);
446: VecAXPY(ys[1],-eta,ctx->Bp[0]);
447: }
448: #endif
449: return(0);
450: }
452: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
453: {
455: PetscMPIInt np,rk,count;
456: PetscScalar *array1,*array2;
457: PetscInt nloc;
460: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
461: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
462: BVGetSizes(pep->V,&nloc,NULL,NULL);
463: if (v) {
464: VecGetArray(v,&array1);
465: VecGetArray(vex,&array2);
466: if (back) {
467: PetscArraycpy(array1,array2,nloc);
468: } else {
469: PetscArraycpy(array2,array1,nloc);
470: }
471: VecRestoreArray(v,&array1);
472: VecRestoreArray(vex,&array2);
473: }
474: if (a) {
475: VecGetArray(vex,&array2);
476: if (back) {
477: PetscArraycpy(a,array2+nloc+off,na);
478: PetscMPIIntCast(na,&count);
479: MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
480: } else {
481: PetscArraycpy(array2+nloc+off,a,na);
482: PetscMPIIntCast(na,&count);
483: MPI_Bcast(array2+nloc+off,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
484: }
485: VecRestoreArray(vex,&array2);
486: }
487: return(0);
488: }
490: /* Computes Phi^hat(lambda) times a vector or its derivative (depends on beval)
491: if no vector is provided returns a matrix
492: */
493: static PetscErrorCode PEPJDEvaluateHatBasis(PEP pep,PetscInt n,PetscScalar *H,PetscInt ldh,PetscScalar *beval,PetscScalar *t,PetscInt idx,PetscScalar *qpp,PetscScalar *qp,PetscScalar *q)
494: {
496: PetscInt j,i;
497: PetscBLASInt n_,ldh_,one=1;
498: PetscReal *a,*b,*g;
499: PetscScalar sone=1.0,zero=0.0;
502: a = pep->pbc; b=a+pep->nmat; g=b+pep->nmat;
503: PetscBLASIntCast(n,&n_);
504: PetscBLASIntCast(ldh,&ldh_);
505: if (idx<1) {
506: PetscArrayzero(q,t?n:n*n);
507: } else if (idx==1) {
508: if (t) {for (j=0;j<n;j++) q[j] = t[j]*beval[idx-1]/a[0];}
509: else {
510: PetscArrayzero(q,n*n);
511: for (j=0;j<n;j++) q[(j+1)*n] = beval[idx-1]/a[0];
512: }
513: } else {
514: if (t) {
515: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n_,&n_,&sone,H,&ldh_,qp,&one,&zero,q,&one));
516: for (j=0;j<n;j++) {
517: q[j] += beval[idx-1]*t[j]-b[idx-1]*qp[j]-g[idx-1]*qpp[j];
518: q[j] /= a[idx-1];
519: }
520: } else {
521: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,H,&ldh_,qp,&n_,&zero,q,&n_));
522: for (j=0;j<n;j++) {
523: q[j+n*j] += beval[idx-1];
524: for (i=0;i<n;i++) {
525: q[i+n*j] += -b[idx-1]*qp[j*n+i]-g[idx-1]*qpp[j*n+i];
526: q[i+n*j] /= a[idx-1];
527: }
528: }
529: }
530: }
531: return(0);
532: }
534: static PetscErrorCode PEPJDComputeResidual(PEP pep,PetscBool derivative,PetscInt sz,Vec *u,PetscScalar *theta,Vec *p,Vec *work)
535: {
536: PEP_JD *pjd = (PEP_JD*)pep->data;
538: PetscMPIInt rk,np,count;
539: Vec tu,tp,w;
540: PetscScalar *dval,*dvali,*array1,*array2,*x2=NULL,*y2,*qj=NULL,*tt=NULL,*xx=NULL,*xxi=NULL,sone=1.0;
541: PetscInt i,j,nconv,nloc;
542: PetscBLASInt n,ld,one=1;
543: #if !defined(PETSC_USE_COMPLEX)
544: Vec tui=NULL,tpi=NULL;
545: PetscScalar *x2i=NULL,*qji=NULL,*qq,*y2i,*arrayi1,*arrayi2;
546: #endif
549: nconv = pjd->nlock;
550: if (!nconv) {
551: PetscMalloc1(2*sz*pep->nmat,&dval);
552: } else {
553: PetscMalloc5(2*pep->nmat,&dval,2*nconv,&xx,nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*pep->nmat,&qj);
554: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
555: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
556: BVGetSizes(pep->V,&nloc,NULL,NULL);
557: VecGetArray(u[0],&array1);
558: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]*PetscSqrtReal(np);
559: VecRestoreArray(u[0],&array1);
560: #if !defined(PETSC_USE_COMPLEX)
561: if (sz==2) {
562: x2i = x2+nconv;
563: VecGetArray(u[1],&arrayi1);
564: for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]*PetscSqrtReal(np);
565: VecRestoreArray(u[1],&arrayi1);
566: }
567: #endif
568: }
569: dvali = dval+pep->nmat;
570: tu = work[0];
571: tp = work[1];
572: w = work[2];
573: VecGetArray(u[0],&array1);
574: VecPlaceArray(tu,array1);
575: VecGetArray(p[0],&array2);
576: VecPlaceArray(tp,array2);
577: VecSet(tp,0.0);
578: #if !defined(PETSC_USE_COMPLEX)
579: if (sz==2) {
580: tui = work[3];
581: tpi = work[4];
582: VecGetArray(u[1],&arrayi1);
583: VecPlaceArray(tui,arrayi1);
584: VecGetArray(p[1],&arrayi2);
585: VecPlaceArray(tpi,arrayi2);
586: VecSet(tpi,0.0);
587: }
588: #endif
589: if (derivative) {
590: PEPEvaluateBasisDerivative(pep,theta[0],theta[1],dval,dvali);
591: } else {
592: PEPEvaluateBasis(pep,theta[0],theta[1],dval,dvali);
593: }
594: for (i=derivative?1:0;i<pep->nmat;i++) {
595: MatMult(pep->A[i],tu,w);
596: VecAXPY(tp,dval[i],w);
597: #if !defined(PETSC_USE_COMPLEX)
598: if (sz==2) {
599: VecAXPY(tpi,dvali[i],w);
600: MatMult(pep->A[i],tui,w);
601: VecAXPY(tpi,dval[i],w);
602: VecAXPY(tp,-dvali[i],w);
603: }
604: #endif
605: }
606: if (nconv) {
607: for (i=0;i<pep->nmat;i++) {
608: 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);
609: }
610: #if !defined(PETSC_USE_COMPLEX)
611: if (sz==2) {
612: qji = qj+nconv*pep->nmat;
613: qq = qji+nconv*pep->nmat;
614: for (i=0;i<pep->nmat;i++) {
615: 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);
616: }
617: for (i=0;i<nconv*pep->nmat;i++) qj[i] -= qji[i];
618: for (i=0;i<pep->nmat;i++) {
619: 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);
620: 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);
621: }
622: for (i=0;i<nconv*pep->nmat;i++) qji[i] += qq[i];
623: for (i=derivative?2:1;i<pep->nmat;i++) {
624: BVMultVec(pjd->AX[i],1.0,1.0,tpi,qji+i*nconv);
625: }
626: }
627: #endif
628: for (i=derivative?2:1;i<pep->nmat;i++) {
629: BVMultVec(pjd->AX[i],1.0,1.0,tp,qj+i*nconv);
630: }
632: /* extended vector part */
633: BVSetActiveColumns(pjd->X,0,nconv);
634: BVDotVec(pjd->X,tu,xx);
635: xxi = xx+nconv;
636: #if !defined(PETSC_USE_COMPLEX)
637: if (sz==2) { BVDotVec(pjd->X,tui,xxi); }
638: #endif
639: if (sz==1) { PetscArrayzero(xxi,nconv); }
640: if (rk==np-1) {
641: PetscBLASIntCast(nconv,&n);
642: PetscBLASIntCast(pjd->ld,&ld);
643: y2 = array2+nloc;
644: PetscArrayzero(y2,nconv);
645: for (j=derivative?1:0;j<pjd->midx;j++) {
646: for (i=0;i<nconv;i++) tt[i] = dval[j]*xx[i]-dvali[j]*xxi[i];
647: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
648: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
649: }
650: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
651: #if !defined(PETSC_USE_COMPLEX)
652: if (sz==2) {
653: y2i = arrayi2+nloc;
654: PetscArrayzero(y2i,nconv);
655: for (j=derivative?1:0;j<pjd->midx;j++) {
656: for (i=0;i<nconv;i++) tt[i] = dval[j]*xxi[i]+dvali[j]*xx[i];
657: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
658: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
659: }
660: for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
661: }
662: #endif
663: }
664: PetscMPIIntCast(nconv,&count);
665: MPI_Bcast(array2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
666: #if !defined(PETSC_USE_COMPLEX)
667: if (sz==2) {
668: MPI_Bcast(arrayi2+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
669: }
670: #endif
671: }
672: if (nconv) {
673: PetscFree5(dval,xx,tt,x2,qj);
674: } else {
675: PetscFree(dval);
676: }
677: VecResetArray(tu);
678: VecRestoreArray(u[0],&array1);
679: VecResetArray(tp);
680: VecRestoreArray(p[0],&array2);
681: #if !defined(PETSC_USE_COMPLEX)
682: if (sz==2) {
683: VecResetArray(tui);
684: VecRestoreArray(u[1],&arrayi1);
685: VecResetArray(tpi);
686: VecRestoreArray(p[1],&arrayi2);
687: }
688: #endif
689: return(0);
690: }
692: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
693: {
694: PEP_JD *pjd = (PEP_JD*)pep->data;
696: PetscScalar *tt,target[2];
697: Vec vg,wg;
698: PetscInt i;
699: PetscReal norm;
702: PetscMalloc1(pjd->ld-1,&tt);
703: if (pep->nini==0) {
704: BVSetRandomColumn(pjd->V,0);
705: for (i=0;i<pjd->ld-1;i++) tt[i] = 0.0;
706: BVGetColumn(pjd->V,0,&vg);
707: PEPJDCopyToExtendedVec(pep,NULL,tt,pjd->ld-1,0,vg,PETSC_FALSE);
708: BVRestoreColumn(pjd->V,0,&vg);
709: BVNormColumn(pjd->V,0,NORM_2,&norm);
710: BVScaleColumn(pjd->V,0,1.0/norm);
711: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
712: BVGetColumn(pjd->V,0,&vg);
713: BVGetColumn(pjd->W,0,&wg);
714: VecSet(wg,0.0);
715: target[0] = pep->target; target[1] = 0.0;
716: PEPJDComputeResidual(pep,PETSC_TRUE,1,&vg,target,&wg,w);
717: BVRestoreColumn(pjd->W,0,&wg);
718: BVRestoreColumn(pjd->V,0,&vg);
719: BVNormColumn(pjd->W,0,NORM_2,&norm);
720: BVScaleColumn(pjd->W,0,1.0/norm);
721: }
722: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for initial vectors not implemented yet");
723: PetscFree(tt);
724: return(0);
725: }
727: static PetscErrorCode PEPJDShellMatMult(Mat P,Vec x,Vec y)
728: {
729: PetscErrorCode ierr;
730: PEP_JD_MATSHELL *matctx;
731: PEP_JD *pjd;
732: PetscInt i,j,nconv,nloc,nmat,ldt,ncv,sz;
733: Vec tx,ty;
734: const Vec *xs,*ys;
735: PetscScalar *array1,*array2,*x2=NULL,*y2,*tt=NULL,*xx=NULL,*xxi,theta[2],sone=1.0,*qj,*val,*vali=NULL;
736: PetscBLASInt n,ld,one=1;
737: PetscMPIInt np;
738: #if !defined(PETSC_USE_COMPLEX)
739: Vec txi=NULL,tyi=NULL;
740: PetscScalar *x2i=NULL,*qji=NULL,*qq,*y2i,*arrayi1,*arrayi2;
741: #endif
744: MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
745: MatShellGetContext(P,(void**)&matctx);
746: pjd = (PEP_JD*)(matctx->pep->data);
747: nconv = pjd->nlock;
748: nmat = matctx->pep->nmat;
749: ncv = matctx->pep->ncv;
750: ldt = pjd->ld;
751: VecCompGetSubVecs(x,&sz,&xs);
752: VecCompGetSubVecs(y,NULL,&ys);
753: theta[0] = matctx->theta[0];
754: theta[1] = (sz==2)?matctx->theta[1]:0.0;
755: if (nconv>0) {
756: PetscMalloc5(nconv,&tt,sz*nconv,&x2,(sz==2?3:1)*nconv*nmat,&qj,2*nconv,&xx,2*nmat,&val);
757: BVGetSizes(matctx->pep->V,&nloc,NULL,NULL);
758: VecGetArray(xs[0],&array1);
759: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]* PetscSqrtReal(np);
760: VecRestoreArray(xs[0],&array1);
761: #if !defined(PETSC_USE_COMPLEX)
762: if (sz==2) {
763: x2i = x2+nconv;
764: VecGetArray(xs[1],&arrayi1);
765: for (i=0;i<nconv;i++) x2i[i] = arrayi1[nloc+i]* PetscSqrtReal(np);
766: VecRestoreArray(xs[1],&arrayi1);
767: }
768: #endif
769: vali = val+nmat;
770: }
771: tx = matctx->work[0];
772: ty = matctx->work[1];
773: VecGetArray(xs[0],&array1);
774: VecPlaceArray(tx,array1);
775: VecGetArray(ys[0],&array2);
776: VecPlaceArray(ty,array2);
777: MatMult(matctx->Pr,tx,ty);
778: #if !defined(PETSC_USE_COMPLEX)
779: if (sz==2) {
780: txi = matctx->work[2];
781: tyi = matctx->work[3];
782: VecGetArray(xs[1],&arrayi1);
783: VecPlaceArray(txi,arrayi1);
784: VecGetArray(ys[1],&arrayi2);
785: VecPlaceArray(tyi,arrayi2);
786: MatMult(matctx->Pr,txi,tyi);
787: if (theta[1]!=0.0) {
788: MatMult(matctx->Pi,txi,matctx->work[4]);
789: VecAXPY(ty,-1.0,matctx->work[4]);
790: MatMult(matctx->Pi,tx,matctx->work[4]);
791: VecAXPY(tyi,1.0,matctx->work[4]);
792: }
793: }
794: #endif
795: if (nconv) {
796: PEPEvaluateBasis(matctx->pep,theta[0],theta[1],val,vali);
797: for (i=0;i<nmat;i++) {
798: 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);
799: }
800: #if !defined(PETSC_USE_COMPLEX)
801: if (sz==2) {
802: qji = qj+nconv*nmat;
803: qq = qji+nconv*nmat;
804: for (i=0;i<nmat;i++) {
805: 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);
806: }
807: for (i=0;i<nconv*nmat;i++) qj[i] -= qji[i];
808: for (i=0;i<nmat;i++) {
809: 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);
810: 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);
811: }
812: for (i=0;i<nconv*nmat;i++) qji[i] += qq[i];
813: for (i=1;i<matctx->pep->nmat;i++) {
814: BVMultVec(pjd->AX[i],1.0,1.0,tyi,qji+i*nconv);
815: }
816: }
817: #endif
818: for (i=1;i<nmat;i++) {
819: BVMultVec(pjd->AX[i],1.0,1.0,ty,qj+i*nconv);
820: }
822: /* extended vector part */
823: BVSetActiveColumns(pjd->X,0,nconv);
824: BVDotVec(pjd->X,tx,xx);
825: xxi = xx+nconv;
826: #if !defined(PETSC_USE_COMPLEX)
827: if (sz==2) { BVDotVec(pjd->X,txi,xxi); }
828: #endif
829: if (sz==1) { PetscArrayzero(xxi,nconv); }
830: PetscBLASIntCast(pjd->nlock,&n);
831: PetscBLASIntCast(ldt,&ld);
832: y2 = array2+nloc;
833: PetscArrayzero(y2,nconv);
834: for (j=0;j<pjd->midx;j++) {
835: for (i=0;i<nconv;i++) tt[i] = val[j]*xx[i]-vali[j]*xxi[i];
836: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qj+j*nconv,&one,&sone,tt,&one));
837: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2,&one));
838: }
839: #if !defined(PETSC_USE_COMPLEX)
840: if (sz==2) {
841: y2i = arrayi2+nloc;
842: PetscArrayzero(y2i,nconv);
843: for (j=0;j<pjd->midx;j++) {
844: for (i=0;i<nconv;i++) tt[i] = val[j]*xxi[i]+vali[j]*xx[i];
845: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,qji+j*nconv,&one,&sone,tt,&one));
846: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&n,&sone,pjd->Tj+j*ld*ld,&ld,tt,&one,&sone,y2i,&one));
847: }
848: for (i=0;i<nconv;i++) arrayi2[nloc+i] /= PetscSqrtReal(np);
849: }
850: #endif
851: for (i=0;i<nconv;i++) array2[nloc+i] /= PetscSqrtReal(np);
852: PetscFree5(tt,x2,qj,xx,val);
853: }
854: VecResetArray(tx);
855: VecRestoreArray(xs[0],&array1);
856: VecResetArray(ty);
857: VecRestoreArray(ys[0],&array2);
858: #if !defined(PETSC_USE_COMPLEX)
859: if (sz==2) {
860: VecResetArray(txi);
861: VecRestoreArray(xs[1],&arrayi1);
862: VecResetArray(tyi);
863: VecRestoreArray(ys[1],&arrayi2);
864: }
865: #endif
866: return(0);
867: }
869: static PetscErrorCode PEPJDSellMatCreateVecs(Mat A,Vec *right,Vec *left)
870: {
871: PetscErrorCode ierr;
872: PEP_JD_MATSHELL *matctx;
873: PEP_JD *pjd;
874: PetscInt kspsf=1,i;
875: Vec v[2];
878: MatShellGetContext(A,(void**)&matctx);
879: pjd = (PEP_JD*)(matctx->pep->data);
880: #if !defined (PETSC_USE_COMPLEX)
881: kspsf = 2;
882: #endif
883: for (i=0;i<kspsf;i++){
884: BVCreateVec(pjd->V,v+i);
885: }
886: if (right) {
887: VecCreateCompWithVecs(v,kspsf,pjd->vtempl,right);
888: }
889: if (left) {
890: VecCreateCompWithVecs(v,kspsf,pjd->vtempl,left);
891: }
892: for (i=0;i<kspsf;i++) {
893: VecDestroy(&v[i]);
894: }
895: return(0);
896: }
898: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
899: {
900: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GETRI) || defined(PETSC_MISSING_LAPACK_GETRF)
902: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GETRI/GETRF - Lapack routines are unavailable");
903: #else
905: PEP_JD *pjd = (PEP_JD*)pep->data;
906: PEP_JD_PCSHELL *pcctx;
907: PetscInt i,j,k,n=pjd->nlock,ld=pjd->ld,deg=pep->nmat-1;
908: PetscScalar *M,*ps,*work,*U,*V,*S,*Sp,*Spp,snone=-1.0,sone=1.0,zero=0.0,*val;
909: PetscReal tol,maxeig=0.0,*sg,*rwork;
910: PetscBLASInt n_,info,ld_,*p,lw_,rk=0;
913: if (n) {
914: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
915: pcctx->theta = theta;
916: pcctx->n = n;
917: M = pcctx->M;
918: PetscBLASIntCast(n,&n_);
919: PetscBLASIntCast(ld,&ld_);
920: if (pjd->midx==1) {
921: PetscArraycpy(M,pjd->XpX,ld*ld);
922: PetscCalloc2(10*n,&work,n,&p);
923: } else {
924: ps = pcctx->ps;
925: PetscCalloc7(2*n*n,&U,3*n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p,deg+1,&val);
926: V = U+n*n;
927: /* pseudo-inverse */
928: for (j=0;j<n;j++) {
929: for (i=0;i<n;i++) S[n*j+i] = -pjd->T[pep->ncv*j+i];
930: S[n*j+j] += theta;
931: }
932: lw_ = 10*n_;
933: #if !defined (PETSC_USE_COMPLEX)
934: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
935: #else
936: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
937: #endif
938: SlepcCheckLapackInfo("gesvd",info);
939: for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
940: tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
941: for (j=0;j<n;j++) {
942: if (sg[j]>tol) {
943: for (i=0;i<n;i++) U[j*n+i] /= sg[j];
944: rk++;
945: } else break;
946: }
947: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));
949: /* compute M */
950: PEPEvaluateBasis(pep,theta,0.0,val,NULL);
951: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&snone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
952: PetscArrayzero(S,2*n*n);
953: Sp = S+n*n;
954: for (j=0;j<n;j++) S[j*(n+1)] = 1.0;
955: for (k=1;k<pjd->midx;k++) {
956: 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];
957: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
958: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
959: Spp = Sp; Sp = S;
960: PEPJDEvaluateHatBasis(pep,n,pjd->T,pep->ncv,val,NULL,k+1,Spp,Sp,S);
961: }
962: }
963: /* inverse */
964: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
965: SlepcCheckLapackInfo("getrf",info);
966: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
967: SlepcCheckLapackInfo("getri",info);
968: if (pjd->midx==1) {
969: PetscFree2(work,p);
970: } else {
971: PetscFree7(U,S,sg,work,rwork,p,val);
972: }
973: }
974: return(0);
975: #endif
976: }
978: static PetscErrorCode PEPJDMatSetUp(PEP pep,PetscInt sz,PetscScalar *theta)
979: {
980: PetscErrorCode ierr;
981: PEP_JD *pjd = (PEP_JD*)pep->data;
982: PEP_JD_MATSHELL *matctx;
983: PEP_JD_PCSHELL *pcctx;
984: MatStructure str;
985: PetscScalar *vals,*valsi;
986: PetscBool skipmat=PETSC_FALSE;
987: PetscInt i;
988: Mat Pr=NULL;
991: if (sz==2 && theta[1]==0.0) sz = 1;
992: MatShellGetContext(pjd->Pshell,(void**)&matctx);
993: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
994: if (matctx->Pr && matctx->theta[0]==theta[0] && ((!matctx->Pi && sz==1) || (sz==2 && matctx->theta[1]==theta[1]))) {
995: if (pcctx->n == pjd->nlock) return(0);
996: skipmat = PETSC_TRUE;
997: }
998: if (!skipmat) {
999: PetscMalloc2(pep->nmat,&vals,pep->nmat,&valsi);
1000: STGetMatStructure(pep->st,&str);
1001: PEPEvaluateBasis(pep,theta[0],theta[1],vals,valsi);
1002: if (!matctx->Pr) {
1003: MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pr);
1004: } else {
1005: MatCopy(pep->A[0],matctx->Pr,str);
1006: }
1007: for (i=1;i<pep->nmat;i++) {
1008: MatAXPY(matctx->Pr,vals[i],pep->A[i],str);
1009: }
1010: if (!pjd->reusepc) {
1011: if (pcctx->PPr && sz==2) {
1012: MatCopy(matctx->Pr,pcctx->PPr,str);
1013: Pr = pcctx->PPr;
1014: } else Pr = matctx->Pr;
1015: }
1016: matctx->theta[0] = theta[0];
1017: #if !defined(PETSC_USE_COMPLEX)
1018: if (sz==2) {
1019: if (!matctx->Pi) {
1020: MatDuplicate(pep->A[0],MAT_COPY_VALUES,&matctx->Pi);
1021: } else {
1022: MatCopy(pep->A[1],matctx->Pi,str);
1023: }
1024: MatScale(matctx->Pi,valsi[1]);
1025: for (i=2;i<pep->nmat;i++) {
1026: MatAXPY(matctx->Pi,valsi[i],pep->A[i],str);
1027: }
1028: matctx->theta[1] = theta[1];
1029: }
1030: #endif
1031: PetscFree2(vals,valsi);
1032: }
1033: if (!pjd->reusepc) {
1034: if (!skipmat) {
1035: PCSetOperators(pcctx->pc,Pr,Pr);
1036: PCSetUp(pcctx->pc);
1037: }
1038: PEPJDUpdateExtendedPC(pep,theta[0]);
1039: }
1040: return(0);
1041: }
1043: static PetscErrorCode PEPJDCreateShellPC(PEP pep,Vec *ww)
1044: {
1045: PEP_JD *pjd = (PEP_JD*)pep->data;
1046: PEP_JD_PCSHELL *pcctx;
1047: PEP_JD_MATSHELL *matctx;
1048: KSP ksp;
1049: PetscInt nloc,mloc,kspsf=1;
1050: Vec v[2];
1051: PetscScalar target[2];
1052: PetscErrorCode ierr;
1053: Mat Pr;
1056: /* Create the reference vector */
1057: BVGetColumn(pjd->V,0,&v[0]);
1058: v[1] = v[0];
1059: #if !defined (PETSC_USE_COMPLEX)
1060: kspsf = 2;
1061: #endif
1062: VecCreateCompWithVecs(v,kspsf,NULL,&pjd->vtempl);
1063: BVRestoreColumn(pjd->V,0,&v[0]);
1064: PetscLogObjectParent((PetscObject)pep,(PetscObject)pjd->vtempl);
1066: /* Replace preconditioner with one containing projectors */
1067: PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
1068: PCSetType(pjd->pcshell,PCSHELL);
1069: PCShellSetName(pjd->pcshell,"PCPEPJD");
1070: PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
1071: PetscNew(&pcctx);
1072: PCShellSetContext(pjd->pcshell,pcctx);
1073: STGetKSP(pep->st,&ksp);
1074: BVCreateVec(pjd->V,&pcctx->Bp[0]);
1075: VecDuplicate(pcctx->Bp[0],&pcctx->Bp[1]);
1076: KSPGetPC(ksp,&pcctx->pc);
1077: PetscObjectReference((PetscObject)pcctx->pc);
1078: MatGetLocalSize(pep->A[0],&mloc,&nloc);
1079: if (pjd->ld>1) {
1080: nloc += pjd->ld-1; mloc += pjd->ld-1;
1081: }
1082: PetscNew(&matctx);
1083: MatCreateShell(PetscObjectComm((PetscObject)pep),kspsf*nloc,kspsf*mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
1084: MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)(void))PEPJDShellMatMult);
1085: MatShellSetOperation(pjd->Pshell,MATOP_CREATE_VECS,(void(*)(void))PEPJDSellMatCreateVecs);
1086: matctx->pep = pep;
1087: target[0] = pep->target; target[1] = 0.0;
1088: PEPJDMatSetUp(pep,1,target);
1089: Pr = matctx->Pr;
1090: pcctx->PPr = NULL;
1091: #if !defined(PETSC_USE_COMPLEX)
1092: if (!pjd->reusepc) {
1093: MatDuplicate(matctx->Pr,MAT_COPY_VALUES,&pcctx->PPr);
1094: Pr = pcctx->PPr;
1095: }
1096: #endif
1097: PCSetOperators(pcctx->pc,Pr,Pr);
1098: PCSetErrorIfFailure(pcctx->pc,PETSC_TRUE);
1099: KSPSetPC(ksp,pjd->pcshell);
1100: if (pjd->reusepc) {
1101: PCSetReusePreconditioner(pcctx->pc,PETSC_TRUE);
1102: KSPSetReusePreconditioner(ksp,PETSC_TRUE);
1103: }
1104: KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
1105: KSPSetUp(ksp);
1106: if (pjd->ld>1) {
1107: PetscMalloc2(pjd->ld*pjd->ld,&pcctx->M,pjd->ld*pjd->ld,&pcctx->ps);
1108: pcctx->pep = pep;
1109: }
1110: matctx->work = ww;
1111: pcctx->work = ww;
1112: return(0);
1113: }
1115: static PetscErrorCode PEPJDEigenvectors(PEP pep)
1116: {
1117: #if defined(SLEPC_MISSING_LAPACK_TREVC)
1119: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
1120: #else
1122: PEP_JD *pjd = (PEP_JD*)pep->data;
1123: PetscBLASInt ld,nconv,info,nc;
1124: PetscScalar *Z,*w;
1125: PetscReal *wr,norm;
1126: PetscInt i;
1127: Mat U;
1128: #if !defined(PETSC_USE_COMPLEX)
1129: Vec v,v1;
1130: #endif
1133: PetscMalloc3(pep->nconv*pep->nconv,&Z,3*pep->ncv,&wr,2*pep->ncv,&w);
1134: PetscBLASIntCast(pep->ncv,&ld);
1135: PetscBLASIntCast(pep->nconv,&nconv);
1136: #if !defined(PETSC_USE_COMPLEX)
1137: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
1138: #else
1139: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
1140: #endif
1141: SlepcCheckLapackInfo("trevc",info);
1142: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
1143: BVSetActiveColumns(pjd->X,0,pep->nconv);
1144: BVMultInPlace(pjd->X,U,0,pep->nconv);
1145: for (i=0;i<pep->nconv;i++) {
1146: #if !defined(PETSC_USE_COMPLEX)
1147: if (pep->eigi[i]!=0.0) { /* first eigenvalue of a complex conjugate pair */
1148: BVGetColumn(pjd->X,i,&v);
1149: BVGetColumn(pjd->X,i+1,&v1);
1150: VecNormalizeComplex(v,v1,PETSC_TRUE,NULL);
1151: BVRestoreColumn(pjd->X,i,&v);
1152: BVRestoreColumn(pjd->X,i+1,&v1);
1153: i++;
1154: } else /* real eigenvalue */
1155: #endif
1156: {
1157: BVNormColumn(pjd->X,i,NORM_2,&norm);
1158: BVScaleColumn(pjd->X,i,1.0/norm);
1159: }
1160: }
1161: MatDestroy(&U);
1162: PetscFree3(Z,wr,w);
1163: return(0);
1164: #endif
1165: }
1167: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv,PetscInt sz)
1168: {
1170: PEP_JD *pjd = (PEP_JD*)pep->data;
1171: PetscInt j,i,*P,ldds,rk=0,nvv=*nv;
1172: Vec v,x,w;
1173: PetscScalar *R,*r,*pX,target[2];
1174: Mat X;
1175: PetscBLASInt sz_,rk_,nv_,info;
1176: PetscMPIInt np;
1179: /* update AX and XpX */
1180: for (i=sz;i>0;i--) {
1181: BVGetColumn(pjd->X,pjd->nlock-i,&x);
1182: for (j=0;j<pep->nmat;j++) {
1183: BVGetColumn(pjd->AX[j],pjd->nlock-i,&v);
1184: MatMult(pep->A[j],x,v);
1185: BVRestoreColumn(pjd->AX[j],pjd->nlock-i,&v);
1186: BVSetActiveColumns(pjd->AX[j],0,pjd->nlock-i+1);
1187: }
1188: BVRestoreColumn(pjd->X,pjd->nlock-i,&x);
1189: BVDotColumn(pjd->X,(pjd->nlock-i),pjd->XpX+(pjd->nlock-i)*(pjd->ld));
1190: pjd->XpX[(pjd->nlock-i)*(1+pjd->ld)] = 1.0;
1191: 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]);
1192: }
1194: /* minimality index */
1195: pjd->midx = PetscMin(pjd->mmidx,pjd->nlock);
1197: /* evaluate the polynomial basis in T */
1198: PetscArrayzero(pjd->Tj,pjd->ld*pjd->ld*pep->nmat);
1199: for (j=0;j<pep->nmat;j++) {
1200: 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);
1201: }
1203: /* Extend search space */
1204: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1205: PetscCalloc3(nvv,&P,nvv*nvv,&R,nvv*sz,&r);
1206: DSGetLeadingDimension(pep->ds,&ldds);
1207: DSGetArray(pep->ds,DS_MAT_X,&pX);
1208: PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
1209: for (j=0;j<sz;j++) {
1210: 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 */
1211: }
1212: PetscBLASIntCast(rk,&rk_);
1213: PetscBLASIntCast(sz,&sz_);
1214: PetscBLASIntCast(nvv,&nv_);
1215: PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
1216: SlepcCheckLapackInfo("trtri",info);
1217: for (i=0;i<sz;i++) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r+i,&sz_));
1218: for (i=0;i<sz*rk;i++) r[i] = PetscConj(r[i])/PetscSqrtReal(np); /* revert */
1219: BVSetActiveColumns(pjd->V,0,nvv);
1220: rk -= sz;
1221: for (j=0;j<rk;j++) {
1222: PetscArraycpy(R+j*nvv,pX+(j+sz)*ldds,nvv);
1223: }
1224: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1225: MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk,R,&X);
1226: BVMultInPlace(pjd->V,X,0,rk);
1227: MatDestroy(&X);
1228: BVSetActiveColumns(pjd->V,0,rk);
1229: for (j=0;j<rk;j++) {
1230: BVGetColumn(pjd->V,j,&v);
1231: PEPJDCopyToExtendedVec(pep,NULL,r+sz*(j+sz),sz,pjd->nlock-sz,v,PETSC_FALSE);
1232: BVRestoreColumn(pjd->V,j,&v);
1233: }
1234: BVOrthogonalize(pjd->V,NULL);
1236: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1237: for (j=0;j<rk;j++) {
1238: /* W = P(target)*V */
1239: BVGetColumn(pjd->W,j,&w);
1240: BVGetColumn(pjd->V,j,&v);
1241: target[0] = pep->target; target[1] = 0.0;
1242: PEPJDComputeResidual(pep,PETSC_FALSE,1,&v,target,&w,pep->work);
1243: BVRestoreColumn(pjd->V,j,&v);
1244: BVRestoreColumn(pjd->W,j,&w);
1245: }
1246: BVSetActiveColumns(pjd->W,0,rk);
1247: BVOrthogonalize(pjd->W,NULL);
1248: }
1249: *nv = rk;
1250: PetscFree3(P,R,r);
1251: return(0);
1252: }
1254: PetscErrorCode PEPJDSystemSetUp(PEP pep,PetscInt sz,PetscScalar *theta,Vec *u,Vec *p,Vec *ww)
1255: {
1257: PEP_JD *pjd = (PEP_JD*)pep->data;
1258: PEP_JD_PCSHELL *pcctx;
1259: #if !defined(PETSC_USE_COMPLEX)
1260: PetscScalar s[2];
1261: #endif
1264: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1265: PEPJDMatSetUp(pep,sz,theta);
1266: pcctx->u[0] = u[0]; pcctx->u[1] = u[1];
1267: /* Compute r'. p is a work space vector */
1268: PEPJDComputeResidual(pep,PETSC_TRUE,sz,u,theta,p,ww);
1269: PEPJDExtendedPCApply(pjd->pcshell,p[0],pcctx->Bp[0]);
1270: VecDot(pcctx->Bp[0],u[0],pcctx->gamma);
1271: #if !defined(PETSC_USE_COMPLEX)
1272: if (sz==2) {
1273: PEPJDExtendedPCApply(pjd->pcshell,p[1],pcctx->Bp[1]);
1274: VecDot(pcctx->Bp[0],u[1],pcctx->gamma+1);
1275: VecMDot(pcctx->Bp[1],2,u,s);
1276: pcctx->gamma[0] += s[1];
1277: pcctx->gamma[1] = -pcctx->gamma[1]+s[0];
1278: }
1279: #endif
1280: if (sz==1) {
1281: VecZeroEntries(pcctx->Bp[1]);
1282: pcctx->gamma[1] = 0.0;
1283: }
1284: return(0);
1285: }
1287: PetscErrorCode PEPSolve_JD(PEP pep)
1288: {
1289: PetscErrorCode ierr;
1290: PEP_JD *pjd = (PEP_JD*)pep->data;
1291: PetscInt k,nv,nvc,ld,minv,dim,bupdated=0,sz=1,kspsf=1,idx,off,maxits,nloc;
1292: PetscMPIInt np,count;
1293: PetscScalar theta[2]={0.0,0.0},ritz[2]={0.0,0.0},*pX,*eig,*eigi,*array;
1294: PetscReal norm,*res,tol=0.0,rtol,abstol, dtol;
1295: PetscBool lindep,ini=PETSC_TRUE;
1296: Vec tc,t[2]={NULL,NULL},u[2]={NULL,NULL},p[2]={NULL,NULL};
1297: Vec rc,rr[2],r[2]={NULL,NULL},*ww=pep->work,v[2];
1298: Mat G,X,Y;
1299: KSP ksp;
1300: PEP_JD_PCSHELL *pcctx;
1301: PEP_JD_MATSHELL *matctx;
1302: #if !defined(PETSC_USE_COMPLEX)
1303: PetscReal norm1;
1304: #endif
1307: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
1308: BVGetSizes(pep->V,&nloc,NULL,NULL);
1309: DSGetLeadingDimension(pep->ds,&ld);
1310: PetscCalloc3(pep->ncv+pep->nev,&eig,pep->ncv+pep->nev,&eigi,pep->ncv+pep->nev,&res);
1311: pjd->nlock = 0;
1312: STGetKSP(pep->st,&ksp);
1313: KSPGetTolerances(ksp,&rtol,&abstol,&dtol,&maxits);
1314: #if !defined (PETSC_USE_COMPLEX)
1315: kspsf = 2;
1316: #endif
1317: PEPJDProcessInitialSpace(pep,ww);
1318: nv = (pep->nini)?pep->nini:1;
1320: /* Replace preconditioner with one containing projectors */
1321: PEPJDCreateShellPC(pep,ww);
1322: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
1324: /* Create auxiliar vectors */
1325: BVCreateVec(pjd->V,&u[0]);
1326: VecDuplicate(u[0],&p[0]);
1327: VecDuplicate(u[0],&r[0]);
1328: #if !defined (PETSC_USE_COMPLEX)
1329: VecDuplicate(u[0],&u[1]);
1330: VecDuplicate(u[0],&p[1]);
1331: VecDuplicate(u[0],&r[1]);
1332: #endif
1334: /* Restart loop */
1335: while (pep->reason == PEP_CONVERGED_ITERATING) {
1336: pep->its++;
1337: DSSetDimensions(pep->ds,nv,0,0,0);
1338: BVSetActiveColumns(pjd->V,bupdated,nv);
1339: PEPJDUpdateTV(pep,bupdated,nv,ww);
1340: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVSetActiveColumns(pjd->W,bupdated,nv); }
1341: for (k=0;k<pep->nmat;k++) {
1342: BVSetActiveColumns(pjd->TV[k],bupdated,nv);
1343: DSGetMat(pep->ds,DSMatExtra[k],&G);
1344: BVMatProject(pjd->TV[k],NULL,pjd->W,G);
1345: DSRestoreMat(pep->ds,DSMatExtra[k],&G);
1346: }
1347: BVSetActiveColumns(pjd->V,0,nv);
1348: BVSetActiveColumns(pjd->W,0,nv);
1350: /* Solve projected problem */
1351: DSSetState(pep->ds,DS_STATE_RAW);
1352: DSSolve(pep->ds,pep->eigr,pep->eigi);
1353: DSSort(pep->ds,pep->eigr,pep->eigi,NULL,NULL,NULL);
1354: DSSynchronize(pep->ds,pep->eigr,pep->eigi);
1355: idx = 0;
1356: do {
1357: ritz[0] = pep->eigr[idx];
1358: #if !defined(PETSC_USE_COMPLEX)
1359: ritz[1] = pep->eigi[idx];
1360: sz = (ritz[1]==0.0)?1:2;
1361: #endif
1362: /* Compute Ritz vector u=V*X(:,1) */
1363: DSGetArray(pep->ds,DS_MAT_X,&pX);
1364: BVSetActiveColumns(pjd->V,0,nv);
1365: BVMultVec(pjd->V,1.0,0.0,u[0],pX+idx*ld);
1366: #if !defined(PETSC_USE_COMPLEX)
1367: if (sz==2) {
1368: BVMultVec(pjd->V,1.0,0.0,u[1],pX+(idx+1)*ld);
1369: }
1370: #endif
1371: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1372: PEPJDComputeResidual(pep,PETSC_FALSE,sz,u,ritz,r,ww);
1373: /* Check convergence */
1374: VecNorm(r[0],NORM_2,&norm);
1375: #if !defined(PETSC_USE_COMPLEX)
1376: if (sz==2) {
1377: VecNorm(r[1],NORM_2,&norm1);
1378: norm = SlepcAbs(norm,norm1);
1379: }
1380: #endif
1381: (*pep->converged)(pep,ritz[0],ritz[1],norm,&pep->errest[pep->nconv],pep->convergedctx);
1382: if (sz==2) pep->errest[pep->nconv+1] = pep->errest[pep->nconv];
1383: if (ini) {
1384: tol = PetscMin(.1,pep->errest[pep->nconv]); ini = PETSC_FALSE;
1385: } else tol = PetscMin(pep->errest[pep->nconv],tol/2);
1386: (*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);
1387: if (pep->errest[pep->nconv]<pep->tol) {
1388: /* Ritz pair converged */
1389: ini = PETSC_TRUE;
1390: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1391: if (pjd->ld>1) {
1392: BVGetColumn(pjd->X,pep->nconv,&v[0]);
1393: PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*pep->nconv,pjd->ld-1,0,u[0],PETSC_TRUE);
1394: BVRestoreColumn(pjd->X,pep->nconv,&v[0]);
1395: BVSetActiveColumns(pjd->X,0,pep->nconv+1);
1396: BVNormColumn(pjd->X,pep->nconv,NORM_2,&norm);
1397: BVScaleColumn(pjd->X,pep->nconv,1.0/norm);
1398: for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*pep->nconv+k] *= PetscSqrtReal(np)/norm;
1399: pjd->T[(pep->ncv+1)*pep->nconv] = ritz[0];
1400: eig[pep->nconv] = ritz[0];
1401: idx++;
1402: #if !defined(PETSC_USE_COMPLEX)
1403: if (sz==2) {
1404: BVGetColumn(pjd->X,pep->nconv+1,&v[0]);
1405: PEPJDCopyToExtendedVec(pep,v[0],pjd->T+pep->ncv*(pep->nconv+1),pjd->ld-1,0,u[1],PETSC_TRUE);
1406: BVRestoreColumn(pjd->X,pep->nconv+1,&v[0]);
1407: BVSetActiveColumns(pjd->X,0,pep->nconv+2);
1408: BVNormColumn(pjd->X,pep->nconv+1,NORM_2,&norm1);
1409: BVScaleColumn(pjd->X,pep->nconv+1,1.0/norm1);
1410: for (k=0;k<pep->nconv;k++) pjd->T[pep->ncv*(pep->nconv+1)+k] *= PetscSqrtReal(np)/norm1;
1411: pjd->T[(pep->ncv+1)*(pep->nconv+1)] = ritz[0];
1412: pjd->T[(pep->ncv+1)*pep->nconv+1] = -ritz[1]*norm1/norm;
1413: pjd->T[(pep->ncv+1)*(pep->nconv+1)-1] = ritz[1]*norm/norm1;
1414: eig[pep->nconv+1] = ritz[0];
1415: eigi[pep->nconv] = ritz[1]; eigi[pep->nconv+1] = -ritz[1];
1416: idx++;
1417: }
1418: #endif
1419: } else {
1420: BVInsertVec(pep->V,pep->nconv,u[0]);
1421: }
1422: pep->nconv += sz;
1423: }
1424: } while (pep->errest[pep->nconv]<pep->tol && pep->nconv<nv);
1426: if (pep->reason==PEP_CONVERGED_ITERATING) {
1427: nvc = nv;
1428: if (idx) {
1429: pjd->nlock +=idx;
1430: PEPJDLockConverged(pep,&nv,idx);
1431: }
1432: if (nv+sz>=pep->ncv-1) {
1433: /* Basis full, force restart */
1434: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
1435: DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
1436: DSGetArray(pep->ds,DS_MAT_X,&pX);
1437: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1438: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
1439: DSGetArray(pep->ds,DS_MAT_Y,&pX);
1440: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
1441: DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
1442: DSGetMat(pep->ds,DS_MAT_X,&X);
1443: BVMultInPlace(pjd->V,X,0,minv);
1444: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1445: DSGetMat(pep->ds,DS_MAT_Y,&Y);
1446: BVMultInPlace(pjd->W,Y,pep->nconv,minv);
1447: DSRestoreMat(pep->ds,DS_MAT_Y,&Y);
1448: }
1449: MatDestroy(&X);
1450: nv = minv;
1451: bupdated = 0;
1452: } else {
1453: if (!idx && pep->errest[pep->nconv]<pjd->fix) {theta[0] = ritz[0]; theta[1] = ritz[1];}
1454: else {theta[0] = pep->target; theta[1] = 0.0;}
1455: /* Update system mat */
1456: PEPJDSystemSetUp(pep,sz,theta,u,p,ww);
1457: /* Solve correction equation to expand basis */
1458: BVGetColumn(pjd->V,nv,&t[0]);
1459: rr[0] = r[0];
1460: if (sz==2) {
1461: BVGetColumn(pjd->V,nv+1,&t[1]);
1462: rr[1] = r[1];
1463: } else {
1464: t[1] = NULL;
1465: rr[1] = NULL;
1466: }
1467: VecCreateCompWithVecs(t,kspsf,pjd->vtempl,&tc);
1468: VecCreateCompWithVecs(rr,kspsf,pjd->vtempl,&rc);
1469: VecCompSetSubVecs(pjd->vtempl,sz,NULL);
1470: tol = PetscMax(rtol,tol/2);
1471: KSPSetTolerances(ksp,tol,abstol,dtol,maxits);
1472: KSPSolve(ksp,rc,tc);
1473: VecDestroy(&tc);
1474: VecDestroy(&rc);
1475: VecGetArray(t[0],&array);
1476: PetscMPIIntCast(pep->nconv,&count);
1477: MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1478: VecRestoreArray(t[0],&array);
1479: BVRestoreColumn(pjd->V,nv,&t[0]);
1480: BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
1481: if (lindep || norm==0.0) {
1482: if (sz==1) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1483: else off = 1;
1484: } else {
1485: off = 0;
1486: BVScaleColumn(pjd->V,nv,1.0/norm);
1487: }
1488: #if !defined(PETSC_USE_COMPLEX)
1489: if (sz==2) {
1490: VecGetArray(t[1],&array);
1491: MPI_Bcast(array+nloc,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
1492: VecRestoreArray(t[1],&array);
1493: BVRestoreColumn(pjd->V,nv+1,&t[1]);
1494: if (off) {
1495: BVCopyColumn(pjd->V,nv+1,nv);
1496: }
1497: BVOrthogonalizeColumn(pjd->V,nv+1-off,NULL,&norm,&lindep);
1498: if (lindep || norm==0.0) {
1499: if (off) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1500: else off = 1;
1501: } else {
1502: BVScaleColumn(pjd->V,nv+1-off,1.0/norm);
1503: }
1504: }
1505: #endif
1506: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) {
1507: BVInsertVec(pjd->W,nv,r[0]);
1508: if (sz==2 && !off) {
1509: BVInsertVec(pjd->W,nv+1,r[1]);
1510: }
1511: BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
1512: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1513: BVScaleColumn(pjd->W,nv,1.0/norm);
1514: if (sz==2 && !off) {
1515: BVOrthogonalizeColumn(pjd->W,nv+1,NULL,&norm,&lindep);
1516: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1517: BVScaleColumn(pjd->W,nv+1,1.0/norm);
1518: }
1519: }
1520: bupdated = idx?0:nv;
1521: nv += sz-off;
1522: }
1523: for (k=0;k<nvc;k++) {
1524: eig[pep->nconv-idx+k] = pep->eigr[k];
1525: #if !defined(PETSC_USE_COMPLEX)
1526: eigi[pep->nconv-idx+k] = pep->eigi[k];
1527: #endif
1528: }
1529: PEPMonitor(pep,pep->its,pep->nconv,eig,eigi,pep->errest,pep->nconv+1);
1530: }
1531: }
1532: if (pjd->ld>1) {
1533: for (k=0;k<pep->nconv;k++) {
1534: pep->eigr[k] = eig[k];
1535: pep->eigi[k] = eigi[k];
1536: }
1537: if (pep->nconv>0) { PEPJDEigenvectors(pep); }
1538: PetscFree2(pcctx->M,pcctx->ps);
1539: }
1540: VecDestroy(&u[0]);
1541: VecDestroy(&r[0]);
1542: VecDestroy(&p[0]);
1543: #if !defined (PETSC_USE_COMPLEX)
1544: VecDestroy(&u[1]);
1545: VecDestroy(&r[1]);
1546: VecDestroy(&p[1]);
1547: #endif
1548: KSPSetTolerances(ksp,rtol,abstol,dtol,maxits);
1549: KSPSetPC(ksp,pcctx->pc);
1550: VecDestroy(&pcctx->Bp[0]);
1551: VecDestroy(&pcctx->Bp[1]);
1552: MatShellGetContext(pjd->Pshell,(void**)&matctx);
1553: MatDestroy(&matctx->Pr);
1554: MatDestroy(&matctx->Pi);
1555: MatDestroy(&pjd->Pshell);
1556: MatDestroy(&pcctx->PPr);
1557: PCDestroy(&pcctx->pc);
1558: PetscFree(pcctx);
1559: PetscFree(matctx);
1560: PCDestroy(&pjd->pcshell);
1561: PetscFree3(eig,eigi,res);
1562: VecDestroy(&pjd->vtempl);
1563: return(0);
1564: }
1566: PetscErrorCode PEPJDSetRestart_JD(PEP pep,PetscReal keep)
1567: {
1568: PEP_JD *pjd = (PEP_JD*)pep->data;
1571: if (keep==PETSC_DEFAULT) pjd->keep = 0.5;
1572: else {
1573: 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]");
1574: pjd->keep = keep;
1575: }
1576: return(0);
1577: }
1579: /*@
1580: PEPJDSetRestart - Sets the restart parameter for the Jacobi-Davidson
1581: method, in particular the proportion of basis vectors that must be kept
1582: after restart.
1584: Logically Collective on pep
1586: Input Parameters:
1587: + pep - the eigenproblem solver context
1588: - keep - the number of vectors to be kept at restart
1590: Options Database Key:
1591: . -pep_jd_restart - Sets the restart parameter
1593: Notes:
1594: Allowed values are in the range [0.1,0.9]. The default is 0.5.
1596: Level: advanced
1598: .seealso: PEPJDGetRestart()
1599: @*/
1600: PetscErrorCode PEPJDSetRestart(PEP pep,PetscReal keep)
1601: {
1607: PetscTryMethod(pep,"PEPJDSetRestart_C",(PEP,PetscReal),(pep,keep));
1608: return(0);
1609: }
1611: PetscErrorCode PEPJDGetRestart_JD(PEP pep,PetscReal *keep)
1612: {
1613: PEP_JD *pjd = (PEP_JD*)pep->data;
1616: *keep = pjd->keep;
1617: return(0);
1618: }
1620: /*@
1621: PEPJDGetRestart - Gets the restart parameter used in the Jacobi-Davidson method.
1623: Not Collective
1625: Input Parameter:
1626: . pep - the eigenproblem solver context
1628: Output Parameter:
1629: . keep - the restart parameter
1631: Level: advanced
1633: .seealso: PEPJDSetRestart()
1634: @*/
1635: PetscErrorCode PEPJDGetRestart(PEP pep,PetscReal *keep)
1636: {
1642: PetscUseMethod(pep,"PEPJDGetRestart_C",(PEP,PetscReal*),(pep,keep));
1643: return(0);
1644: }
1646: PetscErrorCode PEPJDSetFix_JD(PEP pep,PetscReal fix)
1647: {
1648: PEP_JD *pjd = (PEP_JD*)pep->data;
1651: if (fix == PETSC_DEFAULT || fix == PETSC_DECIDE) pjd->fix = 0.01;
1652: else {
1653: if (fix < 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid fix value");
1654: pjd->fix = fix;
1655: }
1656: return(0);
1657: }
1659: /*@
1660: PEPJDSetFix - Sets the threshold for changing the target in the correction
1661: equation.
1663: Logically Collective on pep
1665: Input Parameters:
1666: + pep - the eigenproblem solver context
1667: - fix - threshold for changing the target
1669: Options Database Key:
1670: . -pep_jd_fix - the fix value
1672: Note:
1673: The target in the correction equation is fixed at the first iterations.
1674: When the norm of the residual vector is lower than the fix value,
1675: the target is set to the corresponding eigenvalue.
1677: Level: advanced
1679: .seealso: PEPJDGetFix()
1680: @*/
1681: PetscErrorCode PEPJDSetFix(PEP pep,PetscReal fix)
1682: {
1688: PetscTryMethod(pep,"PEPJDSetFix_C",(PEP,PetscReal),(pep,fix));
1689: return(0);
1690: }
1692: PetscErrorCode PEPJDGetFix_JD(PEP pep,PetscReal *fix)
1693: {
1694: PEP_JD *pjd = (PEP_JD*)pep->data;
1697: *fix = pjd->fix;
1698: return(0);
1699: }
1701: /*@
1702: PEPJDGetFix - Returns the threshold for changing the target in the correction
1703: equation.
1705: Not Collective
1707: Input Parameter:
1708: . pep - the eigenproblem solver context
1710: Output Parameter:
1711: . fix - threshold for changing the target
1713: Note:
1714: The target in the correction equation is fixed at the first iterations.
1715: When the norm of the residual vector is lower than the fix value,
1716: the target is set to the corresponding eigenvalue.
1718: Level: advanced
1720: .seealso: PEPJDSetFix()
1721: @*/
1722: PetscErrorCode PEPJDGetFix(PEP pep,PetscReal *fix)
1723: {
1729: PetscUseMethod(pep,"PEPJDGetFix_C",(PEP,PetscReal*),(pep,fix));
1730: return(0);
1731: }
1733: PetscErrorCode PEPJDSetReusePreconditioner_JD(PEP pep,PetscBool reusepc)
1734: {
1735: PEP_JD *pjd = (PEP_JD*)pep->data;
1738: pjd->reusepc = reusepc;
1739: return(0);
1740: }
1742: /*@
1743: PEPJDSetReusePreconditioner - Sets a flag indicating whether the preconditioner
1744: must be reused or not.
1746: Logically Collective on pep
1748: Input Parameters:
1749: + pep - the eigenproblem solver context
1750: - reusepc - the reuse flag
1752: Options Database Key:
1753: . -pep_jd_reuse_preconditioner - the reuse flag
1755: Note:
1756: The default value is False. If set to True, the preconditioner is built
1757: only at the beginning, using the target value. Otherwise, it may be rebuilt
1758: (depending on the fix parameter) at each iteration from the Ritz value.
1760: Level: advanced
1762: .seealso: PEPJDGetReusePreconditioner(), PEPJDSetFix()
1763: @*/
1764: PetscErrorCode PEPJDSetReusePreconditioner(PEP pep,PetscBool reusepc)
1765: {
1771: PetscTryMethod(pep,"PEPJDSetReusePreconditioner_C",(PEP,PetscBool),(pep,reusepc));
1772: return(0);
1773: }
1775: PetscErrorCode PEPJDGetReusePreconditioner_JD(PEP pep,PetscBool *reusepc)
1776: {
1777: PEP_JD *pjd = (PEP_JD*)pep->data;
1780: *reusepc = pjd->reusepc;
1781: return(0);
1782: }
1784: /*@
1785: PEPJDGetReusePreconditioner - Returns the flag for reusing the preconditioner.
1787: Not Collective
1789: Input Parameter:
1790: . pep - the eigenproblem solver context
1792: Output Parameter:
1793: . reusepc - the reuse flag
1795: Level: advanced
1797: .seealso: PEPJDSetReusePreconditioner()
1798: @*/
1799: PetscErrorCode PEPJDGetReusePreconditioner(PEP pep,PetscBool *reusepc)
1800: {
1806: PetscUseMethod(pep,"PEPJDGetReusePreconditioner_C",(PEP,PetscBool*),(pep,reusepc));
1807: return(0);
1808: }
1810: PetscErrorCode PEPJDSetMinimalityIndex_JD(PEP pep,PetscInt mmidx)
1811: {
1812: PEP_JD *pjd = (PEP_JD*)pep->data;
1815: if (mmidx == PETSC_DEFAULT || mmidx == PETSC_DECIDE) pjd->mmidx = 1;
1816: else {
1817: if (mmidx < 1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid mmidx value");
1818: pjd->mmidx = mmidx;
1819: pep->state = PEP_STATE_INITIAL;
1820: }
1821: return(0);
1822: }
1824: /*@
1825: PEPJDSetMinimalityIndex - Sets the maximum allowed value for the minimality index.
1827: Logically Collective on pep
1829: Input Parameters:
1830: + pep - the eigenproblem solver context
1831: - mmidx - maximum minimality index
1833: Options Database Key:
1834: . -pep_jd_minimality_index - the minimality index value
1836: Note:
1837: The default value is equal to the degree of the polynomial. A smaller value
1838: can be used if the wanted eigenvectors are known to be linearly independent.
1840: Level: advanced
1842: .seealso: PEPJDGetMinimalityIndex()
1843: @*/
1844: PetscErrorCode PEPJDSetMinimalityIndex(PEP pep,PetscInt mmidx)
1845: {
1851: PetscTryMethod(pep,"PEPJDSetMinimalityIndex_C",(PEP,PetscInt),(pep,mmidx));
1852: return(0);
1853: }
1855: PetscErrorCode PEPJDGetMinimalityIndex_JD(PEP pep,PetscInt *mmidx)
1856: {
1857: PEP_JD *pjd = (PEP_JD*)pep->data;
1860: *mmidx = pjd->mmidx;
1861: return(0);
1862: }
1864: /*@
1865: PEPJDGetMinimalityIndex - Returns the maximum allowed value of the minimality
1866: index.
1868: Not Collective
1870: Input Parameter:
1871: . pep - the eigenproblem solver context
1873: Output Parameter:
1874: . mmidx - minimality index
1876: Level: advanced
1878: .seealso: PEPJDSetMinimalityIndex()
1879: @*/
1880: PetscErrorCode PEPJDGetMinimalityIndex(PEP pep,PetscInt *mmidx)
1881: {
1887: PetscUseMethod(pep,"PEPJDGetMinimalityIndex_C",(PEP,PetscInt*),(pep,mmidx));
1888: return(0);
1889: }
1891: PetscErrorCode PEPJDSetProjection_JD(PEP pep,PEPJDProjection proj)
1892: {
1893: PEP_JD *pjd = (PEP_JD*)pep->data;
1896: switch (proj) {
1897: case PEP_JD_PROJECTION_HARMONIC:
1898: case PEP_JD_PROJECTION_ORTHOGONAL:
1899: if (pjd->proj != proj) {
1900: pep->state = PEP_STATE_INITIAL;
1901: pjd->proj = proj;
1902: }
1903: break;
1904: default:
1905: SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'proj' value");
1906: }
1907: return(0);
1908: }
1910: /*@
1911: PEPJDSetProjection - Sets the type of projection to be used in the Jacobi-Davidson solver.
1913: Logically Collective on pep
1915: Input Parameters:
1916: + pep - the eigenproblem solver context
1917: - proj - the type of projection
1919: Options Database Key:
1920: . -pep_jd_projection - the projection type, either orthogonal or harmonic
1922: Level: advanced
1924: .seealso: PEPJDGetProjection()
1925: @*/
1926: PetscErrorCode PEPJDSetProjection(PEP pep,PEPJDProjection proj)
1927: {
1933: PetscTryMethod(pep,"PEPJDSetProjection_C",(PEP,PEPJDProjection),(pep,proj));
1934: return(0);
1935: }
1937: PetscErrorCode PEPJDGetProjection_JD(PEP pep,PEPJDProjection *proj)
1938: {
1939: PEP_JD *pjd = (PEP_JD*)pep->data;
1942: *proj = pjd->proj;
1943: return(0);
1944: }
1946: /*@
1947: PEPJDGetProjection - Returns the type of projection used by the Jacobi-Davidson solver.
1949: Not Collective
1951: Input Parameter:
1952: . pep - the eigenproblem solver context
1954: Output Parameter:
1955: . proj - the type of projection
1957: Level: advanced
1959: .seealso: PEPJDSetProjection()
1960: @*/
1961: PetscErrorCode PEPJDGetProjection(PEP pep,PEPJDProjection *proj)
1962: {
1968: PetscUseMethod(pep,"PEPJDGetProjection_C",(PEP,PEPJDProjection*),(pep,proj));
1969: return(0);
1970: }
1972: PetscErrorCode PEPSetFromOptions_JD(PetscOptionItems *PetscOptionsObject,PEP pep)
1973: {
1974: PetscErrorCode ierr;
1975: PetscBool flg,b1;
1976: PetscReal r1;
1977: PetscInt i1;
1978: PEPJDProjection proj;
1981: PetscOptionsHead(PetscOptionsObject,"PEP JD Options");
1983: PetscOptionsReal("-pep_jd_restart","Proportion of vectors kept after restart","PEPJDSetRestart",0.5,&r1,&flg);
1984: if (flg) { PEPJDSetRestart(pep,r1); }
1986: PetscOptionsReal("-pep_jd_fix","Tolerance for changing the target in the correction equation","PEPJDSetFix",0.01,&r1,&flg);
1987: if (flg) { PEPJDSetFix(pep,r1); }
1989: PetscOptionsBool("-pep_jd_reuse_preconditioner","Whether to reuse the preconditioner","PEPJDSetReusePreconditoiner",PETSC_FALSE,&b1,&flg);
1990: if (flg) { PEPJDSetReusePreconditioner(pep,b1); }
1992: PetscOptionsInt("-pep_jd_minimality_index","Maximum allowed minimality index","PEPJDSetMinimalityIndex",1,&i1,&flg);
1993: if (flg) { PEPJDSetMinimalityIndex(pep,i1); }
1995: PetscOptionsEnum("-pep_jd_projection","Type of projection","PEPJDSetProjection",PEPJDProjectionTypes,(PetscEnum)PEP_JD_PROJECTION_HARMONIC,(PetscEnum*)&proj,&flg);
1996: if (flg) { PEPJDSetProjection(pep,proj); }
1998: PetscOptionsTail();
1999: return(0);
2000: }
2002: PetscErrorCode PEPView_JD(PEP pep,PetscViewer viewer)
2003: {
2005: PEP_JD *pjd = (PEP_JD*)pep->data;
2006: PetscBool isascii;
2009: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
2010: if (isascii) {
2011: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*pjd->keep));
2012: PetscViewerASCIIPrintf(viewer," threshold for changing the target in the correction equation (fix): %g\n",(double)pjd->fix);
2013: PetscViewerASCIIPrintf(viewer," projection type: %s\n",PEPJDProjectionTypes[pjd->proj]);
2014: PetscViewerASCIIPrintf(viewer," maximum allowed minimality index: %d\n",pjd->mmidx);
2015: if (pjd->reusepc) { PetscViewerASCIIPrintf(viewer," reusing the preconditioner\n"); }
2016: }
2017: return(0);
2018: }
2020: PetscErrorCode PEPSetDefaultST_JD(PEP pep)
2021: {
2023: KSP ksp;
2026: if (!((PetscObject)pep->st)->type_name) {
2027: STSetType(pep->st,STPRECOND);
2028: STPrecondSetKSPHasMat(pep->st,PETSC_TRUE);
2029: }
2030: STSetTransform(pep->st,PETSC_FALSE);
2031: STGetKSP(pep->st,&ksp);
2032: if (!((PetscObject)ksp)->type_name) {
2033: KSPSetType(ksp,KSPBCGSL);
2034: KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
2035: }
2036: return(0);
2037: }
2039: PetscErrorCode PEPReset_JD(PEP pep)
2040: {
2042: PEP_JD *pjd = (PEP_JD*)pep->data;
2043: PetscInt i;
2046: for (i=0;i<pep->nmat;i++) {
2047: BVDestroy(pjd->TV+i);
2048: }
2049: if (pjd->proj==PEP_JD_PROJECTION_HARMONIC) { BVDestroy(&pjd->W); }
2050: if (pjd->ld>1) {
2051: BVDestroy(&pjd->V);
2052: for (i=0;i<pep->nmat;i++) {
2053: BVDestroy(pjd->AX+i);
2054: }
2055: BVDestroy(&pjd->N[0]);
2056: BVDestroy(&pjd->N[1]);
2057: PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
2058: }
2059: PetscFree2(pjd->TV,pjd->AX);
2060: return(0);
2061: }
2063: PetscErrorCode PEPDestroy_JD(PEP pep)
2064: {
2068: PetscFree(pep->data);
2069: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
2070: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
2071: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",NULL);
2072: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",NULL);
2073: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",NULL);
2074: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",NULL);
2075: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",NULL);
2076: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",NULL);
2077: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",NULL);
2078: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",NULL);
2079: return(0);
2080: }
2082: SLEPC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
2083: {
2084: PEP_JD *pjd;
2088: PetscNewLog(pep,&pjd);
2089: pep->data = (void*)pjd;
2091: pjd->fix = 0.01;
2092: pjd->mmidx = 0;
2094: pep->ops->solve = PEPSolve_JD;
2095: pep->ops->setup = PEPSetUp_JD;
2096: pep->ops->setfromoptions = PEPSetFromOptions_JD;
2097: pep->ops->destroy = PEPDestroy_JD;
2098: pep->ops->reset = PEPReset_JD;
2099: pep->ops->view = PEPView_JD;
2100: pep->ops->setdefaultst = PEPSetDefaultST_JD;
2102: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
2103: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
2104: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetFix_C",PEPJDSetFix_JD);
2105: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetFix_C",PEPJDGetFix_JD);
2106: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetReusePreconditioner_C",PEPJDSetReusePreconditioner_JD);
2107: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetReusePreconditioner_C",PEPJDGetReusePreconditioner_JD);
2108: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetMinimalityIndex_C",PEPJDSetMinimalityIndex_JD);
2109: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetMinimalityIndex_C",PEPJDGetMinimalityIndex_JD);
2110: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetProjection_C",PEPJDSetProjection_JD);
2111: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetProjection_C",PEPJDGetProjection_JD);
2112: return(0);
2113: }