Actual source code: pjd.c
slepc-3.8.2 2017-12-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, 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: BV V; /* work basis vectors to store the search space */
38: BV W; /* work basis vectors to store the test space */
39: 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) */
40: BV *AX; /* work basis vectors to store A_i*X for locked eigenvectors */
41: BV X; /* locked eigenvectors */
42: PetscScalar *T; /* matrix of the invariant pair */
43: PetscScalar *Tj; /* matrix containing the powers of the invariant pair matrix */
44: PetscScalar *XpX; /* X^H*X */
45: PC pcshell; /* preconditioner including basic precond+projector */
46: Mat Pshell; /* auxiliary shell matrix */
47: PetscInt nconv; /* number of locked vectors in the invariant pair */
48: } PEP_JD;
50: typedef struct {
51: PC pc; /* basic preconditioner */
52: Vec Bp; /* preconditioned residual of derivative polynomial, B\p */
53: Vec u; /* Ritz vector */
54: PetscScalar gamma; /* precomputed scalar u'*B\p */
55: PetscScalar *M;
56: PetscScalar *ps;
57: PetscInt ld;
58: Vec *work;
59: BV X;
60: PetscInt n;
61: } PEP_JD_PCSHELL;
63: typedef struct {
64: Mat P; /* */
65: PEP pep;
66: Vec *work;
67: PetscScalar theta;
68: } PEP_JD_MATSHELL;
70: /*
71: Duplicate and resize auxiliary basis
72: */
73: static PetscErrorCode PEPJDDuplicateBasis(PEP pep,BV *basis)
74: {
75: PetscErrorCode ierr;
76: PetscInt nloc,m;
77: PetscMPIInt rank,nproc;
78: BVType type;
79: BVOrthogType otype;
80: BVOrthogRefineType oref;
81: PetscReal oeta;
82: BVOrthogBlockType oblock;
85: if (pep->nev>1) {
86: BVCreate(PetscObjectComm((PetscObject)pep),basis);
87: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rank);
88: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&nproc);
89: BVGetSizes(pep->V,&nloc,NULL,&m);
90: if (rank==nproc-1) nloc += pep->nev-1;
91: BVSetSizes(*basis,nloc,PETSC_DECIDE,m);
92: BVGetType(pep->V,&type);
93: BVSetType(*basis,type);
94: BVGetOrthogonalization(pep->V,&otype,&oref,&oeta,&oblock);
95: BVSetOrthogonalization(*basis,otype,oref,oeta,oblock);
96: PetscObjectStateIncrease((PetscObject)*basis);
97: } else {
98: BVDuplicate(pep->V,basis);
99: }
100: return(0);
101: }
103: PetscErrorCode PEPSetUp_JD(PEP pep)
104: {
106: PEP_JD *pjd = (PEP_JD*)pep->data;
107: PetscBool isprecond,flg;
108: PetscInt i;
111: pep->lineariz = PETSC_FALSE;
112: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
113: if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
114: if (!pep->which) pep->which = PEP_TARGET_MAGNITUDE;
115: if (pep->which != PEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEPJD only supports which=target_magnitude");;
117: PetscObjectTypeCompare((PetscObject)pep->st,STPRECOND,&isprecond);
118: if (!isprecond) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"JD only works with PRECOND spectral transformation");
120: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
121: STGetTransform(pep->st,&flg);
122: if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");
124: if (!pjd->keep) pjd->keep = 0.5;
126: PEPAllocateSolution(pep,0);
127: PEPSetWorkVecs(pep,5);
128: PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
129: for (i=0;i<pep->nmat;i++) {
130: PEPJDDuplicateBasis(pep,pjd->TV+i);
131: }
132: PEPJDDuplicateBasis(pep,&pjd->W);
133: if (pep->nev>1) {
134: PEPJDDuplicateBasis(pep,&pjd->V);
135: BVSetFromOptions(pjd->V);
136: for (i=0;i<pep->nmat;i++) {
137: BVDuplicateResize(pep->V,pep->nev-1,pjd->AX+i);
138: }
139: BVDuplicateResize(pep->V,pep->nev,&pjd->X);
140: PetscCalloc3((pep->nev)*(pep->nev),&pjd->XpX,pep->nev*pep->nev,&pjd->T,pep->nev*pep->nev*pep->nmat,&pjd->Tj);
141: } else pjd->V = pep->V;
142: DSSetType(pep->ds,DSPEP);
143: DSPEPSetDegree(pep->ds,pep->nmat-1);
144: DSAllocate(pep->ds,pep->ncv);
145: return(0);
146: }
148: /*
149: Updates columns (low to (high-1)) of TV[i]
150: */
151: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
152: {
154: PEP_JD *pjd = (PEP_JD*)pep->data;
155: PetscInt pp,col,i,j,nloc,nconv,deg=pep->nmat-1;
156: Vec v1,v2,t1,t2;
157: PetscScalar *array1,*array2,*x2,*tt,*xx,*y2,zero=0.0,sone=1.0;
158: PetscMPIInt rk,np,count;
159: PetscBLASInt n,ld,one=1;
162: nconv = pjd->nconv;
163: PetscMalloc3(nconv,&tt,nconv,&x2,nconv,&xx);
164: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
165: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
166: BVGetSizes(pep->V,&nloc,NULL,NULL);
167: t1 = w[0];
168: t2 = w[1];
169: for (col=low;col<high;col++) {
170: BVGetColumn(pjd->V,col,&v1);
171: VecGetArray(v1,&array1);
172: if (nconv>0) {
173: if (rk==np-1) { for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]; }
174: PetscMPIIntCast(nconv,&count);
175: MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
176: }
177: VecPlaceArray(t1,array1);
178: for (pp=0;pp<pep->nmat;pp++) {
179: BVGetColumn(pjd->TV[pp],col,&v2);
180: VecGetArray(v2,&array2);
181: VecPlaceArray(t2,array2);
182: MatMult(pep->A[pp],t1,t2);
183: if (nconv) {
184: PetscBLASIntCast(pjd->nconv,&n);
185: PetscBLASIntCast(pep->nev,&ld);
186: for (j=0;j<nconv;j++) tt[j] = x2[j];
187: for (i=pp+1;i<pep->nmat;i++) {
188: BVMultVec(pjd->AX[i],1.0,1.0,t2,tt);
189: if (i!=pep->nmat-1) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,tt,&one));
190: }
191: BVDotVec(pjd->X,t1,xx);
192: if (rk==np-1 && pp<deg) {
193: y2 = array2+nloc;
194: for (j=0;j<nconv;j++) { y2[j] = xx[j]; xx[j] = x2[j]; }
195: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*pp,&ld,y2,&one));
196: for (i=pp+1;i<pep->nmat-1;i++) {
197: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,xx,&one,&zero,tt,&one));
198: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*i,&ld,tt,&one));
199: for (j=0;j<nconv;j++) y2[j] += tt[j];
200: if (i<pep->nmat-2) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,xx,&one));
201: }
202: }
203: }
204: VecResetArray(t2);
205: VecRestoreArray(v2,&array2);
206: BVRestoreColumn(pjd->TV[pp],col,&v2);
207: }
208: VecResetArray(t1);
209: VecRestoreArray(v1,&array1);
210: BVRestoreColumn(pjd->V,col,&v1);
211: }
212: PetscFree3(tt,x2,xx);
213: return(0);
214: }
216: /*
217: RRQR of X. Xin*P=Xou*R. Rank of R is rk
218: */
219: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
220: {
221: #if defined(SLEPC_MISSING_LAPACK_GEQP3) || defined(PETSC_MISSING_LAPACK_ORGQR)
223: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQP3/QRGQR - Lapack routines are unavailable");
224: #else
226: PetscInt i,j,n,r;
227: PetscBLASInt row_,col_,ldx_,*p,lwork,info,n_;
228: PetscScalar *tau,*work;
229: PetscReal tol,*rwork;
232: PetscBLASIntCast(row,&row_);
233: PetscBLASIntCast(col,&col_);
234: PetscBLASIntCast(ldx,&ldx_);
235: n = PetscMin(row,col);
236: PetscBLASIntCast(n,&n_);
237: lwork = 3*col_+1;
238: PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
239: for (i=1;i<col;i++) p[i] = 0;
240: p[0] = 1;
242: /* rank revealing QR */
243: #if defined(PETSC_USE_COMPLEX)
244: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
245: #else
246: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
247: #endif
248: SlepcCheckLapackInfo("geqp3",info);
249: if (P) for (i=0;i<col;i++) P[i] = p[i];
251: /* rank computation */
252: tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
253: r = 1;
254: for (i=1;i<n;i++) {
255: if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
256: else break;
257: }
258: if (rk) *rk=r;
260: /* copy upper triangular matrix if requested */
261: if (R) {
262: for (i=0;i<r;i++) {
263: PetscMemzero(R+i*ldr,r*sizeof(PetscScalar));
264: for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
265: }
266: }
267: PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
268: SlepcCheckLapackInfo("ungqr",info);
269: PetscFree4(p,tau,work,rwork);
270: return(0);
271: #endif
272: }
274: /*
275: Application of extended preconditioner
276: */
277: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
278: {
279: PetscInt i,j,nloc,n,ld;
280: PetscMPIInt rk,np,count;
281: Vec tx,ty;
282: PEP_JD_PCSHELL *ctx;
283: PetscErrorCode ierr;
284: const PetscScalar *array1;
285: PetscScalar *x2=NULL,*t=NULL,*ps,*array2;
286: PetscBLASInt one=1.0,ld_,n_;
289: PCShellGetContext(pc,(void**)&ctx);
290: n = ctx->n;
291: ps = ctx->ps;
292: ld = ctx->ld;
293: if (n) {
294: PetscMalloc2(n,&x2,n,&t);
295: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rk);
296: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
297: if (rk==np-1) {
298: VecGetLocalSize(ctx->work[0],&nloc);
299: VecGetArrayRead(x,&array1);
300: for (i=0;i<n;i++) x2[i] = array1[nloc+i];
301: VecRestoreArrayRead(x,&array1);
302: }
303: PetscMPIIntCast(n,&count);
304: MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pc));
305: }
307: /* y = B\x apply PC */
308: tx = ctx->work[0];
309: ty = ctx->work[1];
310: VecGetArrayRead(x,&array1);
311: VecPlaceArray(tx,array1);
312: VecGetArray(y,&array2);
313: VecPlaceArray(ty,array2);
314: PCApply(ctx->pc,tx,ty);
315: if (n) {
316: for (j=0;j<n;j++) {
317: t[j] = 0.0;
318: for (i=0;i<n;i++) t[j] += ctx->M[i+j*ld]*x2[i];
319: }
320: if (rk==np-1) for (i=0;i<n;i++) array2[nloc+i] = t[i];
321: PetscBLASIntCast(ld,&ld_);
322: PetscBLASIntCast(n,&n_);
323: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n_,ps,&ld_,t,&one));
324: BVMultVec(ctx->X,-1.0,1.0,ty,t);
325: PetscFree2(x2,t);
326: }
327: VecResetArray(tx);
328: VecResetArray(ty);
329: VecRestoreArrayRead(x,&array1);
330: VecRestoreArray(y,&array2);
331: return(0);
332: }
334: /*
335: Application of shell preconditioner:
336: y = B\x - eta*B\p, with eta = (u'*B\x)/(u'*B\p)
337: */
338: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
339: {
341: PetscScalar eta;
342: PEP_JD_PCSHELL *ctx;
345: PCShellGetContext(pc,(void**)&ctx);
347: /* y = B\x apply extended PC */
348: PEPJDExtendedPCApply(pc,x,y);
350: /* Compute eta = u'*y / u'*Bp */
351: VecDot(y,ctx->u,&eta);
352: eta /= ctx->gamma;
354: /* y = y - eta*Bp */
355: VecAXPY(y,-eta,ctx->Bp);
356: return(0);
357: }
359: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
360: {
362: PetscMPIInt np,rk,count;
363: PetscScalar *array1,*array2;
364: PetscInt nloc;
367: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
368: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
369: BVGetSizes(pep->V,&nloc,NULL,NULL);
370: if (v) {
371: VecGetArray(v,&array1);
372: VecGetArray(vex,&array2);
373: if (back) {
374: PetscMemcpy(array1,array2,nloc*sizeof(PetscScalar));
375: } else {
376: PetscMemcpy(array2,array1,nloc*sizeof(PetscScalar));
377: }
378: VecRestoreArray(v,&array1);
379: VecRestoreArray(vex,&array2);
380: }
381: if (a) {
382: if (rk==np-1) {
383: VecGetArray(vex,&array2);
384: if (back) {
385: PetscMemcpy(a,array2+nloc+off,na*sizeof(PetscScalar));
386: } else {
387: PetscMemcpy(array2+nloc+off,a,na*sizeof(PetscScalar));
388: }
389: VecRestoreArray(vex,&array2);
390: }
391: if (back) {
392: PetscMPIIntCast(na,&count);
393: MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
394: }
395: }
396: return(0);
397: }
399: static PetscErrorCode PEPJDComputePResidual(PEP pep,Vec u,PetscScalar theta,Vec p,Vec *work)
400: {
401: PEP_JD *pjd = (PEP_JD*)pep->data;
403: PetscMPIInt rk,np,count;
404: Vec tu,tp,w;
405: PetscScalar *array1,*array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,sone=1.0,zero=0.0;
406: PetscInt i,j,nconv=pjd->nconv,nloc,deg=pep->nmat-1;
407: PetscBLASInt n,ld,one=1;
410: if (nconv>0) {
411: PetscMalloc4(nconv,&xx,nconv,&tt,nconv,&x2,nconv,&q);
412: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
413: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
414: if (rk==np-1) {
415: BVGetSizes(pep->V,&nloc,NULL,NULL);
416: VecGetArray(u,&array1);
417: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
418: VecRestoreArray(u,&array1);
419: }
420: PetscMPIIntCast(nconv,&count);
421: MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
422: }
423: tu = work[0];
424: tp = work[1];
425: w = work[2];
426: VecGetArray(u,&array1);
427: VecPlaceArray(tu,array1);
428: VecGetArray(p,&array2);
429: VecPlaceArray(tp,array2);
430: VecSet(tp,0.0);
431: for (i=1;i<pep->nmat;i++) {
432: MatMult(pep->A[i],tu,w);
433: VecAXPY(tp,fact*(PetscReal)i,w);
434: fact *= theta;
435: }
436: if (nconv) {
437: PetscBLASIntCast(nconv,&n);
438: PetscBLASIntCast(pep->nev,&ld);
439: for (j=0;j<nconv;j++) q[j] = x2[j];
440: fact = theta;
441: for (i=2;i<pep->nmat;i++) {
442: BVMultVec(pjd->AX[i],1.0,1.0,tp,q);
443: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
444: for (j=0;j<nconv;j++) q[j] += (PetscReal)i*fact*x2[j];
445: fact *= theta;
446: }
447: BVSetActiveColumns(pjd->X,0,nconv);
448: BVDotVec(pjd->X,tu,xx);
449: if (rk==np-1) {
450: y2 = array2+nloc;
451: for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
452: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld,&ld,y2,&one));
453: fact = theta;
454: for (j=2;j<deg;j++) {
455: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
456: for (i=0;i<nconv;i++) tt[i] += (PetscReal)j*fact*xx[i];
457: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
458: for (i=0;i<nconv;i++) y2[i] += tt[i];
459: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
460: for (i=0;i<nconv;i++) q[i] += (PetscReal)j*fact*x2[i];
461: fact *= theta;
462: }
463: }
464: PetscFree4(xx,tt,x2,q);
465: }
466: VecResetArray(tu);
467: VecRestoreArray(u,&array1);
468: VecResetArray(tp);
469: VecRestoreArray(p,&array2);
470: return(0);
471: }
473: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
474: {
475: PEP_JD *pjd = (PEP_JD*)pep->data;
477: PetscScalar *tt;
478: Vec vg,wg;
479: PetscInt i;
480: PetscReal norm;
483: PetscMalloc1(pep->nev-1,&tt);
484: if (pep->nini==0) {
485: BVSetRandomColumn(pjd->V,0);
486: for (i=0;i<pep->nev-1;i++) tt[i] = 0.0;
487: BVGetColumn(pjd->V,0,&vg);
488: PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,vg,PETSC_FALSE);
489: BVRestoreColumn(pjd->V,0,&vg);
490: BVNormColumn(pjd->V,0,NORM_2,&norm);
491: BVScaleColumn(pjd->V,0,1.0/norm);
492: BVGetColumn(pjd->V,0,&vg);
493: BVGetColumn(pjd->W,0,&wg);
494: VecSet(wg,0.0);
495: PEPJDComputePResidual(pep,vg,pep->target,wg,w);
496: BVRestoreColumn(pjd->W,0,&wg);
497: BVRestoreColumn(pjd->V,0,&vg);
498: BVNormColumn(pjd->W,0,NORM_2,&norm);
499: BVScaleColumn(pjd->W,0,1.0/norm);
500: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Support for initial vectors not implemented yet");
501: PetscFree(tt);
502: return(0);
503: }
505: static PetscErrorCode PEPJDShellMatMult(Mat P,Vec x,Vec y)
506: {
507: PetscErrorCode ierr;
508: PEP_JD_MATSHELL *matctx;
509: PEP_JD *pjd;
510: PetscMPIInt rk,np,count;
511: PetscInt i,j,nconv,nloc,nmat,ldt,deg;
512: Vec tx,ty;
513: PetscScalar *array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,theta,*yy=NULL,sone=1.0,zero=0.0;
514: PetscBLASInt n,ld,one=1;
515: const PetscScalar *array1;
518: MatShellGetContext(P,(void**)&matctx);
519: pjd = (PEP_JD*)(matctx->pep->data);
520: nconv = pjd->nconv;
521: theta = matctx->theta;
522: nmat = matctx->pep->nmat;
523: deg = nmat-1;
524: ldt = matctx->pep->nev;
525: if (nconv>0) {
526: PetscMalloc5(nconv,&tt,nconv,&x2,nconv,&q,nconv,&xx,nconv,&yy);
527: MPI_Comm_rank(PetscObjectComm((PetscObject)P),&rk);
528: MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
529: if (rk==np-1) {
530: BVGetSizes(matctx->pep->V,&nloc,NULL,NULL);
531: VecGetArrayRead(x,&array1);
532: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
533: VecRestoreArrayRead(x,&array1);
534: }
535: PetscMPIIntCast(nconv,&count);
536: MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)P));
537: }
538: tx = matctx->work[0];
539: ty = matctx->work[1];
540: VecGetArrayRead(x,&array1);
541: VecPlaceArray(tx,array1);
542: VecGetArray(y,&array2);
543: VecPlaceArray(ty,array2);
544: VecSet(ty,0.0);
545: MatMult(matctx->P,tx,ty);
546: if (nconv) {
547: PetscBLASIntCast(pjd->nconv,&n);
548: PetscBLASIntCast(ldt,&ld);
549: for (j=0;j<nconv;j++) q[j] = x2[j];
550: fact = theta;
551: for (i=1;i<nmat;i++) {
552: BVMultVec(pjd->AX[i],1.0,1.0,ty,q);
553: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
554: for (j=0;j<nconv;j++) q[j] += fact*x2[j];
555: fact *= theta;
556: }
557: BVSetActiveColumns(pjd->X,0,nconv);
558: BVDotVec(pjd->X,tx,xx);
559: if (rk==np-1) {
560: y2 = array2+nloc;
561: for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
562: fact = theta;
563: for (j=1;j<deg;j++) {
564: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
565: for (i=0;i<nconv;i++) tt[i] += fact*xx[i];
566: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
567: for (i=0;i<nconv;i++) y2[i] += tt[i];
568: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
569: for (i=0;i<nconv;i++) q[i] += fact*x2[i];
570: fact *= theta;
571: }
572: }
573: PetscFree5(tt,x2,q,xx,yy);
574: }
575: VecResetArray(tx);
576: VecRestoreArrayRead(x,&array1);
577: VecResetArray(ty);
578: VecRestoreArray(y,&array2);
579: return(0);
580: }
582: static PetscErrorCode PEPJDCreateShellPC(PEP pep)
583: {
584: PEP_JD *pjd = (PEP_JD*)pep->data;
585: PEP_JD_PCSHELL *pcctx;
586: PEP_JD_MATSHELL *matctx;
587: KSP ksp;
588: PetscInt nloc,mloc;
589: PetscMPIInt np,rk;
590: PetscErrorCode ierr;
593: PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
594: PCSetType(pjd->pcshell,PCSHELL);
595: PCShellSetName(pjd->pcshell,"PCPEPJD");
596: PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
597: PetscNew(&pcctx);
598: PCShellSetContext(pjd->pcshell,pcctx);
599: STGetKSP(pep->st,&ksp);
600: BVCreateVec(pjd->V,&pcctx->Bp);
601: KSPGetPC(ksp,&pcctx->pc);
602: PetscObjectReference((PetscObject)pcctx->pc);
603: MatGetLocalSize(pep->A[0],&mloc,&nloc);
604: if (pep->nev>1) {
605: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
606: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
607: if (rk==np-1) { nloc += pep->nev-1; mloc += pep->nev-1; }
608: }
609: PetscNew(&matctx);
610: MatCreateShell(PetscObjectComm((PetscObject)pep),nloc,mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
611: MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)())PEPJDShellMatMult);
612: matctx->pep = pep;
613: MatDuplicate(pep->A[0],MAT_DO_NOT_COPY_VALUES,&matctx->P);
614: PCSetOperators(pcctx->pc,matctx->P,matctx->P);
615: KSPSetPC(ksp,pjd->pcshell);
616: KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
617: if (pep->nev>1) {
618: PetscMalloc2(pep->nev*pep->nev,&pcctx->M,pep->nev*pep->nev,&pcctx->ps);
619: pcctx->X = pjd->X;
620: pcctx->ld = pep->nev;
621: }
622: return(0);
623: }
625: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
626: {
627: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GETRI) || defined(PETSC_MISSING_LAPACK_GETRF)
629: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GETRI/GETRF - Lapack routines are unavailable");
630: #else
632: PEP_JD *pjd = (PEP_JD*)pep->data;
633: PEP_JD_PCSHELL *pcctx;
634: PetscInt i,j,k,n=pjd->nconv,ld=pep->nev,deg=pep->nmat-1;
635: PetscScalar fact,*M,*ps,*work,*U,*V,*S,sone=1.0,zero=0.0;
636: PetscReal tol,maxeig=0.0,*sg,*rwork;
637: PetscBLASInt n_,info,ld_,*p,lw_,rk=0;
640: if (n) {
641: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
642: pcctx->n = n;
643: M = pcctx->M;
644: ps = pcctx->ps;
645: /* h, and q are vectors containing diagonal matrices */
646: PetscCalloc7(n*n,&U,n*n,&V,n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p);
647: /* pseudo-inverse */
648: for (j=0;j<n;j++) {
649: for (i=0;i<j;i++) S[n*j+i] = -pjd->T[pep->nev*j+i];
650: S[n*j+j] = theta-pjd->T[pep->nev*j+j];
651: }
652: PetscBLASIntCast(n,&n_);
653: PetscBLASIntCast(ld,&ld_);
654: lw_ = 10*n_;
655: #if !defined (PETSC_USE_COMPLEX)
656: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
657: #else
658: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
659: #endif
660: SlepcCheckLapackInfo("gesvd",info);
661: for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
662: tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
663: for (j=0;j<n;j++) {
664: if (sg[j]>tol) {
665: for (i=0;i<n;i++) U[j*n+i] /= sg[j];
666: rk++;
667: } else break;
668: }
669: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));
671: /* compute M */
672: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
673: fact = theta;
674: PetscMemzero(S,n*n*sizeof(PetscScalar));
675: for (j=0;j<n;j++) S[j*(n+1)] = 1.0; /* q=S */
676: for (k=0;k<deg;k++) {
677: for (j=0;j<n;j++) for (i=0;i<n;i++) V[j*n+i] = S[j*n+i] + M[j*ld+i]*fact;
678: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
679: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
680: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n_,&n_,&sone,pjd->T,&ld_,S,&n_));
681: for (j=0;j<n;j++) S[j*(n+1)] += fact;
682: fact *=theta;
683: }
684: /* inverse */
685: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
686: SlepcCheckLapackInfo("getrf",info);
687: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
688: SlepcCheckLapackInfo("getri",info);
689: PetscFree7(U,V,S,sg,work,rwork,p);
690: }
691: return(0);
692: #endif
693: }
695: static PetscErrorCode PEPJDPCMatSetUp(PEP pep,PetscScalar theta)
696: {
697: PetscErrorCode ierr;
698: PEP_JD *pjd = (PEP_JD*)pep->data;
699: PEP_JD_MATSHELL *matctx;
700: PEP_JD_PCSHELL *pcctx;
701: MatStructure str;
702: PetscScalar t;
703: PetscInt i;
706: MatShellGetContext(pjd->Pshell,(void**)&matctx);
707: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
708: STGetMatStructure(pep->st,&str);
709: MatCopy(pep->A[0],matctx->P,str);
710: t = theta;
711: for (i=1;i<pep->nmat;i++) {
712: if (t!=0.0) { MatAXPY(matctx->P,t,pep->A[i],str); }
713: t *= theta;
714: }
715: PCSetOperators(pcctx->pc,matctx->P,matctx->P);
716: PCSetUp(pcctx->pc);
717: matctx->theta = theta;
718: return(0);
719: }
721: static PetscErrorCode PEPJDEigenvectors(PEP pep)
722: {
723: #if defined(SLEPC_MISSING_LAPACK_TREVC)
725: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
726: #else
728: PEP_JD *pjd = (PEP_JD*)pep->data;
729: PetscBLASInt ld,nconv,info,nc;
730: PetscScalar *Z,*w;
731: PetscReal *wr,norm;
732: PetscInt i;
733: Mat U;
736: PetscMalloc3(pjd->nconv*pjd->nconv,&Z,3*pep->nev,&wr,2*pep->nev,&w);
737: PetscBLASIntCast(pep->nev,&ld);
738: PetscBLASIntCast(pjd->nconv,&nconv);
739: #if !defined(PETSC_USE_COMPLEX)
740: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
741: #else
742: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
743: #endif
744: SlepcCheckLapackInfo("trevc",info);
745: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
746: BVSetActiveColumns(pjd->X,0,pjd->nconv);
747: BVMultInPlace(pjd->X,U,0,pjd->nconv);
748: for (i=0;i<pjd->nconv;i++) {
749: BVNormColumn(pjd->X,i,NORM_2,&norm);
750: BVScaleColumn(pjd->X,i,1.0/norm);
751: }
752: MatDestroy(&U);
753: PetscFree3(Z,wr,w);
754: return(0);
755: #endif
756: }
758: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv)
759: {
760: #if defined(SLEPC_MISSING_LAPACK_TRTRI)
762: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TRTRI - Lapack routine is unavailable");
763: #else
765: PEP_JD *pjd = (PEP_JD*)pep->data;
766: PetscInt j,i,ldds,rk=0,*P,nvv=*nv;
767: Vec v,x;
768: PetscBLASInt n,ld,rk_,nv_,info,one=1;
769: PetscScalar sone=1.0,*Tj,*R,*r,*tt,*pX;
770: Mat X;
773: /* update AX and XpX */
774: BVGetColumn(pjd->X,pjd->nconv-1,&x);
775: for (j=0;j<pep->nmat;j++) {
776: BVGetColumn(pjd->AX[j],pjd->nconv-1,&v);
777: MatMult(pep->A[j],x,v);
778: BVRestoreColumn(pjd->AX[j],pjd->nconv-1,&v);
779: BVSetActiveColumns(pjd->AX[j],0,pjd->nconv);
780: }
781: BVRestoreColumn(pjd->X,pjd->nconv-1,&x);
782: BVDotColumn(pjd->X,(pjd->nconv-1),pjd->XpX+(pjd->nconv-1)*(pep->nev));
783: pjd->XpX[(pjd->nconv-1)*(1+pep->nev)] = 1.0;
784: for (j=0;j<pjd->nconv-1;j++) pjd->XpX[j*(pep->nev)+pjd->nconv-1] = PetscConj(pjd->XpX[(pjd->nconv-1)*(pep->nev)+j]);
786: /* Compute powers of T */
787: PetscBLASIntCast(pjd->nconv,&n);
788: PetscBLASIntCast(pep->nev,&ld);
789: PetscMemzero(pjd->Tj,pep->nev*pep->nev*pep->nmat*sizeof(PetscScalar));
790: Tj = pjd->Tj;
791: for (j=0;j<pjd->nconv;j++) Tj[(pep->nev+1)*j] = 1.0;
792: Tj = pjd->Tj+pep->nev*pep->nev;
793: PetscMemcpy(Tj,pjd->T,pep->nev*pjd->nconv*sizeof(PetscScalar));
794: for (j=2;j<pep->nmat;j++) {
795: PetscMemcpy(Tj+pep->nev*pep->nev,Tj,pep->nev*pjd->nconv*sizeof(PetscScalar));
796: Tj += pep->nev*pep->nev;
797: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n,&n,&sone,pjd->T,&ld,Tj,&ld));
798: }
800: /* Extend search space */
801: PetscCalloc4(nvv,&P,nvv*nvv,&R,nvv,&r,pep->nev-1,&tt);
802: DSGetLeadingDimension(pep->ds,&ldds);
803: DSGetArray(pep->ds,DS_MAT_X,&pX);
804: PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
805: for (i=0;i<rk-1;i++) r[i] = PetscConj(R[nvv*i]*pep->eigr[P[i+1]]); /* first row scaled with permuted diagonal */
806: PetscBLASIntCast(rk,&rk_);
807: PetscBLASIntCast(nvv,&nv_);
808: PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
809: SlepcCheckLapackInfo("trtri",info);
810: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r,&one));
811: for (i=0;i<rk;i++) r[i] = PetscConj(r[i]); /* revert */
812: BVSetActiveColumns(pjd->V,0,nvv);
813: for (j=0;j<rk-1;j++) {
814: PetscMemcpy(R+j*nvv,pX+(j+1)*ldds,nvv*sizeof(PetscScalar));
815: }
816: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
817: MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk-1,R,&X);
818: BVMultInPlace(pjd->V,X,0,rk-1);
819: MatDestroy(&X);
820: BVSetActiveColumns(pjd->V,0,rk-1);
821: for (j=0;j<rk-1;j++) {
822: BVGetColumn(pjd->V,j,&v);
823: PEPJDCopyToExtendedVec(pep,NULL,r+j,1,pjd->nconv-1,v,PETSC_FALSE);
824: BVRestoreColumn(pjd->V,j,&v);
825: }
826: BVOrthogonalize(pjd->V,NULL);
827: for (j=0;j<rk-1;j++) {
828: BVGetColumn(pjd->W,j,&v);
829: PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,v,PETSC_FALSE);
830: BVRestoreColumn(pjd->W,j,&v);
831: }
832: *nv = rk-1;
833: PetscFree4(P,R,r,tt);
834: #endif
835: return(0);
836: }
838: PetscErrorCode PEPSolve_JD(PEP pep)
839: {
840: PetscErrorCode ierr;
841: PEP_JD *pjd = (PEP_JD*)pep->data;
842: PetscInt k,nv,ld,minv,low,high,dim;
843: PetscScalar theta=0.0,*pX,*eig;
844: PetscReal norm,*res;
845: PetscBool lindep,initial=PETSC_FALSE,flglk=PETSC_FALSE,flgre=PETSC_FALSE;
846: Vec t,u,p,r,*ww=pep->work,v;
847: Mat G,X,Y;
848: KSP ksp;
849: PEP_JD_PCSHELL *pcctx;
850: PEP_JD_MATSHELL *matctx;
853: DSGetLeadingDimension(pep->ds,&ld);
854: PetscMalloc2(pep->ncv,&eig,pep->ncv,&res);
855: BVCreateVec(pjd->V,&u);
856: VecDuplicate(u,&p);
857: VecDuplicate(u,&r);
858: STGetKSP(pep->st,&ksp);
860: if (pep->nini) {
861: nv = pep->nini; initial = PETSC_TRUE;
862: } else {
863: theta = pep->target;
864: nv = 1;
865: }
866: PEPJDProcessInitialSpace(pep,ww);
867: BVCopyVec(pjd->V,0,u);
869: /* Restart loop */
870: while (pep->reason == PEP_CONVERGED_ITERATING) {
871: pep->its++;
873: low = (flglk || flgre)? 0: nv-1;
874: high = nv;
875: DSSetDimensions(pep->ds,nv,0,0,0);
876: BVSetActiveColumns(pjd->V,low,high);
877: PEPJDUpdateTV(pep,low,high,ww);
878: BVSetActiveColumns(pjd->W,low,high);
879: for (k=0;k<pep->nmat;k++) {
880: BVSetActiveColumns(pjd->TV[k],low,high);
881: DSGetMat(pep->ds,DSMatExtra[k],&G);
882: BVMatProject(pjd->TV[k],NULL,pjd->W,G);
883: DSRestoreMat(pep->ds,DSMatExtra[k],&G);
884: }
885: BVSetActiveColumns(pjd->V,0,nv);
886: BVSetActiveColumns(pjd->W,0,nv);
888: /* Solve projected problem */
889: if (nv>1 || initial ) {
890: DSSetState(pep->ds,DS_STATE_RAW);
891: DSSolve(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv);
892: DSSort(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv,NULL,NULL,NULL);
893: DSSynchronize(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv);
894: theta = pep->eigr[0];
895: #if !defined(PETSC_USE_COMPLEX)
896: if (PetscAbsScalar(pep->eigi[pep->nconv])!=0.0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"PJD solver not implemented for complex Ritz values in real arithmetic");
897: #endif
899: /* Compute Ritz vector u=V*X(:,1) */
900: DSGetArray(pep->ds,DS_MAT_X,&pX);
901: BVSetActiveColumns(pjd->V,0,nv);
902: BVMultVec(pjd->V,1.0,0.0,u,pX);
903: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
904: }
905: PEPJDUpdateExtendedPC(pep,theta);
907: /* Replace preconditioner with one containing projectors */
908: if (!pjd->pcshell) {
909: PEPJDCreateShellPC(pep);
910: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
911: MatShellGetContext(pjd->Pshell,(void**)&matctx);
912: matctx->work = ww;
913: pcctx->work = ww;
914: }
915: PEPJDPCMatSetUp(pep,theta);
917: /* Compute r and r' */
918: MatMult(pjd->Pshell,u,r);
919: PEPJDComputePResidual(pep,u,theta,p,ww);
920: pcctx->u = u;
922: /* Check convergence */
923: VecNorm(r,NORM_2,&norm);
924: (*pep->converged)(pep,theta,0,norm,&pep->errest[pep->nconv],pep->convergedctx);
925: (*pep->stopping)(pep,pep->its,pep->max_it,(pep->errest[pep->nconv]<pep->tol)?pjd->nconv+1:pjd->nconv,pep->nev,&pep->reason,pep->stoppingctx);
927: if (pep->errest[pep->nconv]<pep->tol) {
929: /* Ritz pair converged */
930: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
931: if (pep->nev>1) {
932: BVGetColumn(pjd->X,pjd->nconv,&v);
933: PEPJDCopyToExtendedVec(pep,v,pjd->T+pep->nev*pjd->nconv,pep->nev-1,0,u,PETSC_TRUE);
934: BVRestoreColumn(pjd->X,pjd->nconv,&v);
935: BVSetActiveColumns(pjd->X,0,pjd->nconv+1);
936: BVNormColumn(pjd->X,pjd->nconv,NORM_2,&norm);
937: BVScaleColumn(pjd->X,pjd->nconv,1.0/norm);
938: for (k=0;k<pjd->nconv;k++) pjd->T[pep->nev*pjd->nconv+k] /= norm;
939: pjd->T[(pep->nev+1)*pjd->nconv] = pep->eigr[0];
940: } else {
941: BVInsertVec(pep->V,pep->nconv,u);
942: }
943: pjd->nconv++;
945: if (pep->reason==PEP_CONVERGED_ITERATING) {
946: PEPJDLockConverged(pep,&nv);
947: BVCopyVec(pjd->V,nv-1,u);
948: if (nv==1) theta = pep->target;
949: }
950: flglk = PETSC_TRUE;
951: } else if (nv==pep->ncv-1) {
953: /* Basis full, force restart */
954: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
955: DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
956: DSGetArray(pep->ds,DS_MAT_X,&pX);
957: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
958: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
959: DSGetArray(pep->ds,DS_MAT_Y,&pX);
960: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
961: DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
962: DSGetMat(pep->ds,DS_MAT_X,&X);
963: BVMultInPlace(pjd->V,X,pep->nconv,minv);
964: MatDestroy(&X);
965: DSGetMat(pep->ds,DS_MAT_Y,&Y);
966: BVMultInPlace(pjd->W,Y,pep->nconv,minv);
967: MatDestroy(&Y);
968: nv = minv;
969: flgre = PETSC_TRUE;
970: } else {
971: /* Solve correction equation to expand basis */
972: PEPJDExtendedPCApply(pjd->pcshell,p,pcctx->Bp);
973: VecDot(pcctx->Bp,u,&pcctx->gamma);
974: BVGetColumn(pjd->V,nv,&t);
975: KSPSolve(ksp,r,t);
976: BVRestoreColumn(pjd->V,nv,&t);
977: BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
978: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
979: BVScaleColumn(pjd->V,nv,1.0/norm);
980: BVInsertVec(pjd->W,nv,r);
981: BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
982: if (lindep) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
983: BVScaleColumn(pjd->W,nv,1.0/norm);
984: nv++;
985: flglk = PETSC_FALSE;
986: flgre = PETSC_FALSE;
987: }
988: for (k=pjd->nconv;k<nv;k++) {
989: eig[k] = pep->eigr[k-pjd->nconv];
990: res[k] = pep->errest[k-pjd->nconv];
991: #if !defined(PETSC_USE_COMPLEX)
992: pep->eigi[k-pjd->nconv] = 0.0;
993: #endif
994: }
995: PEPMonitor(pep,pep->its,pjd->nconv,eig,pep->eigi,res,pjd->nconv+1);
996: }
997: if (pep->nev>1) {
998: if (pjd->nconv>0) { PEPJDEigenvectors(pep); }
999: for (k=0;k<pjd->nconv;k++) {
1000: BVGetColumn(pjd->X,k,&v);
1001: BVInsertVec(pep->V,k,v);
1002: BVRestoreColumn(pjd->X,k,&v);
1003: pep->eigr[k] = pjd->T[(pep->nev+1)*k];
1004: pep->eigi[k] = 0.0;
1005: }
1006: PetscFree2(pcctx->M,pcctx->ps);
1007: }
1008: pep->nconv = pjd->nconv;
1009: KSPSetPC(ksp,pcctx->pc);
1010: MatDestroy(&matctx->P);
1011: VecDestroy(&pcctx->Bp);
1012: MatDestroy(&pjd->Pshell);
1013: PCDestroy(&pcctx->pc);
1014: PetscFree(pcctx);
1015: PetscFree(matctx);
1016: PCDestroy(&pjd->pcshell);
1017: PetscFree2(eig,res);
1018: VecDestroy(&u);
1019: VecDestroy(&r);
1020: VecDestroy(&p);
1021: return(0);
1022: }
1024: PetscErrorCode PEPJDSetRestart_JD(PEP pep,PetscReal keep)
1025: {
1026: PEP_JD *pjd = (PEP_JD*)pep->data;
1029: if (keep==PETSC_DEFAULT) pjd->keep = 0.5;
1030: else {
1031: 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]");
1032: pjd->keep = keep;
1033: }
1034: return(0);
1035: }
1037: /*@
1038: PEPJDSetRestart - Sets the restart parameter for the Jacobi-Davidson
1039: method, in particular the proportion of basis vectors that must be kept
1040: after restart.
1042: Logically Collective on PEP
1044: Input Parameters:
1045: + pep - the eigenproblem solver context
1046: - keep - the number of vectors to be kept at restart
1048: Options Database Key:
1049: . -pep_jd_restart - Sets the restart parameter
1051: Notes:
1052: Allowed values are in the range [0.1,0.9]. The default is 0.5.
1054: Level: advanced
1056: .seealso: PEPJDGetRestart()
1057: @*/
1058: PetscErrorCode PEPJDSetRestart(PEP pep,PetscReal keep)
1059: {
1065: PetscTryMethod(pep,"PEPJDSetRestart_C",(PEP,PetscReal),(pep,keep));
1066: return(0);
1067: }
1069: PetscErrorCode PEPJDGetRestart_JD(PEP pep,PetscReal *keep)
1070: {
1071: PEP_JD *pjd = (PEP_JD*)pep->data;
1074: *keep = pjd->keep;
1075: return(0);
1076: }
1078: /*@
1079: PEPJDGetRestart - Gets the restart parameter used in the Jacobi-Davidson method.
1081: Not Collective
1083: Input Parameter:
1084: . pep - the eigenproblem solver context
1086: Output Parameter:
1087: . keep - the restart parameter
1089: Level: advanced
1091: .seealso: PEPJDSetRestart()
1092: @*/
1093: PetscErrorCode PEPJDGetRestart(PEP pep,PetscReal *keep)
1094: {
1100: PetscUseMethod(pep,"PEPJDGetRestart_C",(PEP,PetscReal*),(pep,keep));
1101: return(0);
1102: }
1104: PetscErrorCode PEPSetFromOptions_JD(PetscOptionItems *PetscOptionsObject,PEP pep)
1105: {
1107: PetscBool flg;
1108: PetscReal r1;
1111: PetscOptionsHead(PetscOptionsObject,"PEP JD Options");
1113: PetscOptionsReal("-pep_jd_restart","Proportion of vectors kept after restart","PEPJDSetRestart",0.5,&r1,&flg);
1114: if (flg) { PEPJDSetRestart(pep,r1); }
1116: PetscOptionsTail();
1117: return(0);
1118: }
1120: PetscErrorCode PEPView_JD(PEP pep,PetscViewer viewer)
1121: {
1123: PEP_JD *pjd = (PEP_JD*)pep->data;
1124: PetscBool isascii;
1127: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
1128: if (isascii) {
1129: PetscViewerASCIIPrintf(viewer," %d%% of basis vectors kept after restart\n",(int)(100*pjd->keep));
1130: }
1131: return(0);
1132: }
1134: PetscErrorCode PEPSetDefaultST_JD(PEP pep)
1135: {
1137: KSP ksp;
1140: if (!((PetscObject)pep->st)->type_name) {
1141: STSetType(pep->st,STPRECOND);
1142: STPrecondSetKSPHasMat(pep->st,PETSC_TRUE);
1143: }
1144: STSetTransform(pep->st,PETSC_FALSE);
1145: STGetKSP(pep->st,&ksp);
1146: if (!((PetscObject)ksp)->type_name) {
1147: KSPSetType(ksp,KSPBCGSL);
1148: KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
1149: }
1150: return(0);
1151: }
1153: PetscErrorCode PEPReset_JD(PEP pep)
1154: {
1156: PEP_JD *pjd = (PEP_JD*)pep->data;
1157: PetscInt i;
1160: for (i=0;i<pep->nmat;i++) {
1161: BVDestroy(pjd->TV+i);
1162: }
1163: BVDestroy(&pjd->W);
1164: if (pep->nev>1) {
1165: BVDestroy(&pjd->V);
1166: for (i=0;i<pep->nmat;i++) {
1167: BVDestroy(pjd->AX+i);
1168: }
1169: BVDestroy(&pjd->X);
1170: PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
1171: }
1172: PetscFree2(pjd->TV,pjd->AX);
1173: return(0);
1174: }
1176: PetscErrorCode PEPDestroy_JD(PEP pep)
1177: {
1181: PetscFree(pep->data);
1182: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
1183: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
1184: return(0);
1185: }
1187: PETSC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
1188: {
1189: PEP_JD *pjd;
1193: PetscNewLog(pep,&pjd);
1194: pep->data = (void*)pjd;
1196: pep->ops->solve = PEPSolve_JD;
1197: pep->ops->setup = PEPSetUp_JD;
1198: pep->ops->setfromoptions = PEPSetFromOptions_JD;
1199: pep->ops->destroy = PEPDestroy_JD;
1200: pep->ops->reset = PEPReset_JD;
1201: pep->ops->view = PEPView_JD;
1202: pep->ops->setdefaultst = PEPSetDefaultST_JD;
1204: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
1205: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
1206: return(0);
1207: }