Actual source code: pjd.c
slepc-3.7.1 2016-05-27
1: /*
3: SLEPc polynomial eigensolver: "jd"
5: Method: Jacobi-Davidson
7: Algorithm:
9: Jacobi-Davidson for polynomial eigenvalue problems.
10: Based on code contributed by the authors of [2] below.
12: References:
14: [1] G.L.G. Sleijpen et al., "Jacobi-Davidson type methods for
15: generalized eigenproblems and polynomial eigenproblems", BIT
16: 36(3):595-633, 1996.
18: [2] Feng-Nan Hwang, Zih-Hao Wei, Tsung-Ming Huang, Weichung Wang,
19: "A Parallel Additive Schwarz Preconditioned Jacobi-Davidson
20: Algorithm for Polynomial Eigenvalue Problems in Quantum Dot
21: Simulation", J. Comput. Phys. 229(8):2932-2947, 2010.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: SLEPc - Scalable Library for Eigenvalue Problem Computations
25: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
27: This file is part of SLEPc.
29: SLEPc is free software: you can redistribute it and/or modify it under the
30: terms of version 3 of the GNU Lesser General Public License as published by
31: the Free Software Foundation.
33: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
34: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
35: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
36: more details.
38: You should have received a copy of the GNU Lesser General Public License
39: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
40: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41: */
43: #include <slepc/private/pepimpl.h> /*I "slepcpep.h" I*/
44: #include pjdp.h
45: #include <slepcblaslapack.h>
49: /*
50: Duplicate and resize auxiliary basis
51: */
52: static PetscErrorCode PEPJDDuplicateBasis(PEP pep,BV *basis)
53: {
54: PetscErrorCode ierr;
55: PetscInt nloc,m;
56: PetscMPIInt rank,nproc;
57: BVType type;
58: BVOrthogType otype;
59: BVOrthogRefineType oref;
60: PetscReal oeta;
61: BVOrthogBlockType oblock;
64: if (pep->nev>1) {
65: BVCreate(PetscObjectComm((PetscObject)pep),basis);
66: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rank);
67: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&nproc);
68: BVGetSizes(pep->V,&nloc,NULL,&m);
69: if (rank==nproc-1) nloc += pep->nev-1;
70: BVSetSizes(*basis,nloc,PETSC_DECIDE,m);
71: BVGetType(pep->V,&type);
72: BVSetType(*basis,type);
73: BVGetOrthogonalization(pep->V,&otype,&oref,&oeta,&oblock);
74: BVSetOrthogonalization(*basis,otype,oref,oeta,oblock);
75: PetscObjectStateIncrease((PetscObject)*basis);
76: } else {
77: BVDuplicate(pep->V,basis);
78: }
79: return(0);
80: }
84: PetscErrorCode PEPSetUp_JD(PEP pep)
85: {
87: PEP_JD *pjd = (PEP_JD*)pep->data;
88: PetscBool isprecond,flg;
89: PetscInt i;
90: KSP ksp;
93: pep->lineariz = PETSC_FALSE;
94: PEPSetDimensions_Default(pep,pep->nev,&pep->ncv,&pep->mpd);
95: if (!pep->max_it) pep->max_it = PetscMax(100,2*pep->n/pep->ncv);
96: if (!pep->which) pep->which = PEP_TARGET_MAGNITUDE;
97: if (pep->which != PEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEPJD only supports which=target_magnitude");;
99: /* Set STPRECOND as the default ST */
100: if (!((PetscObject)pep->st)->type_name) {
101: STSetType(pep->st,STPRECOND);
102: }
103: PetscObjectTypeCompare((PetscObject)pep->st,STPRECOND,&isprecond);
104: if (!isprecond) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"JD only works with PRECOND spectral transformation");
106: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver not implemented for non-monomial bases");
107: STGetTransform(pep->st,&flg);
108: if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Solver requires the ST transformation flag unset, see STSetTransform()");
110: /* Set the default options of the KSP */
111: STGetKSP(pep->st,&ksp);
112: if (!((PetscObject)ksp)->type_name) {
113: KSPSetType(ksp,KSPBCGSL);
114: KSPSetTolerances(ksp,1e-5,PETSC_DEFAULT,PETSC_DEFAULT,100);
115: }
117: if (!pjd->keep) pjd->keep = 0.5;
119: PEPAllocateSolution(pep,0);
120: PEPSetWorkVecs(pep,5);
121: PetscMalloc2(pep->nmat,&pjd->TV,pep->nmat,&pjd->AX);
122: for (i=0;i<pep->nmat;i++) {
123: PEPJDDuplicateBasis(pep,pjd->TV+i);
124: }
125: PEPJDDuplicateBasis(pep,&pjd->W);
126: if (pep->nev>1) {
127: PEPJDDuplicateBasis(pep,&pjd->V);
128: BVSetFromOptions(pjd->V);
129: for (i=0;i<pep->nmat;i++) {
130: BVDuplicateResize(pep->V,pep->nev-1,pjd->AX+i);
131: }
132: BVDuplicateResize(pep->V,pep->nev,&pjd->X);
133: PetscCalloc3((pep->nev)*(pep->nev),&pjd->XpX,pep->nev*pep->nev,&pjd->T,pep->nev*pep->nev*pep->nmat,&pjd->Tj);
134: } else pjd->V = pep->V;
135: DSSetType(pep->ds,DSPEP);
136: DSPEPSetDegree(pep->ds,pep->nmat-1);
137: DSAllocate(pep->ds,pep->ncv);
138: return(0);
139: }
143: /*
144: Updates columns (low to (high-1)) of TV[i]
145: */
146: static PetscErrorCode PEPJDUpdateTV(PEP pep,PetscInt low,PetscInt high,Vec *w)
147: {
149: PEP_JD *pjd = (PEP_JD*)pep->data;
150: PetscInt pp,col,i,j,nloc,nconv,deg=pep->nmat-1;
151: Vec v1,v2,t1,t2;
152: PetscScalar *array1,*array2,*x2,*tt,*xx,*y2,zero=0.0,sone=1.0;
153: PetscMPIInt rk,np,count;
154: PetscBLASInt n,ld,one=1;
157: nconv = pjd->nconv;
158: PetscMalloc3(nconv,&tt,nconv,&x2,nconv,&xx);
159: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
160: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
161: BVGetSizes(pep->V,&nloc,NULL,NULL);
162: t1 = w[0];
163: t2 = w[1];
164: for (col=low;col<high;col++) {
165: BVGetColumn(pjd->V,col,&v1);
166: VecGetArray(v1,&array1);
167: if (nconv>0) {
168: if (rk==np-1) { for (i=0;i<nconv;i++) x2[i] = array1[nloc+i]; }
169: PetscMPIIntCast(nconv,&count);
170: MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
171: }
172: VecPlaceArray(t1,array1);
173: for (pp=0;pp<pep->nmat;pp++) {
174: BVGetColumn(pjd->TV[pp],col,&v2);
175: VecGetArray(v2,&array2);
176: VecPlaceArray(t2,array2);
177: MatMult(pep->A[pp],t1,t2);
178: if (nconv) {
179: PetscBLASIntCast(pjd->nconv,&n);
180: PetscBLASIntCast(pep->nev,&ld);
181: for (j=0;j<nconv;j++) tt[j] = x2[j];
182: for (i=pp+1;i<pep->nmat;i++) {
183: BVMultVec(pjd->AX[i],1.0,1.0,t2,tt);
184: if (i!=pep->nmat-1) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,tt,&one));
185: }
186: BVDotVec(pjd->X,t1,xx);
187: if (rk==np-1 && pp<deg) {
188: y2 = array2+nloc;
189: for (j=0;j<nconv;j++) { y2[j] = xx[j]; xx[j] = x2[j]; }
190: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*pp,&ld,y2,&one));
191: for (i=pp+1;i<pep->nmat-1;i++) {
192: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,xx,&one,&zero,tt,&one));
193: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*i,&ld,tt,&one));
194: for (j=0;j<nconv;j++) y2[j] += tt[j];
195: if (i<pep->nmat-2) PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,xx,&one));
196: }
197: }
198: }
199: VecResetArray(t2);
200: VecRestoreArray(v2,&array2);
201: BVRestoreColumn(pjd->TV[pp],col,&v2);
202: }
203: VecResetArray(t1);
204: VecRestoreArray(v1,&array1);
205: BVRestoreColumn(pjd->V,col,&v1);
206: }
207: PetscFree3(tt,x2,xx);
208: return(0);
209: }
213: /*
214: RRQR of X. Xin*P=Xou*R. Rank of R is rk
215: */
216: static PetscErrorCode PEPJDOrthogonalize(PetscInt row,PetscInt col,PetscScalar *X,PetscInt ldx,PetscInt *rk,PetscInt *P,PetscScalar *R,PetscInt ldr)
217: {
218: #if defined(SLEPC_MISSING_LAPACK_GEQP3)
220: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQP3 - Lapack routine is unavailable");
221: #else
223: PetscInt i,j,n,r;
224: PetscBLASInt row_,col_,ldx_,*p,lwork,info,n_;
225: PetscScalar *tau,*work;
226: PetscReal tol,*rwork;
229: PetscBLASIntCast(row,&row_);
230: PetscBLASIntCast(col,&col_);
231: PetscBLASIntCast(ldx,&ldx_);
232: n = PetscMin(row,col);
233: PetscBLASIntCast(n,&n_);
234: lwork = 3*col_+1;
235: PetscMalloc4(col,&p,n,&tau,lwork,&work,2*col,&rwork);
236: for (i=1;i<col;i++) p[i] = 0;
237: p[0] = 1;
239: /* rank revealing QR */
240: #if defined(PETSC_USE_COMPLEX)
241: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,rwork,&info));
242: #else
243: PetscStackCallBLAS("LAPACKgeqp3",LAPACKgeqp3_(&row_,&col_,X,&ldx_,p,tau,work,&lwork,&info));
244: #endif
245: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQP3 %d",info);
246: if (P) for (i=0;i<col;i++) P[i] = p[i];
248: /* rank computation */
249: tol = PetscMax(row,col)*PETSC_MACHINE_EPSILON*PetscAbsScalar(X[0]);
250: r = 1;
251: for (i=1;i<n;i++) {
252: if (PetscAbsScalar(X[i+ldx*i])>tol) r++;
253: else break;
254: }
255: if (rk) *rk=r;
257: /* copy upper triangular matrix if requested */
258: if (R) {
259: for (i=0;i<r;i++) {
260: PetscMemzero(R+i*ldr,r*sizeof(PetscScalar));
261: for (j=0;j<=i;j++) R[i*ldr+j] = X[i*ldx+j];
262: }
263: }
264: PetscStackCallBLAS("LAPACKungqr",LAPACKungqr_(&row_,&n_,&n_,X,&ldx_,tau,work,&lwork,&info));
265: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
266: PetscFree4(p,tau,work,rwork);
267: return(0);
268: #endif
269: }
273: /*
274: Application of extended preconditioner
275: */
276: static PetscErrorCode PEPJDExtendedPCApply(PC pc,Vec x,Vec y)
277: {
278: PetscInt i,j,nloc,n,ld;
279: PetscMPIInt rk,np,count;
280: Vec tx,ty;
281: PEP_JD_PCSHELL *ctx;
282: PetscErrorCode ierr;
283: const PetscScalar *array1;
284: PetscScalar *x2=NULL,*t=NULL,*ps,*array2;
285: PetscBLASInt one=1.0,ld_,n_;
288: PCShellGetContext(pc,(void**)&ctx);
289: n = ctx->n;
290: ps = ctx->ps;
291: ld = ctx->ld;
292: if (n) {
293: PetscMalloc2(n,&x2,n,&t);
294: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rk);
295: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&np);
296: if (rk==np-1) {
297: VecGetLocalSize(ctx->work[0],&nloc);
298: VecGetArrayRead(x,&array1);
299: for (i=0;i<n;i++) x2[i] = array1[nloc+i];
300: VecRestoreArrayRead(x,&array1);
301: }
302: PetscMPIIntCast(n,&count);
303: MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pc));
304: }
306: /* y = B\x apply PC */
307: tx = ctx->work[0];
308: ty = ctx->work[1];
309: VecGetArrayRead(x,&array1);
310: VecPlaceArray(tx,array1);
311: VecGetArray(y,&array2);
312: VecPlaceArray(ty,array2);
313: PCApply(ctx->pc,tx,ty);
314: if (n) {
315: for (j=0;j<n;j++) {
316: t[j] = 0.0;
317: for (i=0;i<n;i++) t[j] += ctx->M[i+j*ld]*x2[i];
318: }
319: if (rk==np-1) for (i=0;i<n;i++) array2[nloc+i] = t[i];
320: PetscBLASIntCast(ld,&ld_);
321: PetscBLASIntCast(n,&n_);
322: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n_,ps,&ld_,t,&one));
323: BVMultVec(ctx->X,-1.0,1.0,ty,t);
324: PetscFree2(x2,t);
325: }
326: VecResetArray(tx);
327: VecResetArray(ty);
328: VecRestoreArrayRead(x,&array1);
329: VecRestoreArray(y,&array2);
330: return(0);
331: }
335: /*
336: Application of shell preconditioner:
337: y = B\x - eta*B\p, with eta = (u'*B\x)/(u'*B\p)
338: */
339: static PetscErrorCode PCShellApply_PEPJD(PC pc,Vec x,Vec y)
340: {
342: PetscScalar eta;
343: PEP_JD_PCSHELL *ctx;
346: PCShellGetContext(pc,(void**)&ctx);
348: /* y = B\x apply extended PC */
349: PEPJDExtendedPCApply(pc,x,y);
351: /* Compute eta = u'*y / u'*Bp */
352: VecDot(y,ctx->u,&eta);
353: eta /= ctx->gamma;
354:
355: /* y = y - eta*Bp */
356: VecAXPY(y,-eta,ctx->Bp);
357: return(0);
358: }
362: static PetscErrorCode PEPJDCopyToExtendedVec(PEP pep,Vec v,PetscScalar *a,PetscInt na,PetscInt off,Vec vex,PetscBool back)
363: {
365: PetscMPIInt np,rk,count;
366: PetscScalar *array1,*array2;
367: PetscInt nloc;
370: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
371: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
372: BVGetSizes(pep->V,&nloc,NULL,NULL);
373: if (v) {
374: VecGetArray(v,&array1);
375: VecGetArray(vex,&array2);
376: if (back) {
377: PetscMemcpy(array1,array2,nloc*sizeof(PetscScalar));
378: } else {
379: PetscMemcpy(array2,array1,nloc*sizeof(PetscScalar));
380: }
381: VecRestoreArray(v,&array1);
382: VecRestoreArray(vex,&array2);
383: }
384: if (a) {
385: if (rk==np-1) {
386: VecGetArray(vex,&array2);
387: if (back) {
388: PetscMemcpy(a,array2+nloc+off,na*sizeof(PetscScalar));
389: } else {
390: PetscMemcpy(array2+nloc+off,a,na*sizeof(PetscScalar));
391: }
392: VecRestoreArray(vex,&array2);
393: }
394: if (back) {
395: PetscMPIIntCast(na,&count);
396: MPI_Bcast(a,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
397: }
398: }
399: return(0);
400: }
404: static PetscErrorCode PEPJDComputePResidual(PEP pep,Vec u,PetscScalar theta,Vec p,Vec *work)
405: {
406: PEP_JD *pjd = (PEP_JD*)pep->data;
408: PetscMPIInt rk,np,count;
409: Vec tu,tp,w;
410: PetscScalar *array1,*array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,sone=1.0,zero=0.0;
411: PetscInt i,j,nconv=pjd->nconv,nloc,deg=pep->nmat-1;
412: PetscBLASInt n,ld,one=1;
415: if (nconv>0) {
416: PetscMalloc4(nconv,&xx,nconv,&tt,nconv,&x2,nconv,&q);
417: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
418: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
419: if (rk==np-1) {
420: BVGetSizes(pep->V,&nloc,NULL,NULL);
421: VecGetArray(u,&array1);
422: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
423: VecRestoreArray(u,&array1);
424: }
425: PetscMPIIntCast(nconv,&count);
426: MPI_Bcast(x2,count,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)pep));
427: }
428: tu = work[0];
429: tp = work[1];
430: w = work[2];
431: VecGetArray(u,&array1);
432: VecPlaceArray(tu,array1);
433: VecGetArray(p,&array2);
434: VecPlaceArray(tp,array2);
435: VecSet(tp,0.0);
436: for (i=1;i<pep->nmat;i++) {
437: MatMult(pep->A[i],tu,w);
438: VecAXPY(tp,fact*(PetscReal)i,w);
439: fact *= theta;
440: }
441: if (nconv) {
442: PetscBLASIntCast(nconv,&n);
443: PetscBLASIntCast(pep->nev,&ld);
444: for (j=0;j<nconv;j++) q[j] = x2[j];
445: fact = theta;
446: for (i=2;i<pep->nmat;i++) {
447: BVMultVec(pjd->AX[i],1.0,1.0,tp,q);
448: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
449: for (j=0;j<nconv;j++) q[j] += (PetscReal)i*fact*x2[j];
450: fact *= theta;
451: }
452: BVSetActiveColumns(pjd->X,0,nconv);
453: BVDotVec(pjd->X,tu,xx);
454: if (rk==np-1) {
455: y2 = array2+nloc;
456: for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
457: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld,&ld,y2,&one));
458: fact = theta;
459: for (j=2;j<deg;j++) {
460: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
461: for (i=0;i<nconv;i++) tt[i] += (PetscReal)j*fact*xx[i];
462: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
463: for (i=0;i<nconv;i++) y2[i] += tt[i];
464: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
465: for (i=0;i<nconv;i++) q[i] += (PetscReal)j*fact*x2[i];
466: fact *= theta;
467: }
468: }
469: PetscFree4(xx,tt,x2,q);
470: }
471: VecResetArray(tu);
472: VecRestoreArray(u,&array1);
473: VecResetArray(tp);
474: VecRestoreArray(p,&array2);
475: return(0);
476: }
480: static PetscErrorCode PEPJDProcessInitialSpace(PEP pep,Vec *w)
481: {
482: PEP_JD *pjd = (PEP_JD*)pep->data;
484: PetscScalar *tt;
485: Vec vg,wg;
486: PetscInt i;
487: PetscReal norm;
490: PetscMalloc1(pep->nev-1,&tt);
491: if (pep->nini==0) {
492: BVSetRandomColumn(pjd->V,0);
493: for (i=0;i<pep->nev-1;i++) tt[i] = 0.0;
494: BVGetColumn(pjd->V,0,&vg);
495: PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,vg,PETSC_FALSE);
496: BVRestoreColumn(pjd->V,0,&vg);
497: BVNormColumn(pjd->V,0,NORM_2,&norm);
498: BVScaleColumn(pjd->V,0,1.0/norm);
499: BVGetColumn(pjd->V,0,&vg);
500: BVGetColumn(pjd->W,0,&wg);
501: VecSet(wg,0.0);
502: PEPJDComputePResidual(pep,vg,pep->target,wg,w);
503: BVRestoreColumn(pjd->W,0,&wg);
504: BVRestoreColumn(pjd->V,0,&vg);
505: BVNormColumn(pjd->W,0,NORM_2,&norm);
506: BVScaleColumn(pjd->W,0,1.0/norm);
507: } else {
508: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TO DO");
509: }
510: PetscFree(tt);
511: return(0);
512: }
516: static PetscErrorCode PEPJDShellMatMult(Mat P,Vec x,Vec y)
517: {
518: PetscErrorCode ierr;
519: PEP_JD_MATSHELL *matctx;
520: PEP_JD *pjd;
521: PetscMPIInt rk,np,count;
522: PetscInt i,j,nconv,nloc,nmat,ldt,deg;
523: Vec tx,ty;
524: PetscScalar *array2,*x2=NULL,*y2,fact=1.0,*q=NULL,*tt=NULL,*xx=NULL,theta,*yy=NULL,sone=1.0,zero=0.0;
525: PetscBLASInt n,ld,one=1;
526: const PetscScalar *array1;
529: MatShellGetContext(P,(void**)&matctx);
530: pjd = (PEP_JD*)(matctx->pep->data);
531: nconv = pjd->nconv;
532: theta = matctx->theta;
533: nmat = matctx->pep->nmat;
534: deg = nmat-1;
535: ldt = matctx->pep->nev;
536: if (nconv>0) {
537: PetscMalloc5(nconv,&tt,nconv,&x2,nconv,&q,nconv,&xx,nconv,&yy);
538: MPI_Comm_rank(PetscObjectComm((PetscObject)P),&rk);
539: MPI_Comm_size(PetscObjectComm((PetscObject)P),&np);
540: if (rk==np-1) {
541: BVGetSizes(matctx->pep->V,&nloc,NULL,NULL);
542: VecGetArrayRead(x,&array1);
543: for (i=0;i<nconv;i++) x2[i] = array1[nloc+i];
544: VecRestoreArrayRead(x,&array1);
545: }
546: PetscMPIIntCast(nconv,&count);
547: MPI_Bcast(x2,nconv,MPIU_SCALAR,np-1,PetscObjectComm((PetscObject)P));
548: }
549: tx = matctx->work[0];
550: ty = matctx->work[1];
551: VecGetArrayRead(x,&array1);
552: VecPlaceArray(tx,array1);
553: VecGetArray(y,&array2);
554: VecPlaceArray(ty,array2);
555: VecSet(ty,0.0);
556: MatMult(matctx->P,tx,ty);
557: if (nconv) {
558: PetscBLASIntCast(pjd->nconv,&n);
559: PetscBLASIntCast(ldt,&ld);
560: for (j=0;j<nconv;j++) q[j] = x2[j];
561: fact = theta;
562: for (i=1;i<nmat;i++) {
563: BVMultVec(pjd->AX[i],1.0,1.0,ty,q);
564: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
565: for (j=0;j<nconv;j++) q[j] += fact*x2[j];
566: fact *= theta;
567: }
568: BVSetActiveColumns(pjd->X,0,nconv);
569: BVDotVec(pjd->X,tx,xx);
570: if (rk==np-1) {
571: y2 = array2+nloc;
572: for (i=0;i<nconv;i++) { q[i] = x2[i]; y2[i] = xx[i]; }
573: fact = theta;
574: for (j=1;j<deg;j++) {
575: PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&n,&sone,pjd->XpX,&ld,q,&one,&zero,tt,&one));
576: for (i=0;i<nconv;i++) tt[i] += fact*xx[i];
577: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&n,pjd->Tj+ld*ld*j,&ld,tt,&one));
578: for (i=0;i<nconv;i++) y2[i] += tt[i];
579: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","N","N",&n,pjd->T,&ld,q,&one));
580: for (i=0;i<nconv;i++) q[i] += fact*x2[i];
581: fact *= theta;
582: }
583: }
584: PetscFree5(tt,x2,q,xx,yy);
585: }
586: VecResetArray(tx);
587: VecRestoreArrayRead(x,&array1);
588: VecResetArray(ty);
589: VecRestoreArray(y,&array2);
590: return(0);
591: }
595: static PetscErrorCode PEPJDCreateShellPC(PEP pep)
596: {
597: PEP_JD *pjd = (PEP_JD*)pep->data;
598: PEP_JD_PCSHELL *pcctx;
599: PEP_JD_MATSHELL *matctx;
600: KSP ksp;
601: PetscInt nloc,mloc;
602: PetscMPIInt np,rk;
603: PetscErrorCode ierr;
606: PCCreate(PetscObjectComm((PetscObject)pep),&pjd->pcshell);
607: PCSetType(pjd->pcshell,PCSHELL);
608: PCShellSetName(pjd->pcshell,"PCPEPJD");
609: PCShellSetApply(pjd->pcshell,PCShellApply_PEPJD);
610: PetscNew(&pcctx);
611: PCShellSetContext(pjd->pcshell,pcctx);
612: STGetKSP(pep->st,&ksp);
613: BVCreateVec(pjd->V,&pcctx->Bp);
614: KSPGetPC(ksp,&pcctx->pc);
615: PetscObjectReference((PetscObject)pcctx->pc);
616: MatGetLocalSize(pep->A[0],&mloc,&nloc);
617: if (pep->nev>1) {
618: MPI_Comm_rank(PetscObjectComm((PetscObject)pep),&rk);
619: MPI_Comm_size(PetscObjectComm((PetscObject)pep),&np);
620: if (rk==np-1) { nloc += pep->nev-1; mloc += pep->nev-1; }
621: }
622: PetscNew(&matctx);
623: MatCreateShell(PetscObjectComm((PetscObject)pep),nloc,mloc,PETSC_DETERMINE,PETSC_DETERMINE,matctx,&pjd->Pshell);
624: MatShellSetOperation(pjd->Pshell,MATOP_MULT,(void(*)())PEPJDShellMatMult);
625: matctx->pep = pep;
626: MatDuplicate(pep->A[0],MAT_DO_NOT_COPY_VALUES,&matctx->P);
627: PCSetOperators(pcctx->pc,matctx->P,matctx->P);
628: KSPSetPC(ksp,pjd->pcshell);
629: KSPSetOperators(ksp,pjd->Pshell,pjd->Pshell);
630: if (pep->nev>1) {
631: PetscMalloc2(pep->nev*pep->nev,&pcctx->M,pep->nev*pep->nev,&pcctx->ps);
632: pcctx->X = pjd->X;
633: pcctx->ld = pep->nev;
634: }
635: return(0);
636: }
640: static PetscErrorCode PEPJDUpdateExtendedPC(PEP pep,PetscScalar theta)
641: {
643: PEP_JD *pjd = (PEP_JD*)pep->data;
644: PEP_JD_PCSHELL *pcctx;
645: PetscInt i,j,k,n=pjd->nconv,ld=pep->nev,deg=pep->nmat-1;
646: PetscScalar fact,*M,*ps,*work,*U,*V,*S,sone=1.0,zero=0.0;
647: PetscReal tol,maxeig=0.0,*sg,*rwork;
648: PetscBLASInt n_,info,ld_,*p,lw_,rk=0;
651: #if defined(PETSC_MISSING_LAPACK_GESVD) || defined(PETSC_MISSING_LAPACK_GETRI) || defined(PETSC_MISSING_LAPACK_GETRF)
653: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GESVD/GETRI/GETRF - Lapack routine is unavailable");
654: #else
655: if (n) {
656: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
657: pcctx->n = n;
658: M = pcctx->M;
659: ps = pcctx->ps;
660: /* h, and q are vectors containing diagonal matrices */
661: PetscCalloc7(n*n,&U,n*n,&V,n*n,&S,n,&sg,10*n,&work,5*n,&rwork,n,&p);
662: /* pseudo-inverse */
663: for (j=0;j<n;j++) {
664: for (i=0;i<j;i++) S[n*j+i] = -pjd->T[pep->nev*j+i];
665: S[n*j+j] = theta-pjd->T[pep->nev*j+j];
666: }
667: PetscBLASIntCast(n,&n_);
668: PetscBLASIntCast(ld,&ld_);
669: lw_ = 10*n_;
670: #if !defined (PETSC_USE_COMPLEX)
671: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,&info));
672: #else
673: PetscStackCallBLAS("LAPACKgesvd",LAPACKgesvd_("S","S",&n_,&n_,S,&n_,sg,U,&n_,V,&n_,work,&lw_,rwork,&info));
674: #endif
675: for (i=0;i<n;i++) maxeig = PetscMax(maxeig,sg[i]);
676: tol = 10*PETSC_MACHINE_EPSILON*n*maxeig;
677: for (j=0;j<n;j++) {
678: if (sg[j]>tol) {
679: for (i=0;i<n;i++) U[j*n+i] /= sg[j];
680: rk++;
681: } else break;
682: }
683: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&rk,&sone,U,&n_,V,&n_,&zero,ps,&ld_));
685: /* compute M */
686: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,ps,&ld_,&zero,M,&ld_));
687: fact = theta;
688: PetscMemzero(S,n*n*sizeof(PetscScalar));
689: for (j=0;j<n;j++) S[j*(n+1)] = 1.0; /* q=S */
690: for (k=0;k<deg;k++) {
691: 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;
692: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&n_,&n_,&n_,&sone,pjd->XpX,&ld_,V,&n_,&zero,U,&n_));
693: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&n_,&n_,&n_,&sone,pjd->Tj+k*ld*ld,&ld_,U,&n_,&sone,M,&ld_));
694: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n_,&n_,&sone,pjd->T,&ld_,S,&n_));
695: for (j=0;j<n;j++) S[j*(n+1)] += fact;
696: fact *=theta;
697: }
698: /* inverse */
699: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&n_,&n_,M,&ld_,p,&info));
700: PetscStackCallBLAS("LAPACKgetri",LAPACKgetri_(&n_,M,&ld_,p,work,&n_,&info));
701: PetscFree7(U,V,S,sg,work,rwork,p);
702: }
703: return(0);
704: #endif
705: }
709: static PetscErrorCode PEPJDPCMatSetUp(PEP pep,PetscScalar theta)
710: {
711: PetscErrorCode ierr;
712: PEP_JD *pjd = (PEP_JD*)pep->data;
713: PEP_JD_MATSHELL *matctx;
714: PEP_JD_PCSHELL *pcctx;
715: MatStructure str;
716: PetscScalar t;
717: PetscInt i;
720: MatShellGetContext(pjd->Pshell,(void**)&matctx);
721: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
722: STGetMatStructure(pep->st,&str);
723: MatCopy(pep->A[0],matctx->P,str);
724: t = theta;
725: for (i=1;i<pep->nmat;i++) {
726: if (t!=0.0) { MatAXPY(matctx->P,t,pep->A[i],str); }
727: t *= theta;
728: }
729: PCSetOperators(pcctx->pc,matctx->P,matctx->P);
730: PCSetUp(pcctx->pc);
731: matctx->theta = theta;
732: return(0);
733: }
737: static PetscErrorCode PEPJDEigenvectors(PEP pep)
738: {
740: PEP_JD *pjd = (PEP_JD*)pep->data;
741: PetscBLASInt ld,nconv,info,nc;
742: PetscScalar *Z,*w;
743: PetscReal *wr,norm;
744: PetscInt i;
745: Mat U;
746:
748: #if defined(SLEPC_MISSING_LAPACK_TREVC)
750: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TREVC - Lapack routine is unavailable");
751: #else
752: PetscMalloc3(pjd->nconv*pjd->nconv,&Z,3*pep->nev,&wr,2*pep->nev,&w);
753: PetscBLASIntCast(pep->nev,&ld);
754: PetscBLASIntCast(pjd->nconv,&nconv);
755: #if !defined(PETSC_USE_COMPLEX)
756: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,wr,&info));
757: #else
758: PetscStackCallBLAS("LAPACKtrevc",LAPACKtrevc_("R","A",NULL,&nconv,pjd->T,&ld,NULL,&nconv,Z,&nconv,&nconv,&nc,w,wr,&info));
759: #endif
760: MatCreateSeqDense(PETSC_COMM_SELF,nconv,nconv,Z,&U);
761: BVSetActiveColumns(pjd->X,0,pjd->nconv);
762: BVMultInPlace(pjd->X,U,0,pjd->nconv);
763: for (i=0;i<pjd->nconv;i++) {
764: BVNormColumn(pjd->X,i,NORM_2,&norm);
765: BVScaleColumn(pjd->X,i,1.0/norm);
766: }
767: MatDestroy(&U);
768: PetscFree3(Z,wr,w);
769: return(0);
770: #endif
771: }
775: static PetscErrorCode PEPJDLockConverged(PEP pep,PetscInt *nv)
776: {
777: PetscErrorCode ierr;
778: PEP_JD *pjd = (PEP_JD*)pep->data;
779: PetscInt j,i,ldds,rk,*P,nvv=*nv;
780: Vec v,x;
781: PetscBLASInt n,ld,rk_,nv_,info,one=1;
782: PetscScalar sone=1.0,*Tj,*R,*r,*tt,*pX;
783: Mat X;
786: #if defined(SLEPC_MISSING_LAPACK_TRTRI)
788: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"TRTRI - Lapack routine is unavailable");
789: #else
790: /* update AX and XpX */
791: BVGetColumn(pjd->X,pjd->nconv-1,&x);
792: for (j=0;j<pep->nmat;j++) {
793: BVGetColumn(pjd->AX[j],pjd->nconv-1,&v);
794: MatMult(pep->A[j],x,v);
795: BVRestoreColumn(pjd->AX[j],pjd->nconv-1,&v);
796: BVSetActiveColumns(pjd->AX[j],0,pjd->nconv);
797: }
798: BVRestoreColumn(pjd->X,pjd->nconv-1,&x);
799: BVDotColumn(pjd->X,(pjd->nconv-1),pjd->XpX+(pjd->nconv-1)*(pep->nev));
800: pjd->XpX[(pjd->nconv-1)*(1+pep->nev)] = 1.0;
801: 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]);
802:
803: /* Compute powers of T */
804: PetscBLASIntCast(pjd->nconv,&n);
805: PetscBLASIntCast(pep->nev,&ld);
806: PetscMemzero(pjd->Tj,pep->nev*pep->nev*pep->nmat*sizeof(PetscScalar));
807: Tj = pjd->Tj;
808: for (j=0;j<pjd->nconv;j++) Tj[(pep->nev+1)*j] = 1.0;
809: Tj = pjd->Tj+pep->nev*pep->nev;
810: PetscMemcpy(Tj,pjd->T,pep->nev*pjd->nconv*sizeof(PetscScalar));
811: for (j=2;j<pep->nmat;j++) {
812: PetscMemcpy(Tj+pep->nev*pep->nev,Tj,pep->nev*pjd->nconv*sizeof(PetscScalar));
813: Tj += pep->nev*pep->nev;
814: PetscStackCallBLAS("BLAStrmm",BLAStrmm_("L","U","N","N",&n,&n,&sone,pjd->T,&ld,Tj,&ld));
815: }
817: /* Extend search space */
818: PetscCalloc4(nvv,&P,nvv*nvv,&R,nvv,&r,pep->nev-1,&tt);
819: DSGetLeadingDimension(pep->ds,&ldds);
820: DSGetArray(pep->ds,DS_MAT_X,&pX);
821: PEPJDOrthogonalize(nvv,nvv,pX,ldds,&rk,P,R,nvv);
822: for (i=0;i<rk-1;i++) r[i] = PetscConj(R[nvv*i]*pep->eigr[P[i+1]]); /* first row scaled with permuted diagonal */
823: PetscBLASIntCast(rk,&rk_);
824: PetscBLASIntCast(nvv,&nv_);
825: PetscStackCallBLAS("LAPACKtrtri",LAPACKtrtri_("U","N",&rk_,R,&nv_,&info));
826: if (info) SETERRQ1(PETSC_COMM_SELF,1,"Error in xTRTRI, info=%D",(PetscInt)info);
827: PetscStackCallBLAS("BLAStrmv",BLAStrmv_("U","C","N",&rk_,R,&nv_,r,&one));
828: for (i=0;i<rk;i++) r[i] = PetscConj(r[i]); /* revert */
829: BVSetActiveColumns(pjd->V,0,nvv);
830: for (j=0;j<rk-1;j++) {
831: PetscMemcpy(R+j*nvv,pX+(j+1)*ldds,nvv*sizeof(PetscScalar));
832: }
833: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
834: MatCreateSeqDense(PETSC_COMM_SELF,nvv,rk-1,R,&X);
835: BVMultInPlace(pjd->V,X,0,rk-1);
836: MatDestroy(&X);
837: BVSetActiveColumns(pjd->V,0,rk-1);
838: for (j=0;j<rk-1;j++) {
839: BVGetColumn(pjd->V,j,&v);
840: PEPJDCopyToExtendedVec(pep,NULL,r+j,1,pjd->nconv-1,v,PETSC_FALSE);
841: BVRestoreColumn(pjd->V,j,&v);
842: }
843: BVOrthogonalize(pjd->V,NULL);
844: for (j=0;j<rk-1;j++) {
845: BVGetColumn(pjd->W,j,&v);
846: PEPJDCopyToExtendedVec(pep,NULL,tt,pep->nev-1,0,v,PETSC_FALSE);
847: BVRestoreColumn(pjd->W,j,&v);
848: }
849: *nv = rk-1;
850: PetscFree4(P,R,r,tt);
851: #endif
852: return(0);
853: }
857: PetscErrorCode PEPSolve_JD(PEP pep)
858: {
859: PetscErrorCode ierr;
860: PEP_JD *pjd = (PEP_JD*)pep->data;
861: PetscInt k,nv,ld,minv,low,high,dim;
862: PetscScalar theta=0.0,*pX,*eig;
863: PetscReal norm,*res;
864: PetscBool lindep,initial=PETSC_FALSE,flglk=PETSC_FALSE,flgre=PETSC_FALSE;
865: Vec t,u,p,r,*ww=pep->work,v;
866: Mat G,X,Y;
867: KSP ksp;
868: PEP_JD_PCSHELL *pcctx;
869: PEP_JD_MATSHELL *matctx;
872: DSGetLeadingDimension(pep->ds,&ld);
873: PetscMalloc2(pep->ncv,&eig,pep->ncv,&res);
874: BVCreateVec(pjd->V,&u);
875: VecDuplicate(u,&p);
876: VecDuplicate(u,&r);
877: STGetKSP(pep->st,&ksp);
879: if (pep->nini) {
880: nv = pep->nini; initial = PETSC_TRUE;
881: } else {
882: theta = pep->target;
883: nv = 1;
884: }
885: PEPJDProcessInitialSpace(pep,ww);
886: BVCopyVec(pjd->V,0,u);
888: /* Restart loop */
889: while (pep->reason == PEP_CONVERGED_ITERATING) {
890: pep->its++;
892: low = (flglk || flgre)? 0: nv-1;
893: high = nv;
894: DSSetDimensions(pep->ds,nv,0,0,0);
895: BVSetActiveColumns(pjd->V,low,high);
896: PEPJDUpdateTV(pep,low,high,ww);
897: BVSetActiveColumns(pjd->W,low,high);
898: for (k=0;k<pep->nmat;k++) {
899: BVSetActiveColumns(pjd->TV[k],low,high);
900: DSGetMat(pep->ds,DSMatExtra[k],&G);
901: BVMatProject(pjd->TV[k],NULL,pjd->W,G);
902: DSRestoreMat(pep->ds,DSMatExtra[k],&G);
903: }
904: BVSetActiveColumns(pjd->V,0,nv);
905: BVSetActiveColumns(pjd->W,0,nv);
907: /* Solve projected problem */
908: if (nv>1 || initial ) {
909: DSSetState(pep->ds,DS_STATE_RAW);
910: DSSolve(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv);
911: DSSort(pep->ds,pep->eigr+pep->nconv,pep->eigi+pep->nconv,NULL,NULL,NULL);
912: theta = pep->eigr[0];
913: #if !defined(PETSC_USE_COMPLEX)
914: 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");
915: #endif
917: /* Compute Ritz vector u=V*X(:,1) */
918: DSGetArray(pep->ds,DS_MAT_X,&pX);
919: BVSetActiveColumns(pjd->V,0,nv);
920: BVMultVec(pjd->V,1.0,0.0,u,pX);
921: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
922: }
923: PEPJDUpdateExtendedPC(pep,theta);
925: /* Replace preconditioner with one containing projectors */
926: if (!pjd->pcshell) {
927: PEPJDCreateShellPC(pep);
928: PCShellGetContext(pjd->pcshell,(void**)&pcctx);
929: MatShellGetContext(pjd->Pshell,(void**)&matctx);
930: matctx->work = ww;
931: pcctx->work = ww;
932: }
933: PEPJDPCMatSetUp(pep,theta);
934:
935: /* Compute r and r' */
936: MatMult(pjd->Pshell,u,r);
937: PEPJDComputePResidual(pep,u,theta,p,ww);
938: pcctx->u = u;
940: /* Check convergence */
941: VecNorm(r,NORM_2,&norm);
942: (*pep->converged)(pep,theta,0,norm,&pep->errest[pep->nconv],pep->convergedctx);
943: (*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);
945: if (pep->errest[pep->nconv]<pep->tol) {
947: /* Ritz pair converged */
948: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
949: if (pep->nev>1) {
950: BVGetColumn(pjd->X,pjd->nconv,&v);
951: PEPJDCopyToExtendedVec(pep,v,pjd->T+pep->nev*pjd->nconv,pep->nev-1,0,u,PETSC_TRUE);
952: BVRestoreColumn(pjd->X,pjd->nconv,&v);
953: BVSetActiveColumns(pjd->X,0,pjd->nconv+1);
954: BVNormColumn(pjd->X,pjd->nconv,NORM_2,&norm);
955: BVScaleColumn(pjd->X,pjd->nconv,1.0/norm);
956: for (k=0;k<pjd->nconv;k++) pjd->T[pep->nev*pjd->nconv+k] /= norm;
957: pjd->T[(pep->nev+1)*pjd->nconv] = pep->eigr[0];
958: } else {
959: BVInsertVec(pep->V,pep->nconv,u);
960: }
961: pjd->nconv++;
963: if (pep->reason==PEP_CONVERGED_ITERATING) {
964: PEPJDLockConverged(pep,&nv);
965: BVCopyVec(pjd->V,nv-1,u);
966: if (nv==1) theta = pep->target;
967: }
968: flglk = PETSC_TRUE;
969: } else if (nv==pep->ncv-1) {
971: /* Basis full, force restart */
972: minv = PetscMin(nv,(PetscInt)(pjd->keep*pep->ncv));
973: DSGetDimensions(pep->ds,&dim,NULL,NULL,NULL,NULL);
974: DSGetArray(pep->ds,DS_MAT_X,&pX);
975: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
976: DSRestoreArray(pep->ds,DS_MAT_X,&pX);
977: DSGetArray(pep->ds,DS_MAT_Y,&pX);
978: PEPJDOrthogonalize(dim,minv,pX,ld,&minv,NULL,NULL,ld);
979: DSRestoreArray(pep->ds,DS_MAT_Y,&pX);
980: DSGetMat(pep->ds,DS_MAT_X,&X);
981: BVMultInPlace(pjd->V,X,pep->nconv,minv);
982: DSRestoreMat(pep->ds,DS_MAT_X,&X);
983: DSGetMat(pep->ds,DS_MAT_Y,&Y);
984: BVMultInPlace(pjd->W,Y,pep->nconv,minv);
985: DSRestoreMat(pep->ds,DS_MAT_Y,&Y);
986: nv = minv;
987: flgre = PETSC_TRUE;
988: } else {
989: /* Solve correction equation to expand basis */
990: PEPJDExtendedPCApply(pjd->pcshell,p,pcctx->Bp);
991: VecDot(pcctx->Bp,u,&pcctx->gamma);
992: BVGetColumn(pjd->V,nv,&t);
993: KSPSolve(ksp,r,t);
994: BVRestoreColumn(pjd->V,nv,&t);
995: BVOrthogonalizeColumn(pjd->V,nv,NULL,&norm,&lindep);
996: if (lindep || norm==0.0) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
997: BVScaleColumn(pjd->V,nv,1.0/norm);
998: BVInsertVec(pjd->W,nv,r);
999: BVOrthogonalizeColumn(pjd->W,nv,NULL,&norm,&lindep);
1000: if (lindep) SETERRQ(PETSC_COMM_SELF,1,"Linearly dependent continuation vector");
1001: BVScaleColumn(pjd->W,nv,1.0/norm);
1002: nv++;
1003: flglk = PETSC_FALSE;
1004: flgre = PETSC_FALSE;
1005: }
1006: for (k=pjd->nconv;k<nv;k++) {
1007: eig[k] = pep->eigr[k-pjd->nconv];
1008: res[k] = pep->errest[k-pjd->nconv];
1009: #if !defined(PETSC_USE_COMPLEX)
1010: pep->eigi[k-pjd->nconv] = 0.0;
1011: #endif
1012: }
1013: PEPMonitor(pep,pep->its,pjd->nconv,eig,pep->eigi,res,pjd->nconv+1);
1014: }
1015: if (pep->nev>1) {
1016: if (pjd->nconv>0) { PEPJDEigenvectors(pep); }
1017: for (k=0;k<pjd->nconv;k++) {
1018: BVGetColumn(pjd->X,k,&v);
1019: BVInsertVec(pep->V,k,v);
1020: BVRestoreColumn(pjd->X,k,&v);
1021: pep->eigr[k] = pjd->T[(pep->nev+1)*k];
1022: pep->eigi[k] = 0.0;
1023: }
1024: PetscFree2(pcctx->M,pcctx->ps);
1025: }
1026: pep->nconv = pjd->nconv;
1027: KSPSetPC(ksp,pcctx->pc);
1028: MatDestroy(&matctx->P);
1029: VecDestroy(&pcctx->Bp);
1030: MatDestroy(&pjd->Pshell);
1031: PCDestroy(&pcctx->pc);
1032: PetscFree(pcctx);
1033: PetscFree(matctx);
1034: PCDestroy(&pjd->pcshell);
1035: PetscFree2(eig,res);
1036: VecDestroy(&u);
1037: VecDestroy(&r);
1038: VecDestroy(&p);
1039: return(0);
1040: }
1044: PetscErrorCode PEPReset_JD(PEP pep)
1045: {
1047: PEP_JD *pjd = (PEP_JD*)pep->data;
1048: PetscInt i;
1051: for (i=0;i<pep->nmat;i++) {
1052: BVDestroy(pjd->TV+i);
1053: }
1054: BVDestroy(&pjd->W);
1055: if (pep->nev>1) {
1056: BVDestroy(&pjd->V);
1057: for (i=0;i<pep->nmat;i++) {
1058: BVDestroy(pjd->AX+i);
1059: }
1060: BVDestroy(&pjd->X);
1061: PetscFree3(pjd->XpX,pjd->T,pjd->Tj);
1062: }
1063: PetscFree2(pjd->TV,pjd->AX);
1064: return(0);
1065: }
1069: PetscErrorCode PEPDestroy_JD(PEP pep)
1070: {
1074: PetscFree(pep->data);
1075: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",NULL);
1076: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",NULL);
1077: return(0);
1078: }
1082: PETSC_EXTERN PetscErrorCode PEPCreate_JD(PEP pep)
1083: {
1084: PEP_JD *pjd;
1088: PetscNewLog(pep,&pjd);
1089: pep->data = (void*)pjd;
1091: pjd->keep = 0;
1092: pep->ops->solve = PEPSolve_JD;
1093: pep->ops->setup = PEPSetUp_JD;
1094: pep->ops->setfromoptions = PEPSetFromOptions_JD;
1095: pep->ops->reset = PEPReset_JD;
1096: pep->ops->destroy = PEPDestroy_JD;
1097: pep->ops->view = PEPView_JD;
1098: PetscObjectComposeFunction((PetscObject)pep,"PEPJDSetRestart_C",PEPJDSetRestart_JD);
1099: PetscObjectComposeFunction((PetscObject)pep,"PEPJDGetRestart_C",PEPJDGetRestart_JD);
1100: return(0);
1101: }