Actual source code: dvdimprovex.c
slepc-3.18.0 2022-10-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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 eigensolver: "davidson"
13: Step: improve the eigenvectors X
14: */
16: #include "davidson.h"
17: #include <slepcblaslapack.h>
19: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r ****/
21: typedef struct {
22: PetscInt size_X;
23: KSP ksp; /* correction equation solver */
24: Vec friends; /* reference vector for composite vectors */
25: PetscScalar theta[4],thetai[2]; /* the shifts used in the correction eq. */
26: PetscInt maxits; /* maximum number of iterations */
27: PetscInt r_s,r_e; /* the selected eigenpairs to improve */
28: PetscInt ksp_max_size; /* the ksp maximum subvectors size */
29: PetscReal tol; /* the maximum solution tolerance */
30: PetscReal lastTol; /* last tol for dynamic stopping criterion */
31: PetscReal fix; /* tolerance for using the approx. eigenvalue */
32: PetscBool dynamic; /* if the dynamic stopping criterion is applied */
33: dvdDashboard *d; /* the current dvdDashboard reference */
34: PC old_pc; /* old pc in ksp */
35: BV KZ; /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
36: BV U; /* new X vectors */
37: PetscScalar *XKZ; /* X'*KZ */
38: PetscScalar *iXKZ; /* inverse of XKZ */
39: PetscInt ldXKZ; /* leading dimension of XKZ */
40: PetscInt size_iXKZ; /* size of iXKZ */
41: PetscInt ldiXKZ; /* leading dimension of iXKZ */
42: PetscInt size_cX; /* last value of d->size_cX */
43: PetscInt old_size_X; /* last number of improved vectors */
44: PetscBLASInt *iXKZPivots; /* array of pivots */
45: } dvdImprovex_jd;
47: /*
48: Compute (I - KZ*iXKZ*X')*V where,
49: V, the vectors to apply the projector,
50: cV, the number of vectors in V,
51: */
52: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV)
53: {
54: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
55: PetscInt i,ldh,k,l;
56: PetscScalar *h;
57: PetscBLASInt cV_,n,info,ld;
58: #if defined(PETSC_USE_COMPLEX)
59: PetscInt j;
60: #endif
62: PetscAssert(cV<=2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
64: /* h <- X'*V */
65: PetscMalloc1(data->size_iXKZ*cV,&h);
66: ldh = data->size_iXKZ;
67: BVGetActiveColumns(data->U,&l,&k);
68: PetscAssert(ldh==k,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
69: BVSetActiveColumns(data->U,0,k);
70: for (i=0;i<cV;i++) {
71: BVDotVec(data->U,V[i],&h[ldh*i]);
72: #if defined(PETSC_USE_COMPLEX)
73: for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
74: #endif
75: }
76: BVSetActiveColumns(data->U,l,k);
78: /* h <- iXKZ\h */
79: PetscBLASIntCast(cV,&cV_);
80: PetscBLASIntCast(data->size_iXKZ,&n);
81: PetscBLASIntCast(data->ldiXKZ,&ld);
82: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
83: PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
84: PetscFPTrapPop();
85: SlepcCheckLapackInfo("getrs",info);
87: /* V <- V - KZ*h */
88: BVSetActiveColumns(data->KZ,0,k);
89: for (i=0;i<cV;i++) BVMultVec(data->KZ,-1.0,1.0,V[i],&h[ldh*i]);
90: BVSetActiveColumns(data->KZ,l,k);
91: PetscFree(h);
92: return 0;
93: }
95: /*
96: Compute (I - X*iXKZ*KZ')*V where,
97: V, the vectors to apply the projector,
98: cV, the number of vectors in V,
99: */
100: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV)
101: {
102: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
103: PetscInt i,ldh,k,l;
104: PetscScalar *h;
105: PetscBLASInt cV_, n, info, ld;
106: #if defined(PETSC_USE_COMPLEX)
107: PetscInt j;
108: #endif
110: PetscAssert(cV<=2,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
112: /* h <- KZ'*V */
113: PetscMalloc1(data->size_iXKZ*cV,&h);
114: ldh = data->size_iXKZ;
115: BVGetActiveColumns(data->U,&l,&k);
116: PetscAssert(ldh==k,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
117: BVSetActiveColumns(data->KZ,0,k);
118: for (i=0;i<cV;i++) {
119: BVDotVec(data->KZ,V[i],&h[ldh*i]);
120: #if defined(PETSC_USE_COMPLEX)
121: for (j=0;j<k;j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
122: #endif
123: }
124: BVSetActiveColumns(data->KZ,l,k);
126: /* h <- iXKZ\h */
127: PetscBLASIntCast(cV,&cV_);
128: PetscBLASIntCast(data->size_iXKZ,&n);
129: PetscBLASIntCast(data->ldiXKZ,&ld);
130: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
131: PetscCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
132: PetscFPTrapPop();
133: SlepcCheckLapackInfo("getrs",info);
135: /* V <- V - U*h */
136: BVSetActiveColumns(data->U,0,k);
137: for (i=0;i<cV;i++) BVMultVec(data->U,-1.0,1.0,V[i],&h[ldh*i]);
138: BVSetActiveColumns(data->U,l,k);
139: PetscFree(h);
140: return 0;
141: }
143: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
144: {
145: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
147: VecDestroy(&data->friends);
149: /* Restore the pc of ksp */
150: if (data->old_pc) {
151: KSPSetPC(data->ksp, data->old_pc);
152: PCDestroy(&data->old_pc);
153: }
154: return 0;
155: }
157: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
158: {
159: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
161: /* Free local data and objects */
162: PetscFree(data->XKZ);
163: PetscFree(data->iXKZ);
164: PetscFree(data->iXKZPivots);
165: BVDestroy(&data->KZ);
166: BVDestroy(&data->U);
167: PetscFree(data);
168: return 0;
169: }
171: /*
172: y <- theta[1]A*x - theta[0]*B*x
173: auxV, two auxiliary vectors
174: */
175: static inline PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y)
176: {
177: PetscInt n,i;
178: const Vec *Bx;
179: Vec *auxV;
181: n = data->r_e - data->r_s;
182: for (i=0;i<n;i++) MatMult(data->d->A,x[i],y[i]);
184: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
185: for (i=0;i<n;i++) {
186: #if !defined(PETSC_USE_COMPLEX)
187: if (PetscUnlikely(data->d->eigi[data->r_s+i] != 0.0)) {
188: if (data->d->B) {
189: MatMult(data->d->B,x[i],auxV[0]);
190: MatMult(data->d->B,x[i+1],auxV[1]);
191: Bx = auxV;
192: } else Bx = &x[i];
194: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i + ti_i*Bx_i+1;
195: y_i+1 t_2i+1*A*x_i+1 - ti_i*Bx_i - t_2i*Bx_i+1 ] */
196: VecAXPBYPCZ(y[i],-data->theta[2*i],data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
197: VecAXPBYPCZ(y[i+1],-data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
198: i++;
199: } else
200: #endif
201: {
202: if (data->d->B) {
203: MatMult(data->d->B,x[i],auxV[0]);
204: Bx = auxV;
205: } else Bx = &x[i];
206: VecAXPBY(y[i],-data->theta[i*2],data->theta[i*2+1],Bx[0]);
207: }
208: }
209: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
210: return 0;
211: }
213: /*
214: y <- theta[1]'*A'*x - theta[0]'*B'*x
215: */
216: static inline PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y)
217: {
218: PetscInt n,i;
219: const Vec *Bx;
220: Vec *auxV;
222: n = data->r_e - data->r_s;
223: for (i=0;i<n;i++) MatMultTranspose(data->d->A,x[i],y[i]);
225: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
226: for (i=0;i<n;i++) {
227: #if !defined(PETSC_USE_COMPLEX)
228: if (data->d->eigi[data->r_s+i] != 0.0) {
229: if (data->d->B) {
230: MatMultTranspose(data->d->B,x[i],auxV[0]);
231: MatMultTranspose(data->d->B,x[i+1],auxV[1]);
232: Bx = auxV;
233: } else Bx = &x[i];
235: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i - ti_i*Bx_i+1;
236: y_i+1 t_2i+1*A*x_i+1 + ti_i*Bx_i - t_2i*Bx_i+1 ] */
237: VecAXPBYPCZ(y[i],-data->theta[2*i],-data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
238: VecAXPBYPCZ(y[i+1],data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
239: i++;
240: } else
241: #endif
242: {
243: if (data->d->B) {
244: MatMultTranspose(data->d->B,x[i],auxV[0]);
245: Bx = auxV;
246: } else Bx = &x[i];
247: VecAXPBY(y[i],PetscConj(-data->theta[i*2]),PetscConj(data->theta[i*2+1]),Bx[0]);
248: }
249: }
250: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
251: return 0;
252: }
254: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w)
255: {
256: dvdImprovex_jd *data;
257: PetscInt n,i;
258: const Vec *inx,*outx,*wx;
259: Vec *auxV;
260: Mat A;
262: PCGetOperators(pc,&A,NULL);
263: MatShellGetContext(A,&data);
264: VecCompGetSubVecs(in,NULL,&inx);
265: VecCompGetSubVecs(out,NULL,&outx);
266: VecCompGetSubVecs(w,NULL,&wx);
267: n = data->r_e - data->r_s;
268: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
269: switch (side) {
270: case PC_LEFT:
271: /* aux <- theta[1]A*in - theta[0]*B*in */
272: dvd_aux_matmult(data,inx,auxV);
274: /* out <- K * aux */
275: for (i=0;i<n;i++) data->d->improvex_precond(data->d,data->r_s+i,auxV[i],outx[i]);
276: break;
277: case PC_RIGHT:
278: /* aux <- K * in */
279: for (i=0;i<n;i++) data->d->improvex_precond(data->d,data->r_s+i,inx[i],auxV[i]);
281: /* out <- theta[1]A*auxV - theta[0]*B*auxV */
282: dvd_aux_matmult(data,auxV,outx);
283: break;
284: case PC_SYMMETRIC:
285: /* aux <- K^{1/2} * in */
286: for (i=0;i<n;i++) PCApplySymmetricRight(data->old_pc,inx[i],auxV[i]);
288: /* wx <- theta[1]A*auxV - theta[0]*B*auxV */
289: dvd_aux_matmult(data,auxV,wx);
291: /* aux <- K^{1/2} * in */
292: for (i=0;i<n;i++) PCApplySymmetricLeft(data->old_pc,wx[i],outx[i]);
293: break;
294: default:
295: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Unsupported KSP side");
296: }
297: /* out <- out - v*(u'*out) */
298: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
299: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
300: return 0;
301: }
303: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out)
304: {
305: dvdImprovex_jd *data;
306: PetscInt n,i;
307: const Vec *inx, *outx;
308: Mat A;
310: PCGetOperators(pc,&A,NULL);
311: MatShellGetContext(A,&data);
312: VecCompGetSubVecs(in,NULL,&inx);
313: VecCompGetSubVecs(out,NULL,&outx);
314: n = data->r_e - data->r_s;
315: /* out <- K * in */
316: for (i=0;i<n;i++) data->d->improvex_precond(data->d,data->r_s+i,inx[i],outx[i]);
317: /* out <- out - v*(u'*out) */
318: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
319: return 0;
320: }
322: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out)
323: {
324: dvdImprovex_jd *data;
325: PetscInt n,i;
326: const Vec *inx, *outx;
327: Vec *auxV;
328: Mat A;
330: PCGetOperators(pc,&A,NULL);
331: MatShellGetContext(A,&data);
332: VecCompGetSubVecs(in,NULL,&inx);
333: VecCompGetSubVecs(out,NULL,&outx);
334: n = data->r_e - data->r_s;
335: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
336: /* auxV <- in */
337: for (i=0;i<n;i++) VecCopy(inx[i],auxV[i]);
338: /* auxV <- auxV - u*(v'*auxV) */
339: dvd_improvex_applytrans_proj(data->d,auxV,n);
340: /* out <- K' * aux */
341: for (i=0;i<n;i++) PCApplyTranspose(data->old_pc,auxV[i],outx[i]);
342: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
343: return 0;
344: }
346: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out)
347: {
348: dvdImprovex_jd *data;
349: PetscInt n;
350: const Vec *inx, *outx;
351: PCSide side;
353: MatShellGetContext(A,&data);
354: VecCompGetSubVecs(in,NULL,&inx);
355: VecCompGetSubVecs(out,NULL,&outx);
356: n = data->r_e - data->r_s;
357: /* out <- theta[1]A*in - theta[0]*B*in */
358: dvd_aux_matmult(data,inx,outx);
359: KSPGetPCSide(data->ksp,&side);
360: if (side == PC_RIGHT) {
361: /* out <- out - v*(u'*out) */
362: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
363: }
364: return 0;
365: }
367: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out)
368: {
369: dvdImprovex_jd *data;
370: PetscInt n,i;
371: const Vec *inx,*outx,*r;
372: Vec *auxV;
373: PCSide side;
375: MatShellGetContext(A,&data);
376: VecCompGetSubVecs(in,NULL,&inx);
377: VecCompGetSubVecs(out,NULL,&outx);
378: n = data->r_e - data->r_s;
379: KSPGetPCSide(data->ksp,&side);
380: if (side == PC_RIGHT) {
381: /* auxV <- in */
382: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
383: for (i=0;i<n;i++) VecCopy(inx[i],auxV[i]);
384: /* auxV <- auxV - v*(u'*auxV) */
385: dvd_improvex_applytrans_proj(data->d,auxV,n);
386: r = auxV;
387: } else r = inx;
388: /* out <- theta[1]A*r - theta[0]*B*r */
389: dvd_aux_matmulttrans(data,r,outx);
390: if (side == PC_RIGHT) SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
391: return 0;
392: }
394: static PetscErrorCode MatCreateVecs_dvd_jd(Mat A,Vec *right,Vec *left)
395: {
396: Vec *r,*l;
397: dvdImprovex_jd *data;
398: PetscInt n,i;
400: MatShellGetContext(A,&data);
401: n = data->ksp_max_size;
402: if (right) PetscMalloc1(n,&r);
403: if (left) PetscMalloc1(n,&l);
404: for (i=0;i<n;i++) MatCreateVecs(data->d->A,right?&r[i]:NULL,left?&l[i]:NULL);
405: if (right) {
406: VecCreateCompWithVecs(r,n,data->friends,right);
407: for (i=0;i<n;i++) VecDestroy(&r[i]);
408: }
409: if (left) {
410: VecCreateCompWithVecs(l,n,data->friends,left);
411: for (i=0;i<n;i++) VecDestroy(&l[i]);
412: }
414: if (right) PetscFree(r);
415: if (left) PetscFree(l);
416: return 0;
417: }
419: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
420: {
421: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
422: PetscInt rA, cA, rlA, clA;
423: Mat A;
424: PetscBool t;
425: PC pc;
426: Vec v0[2];
428: data->size_cX = data->old_size_X = 0;
429: data->lastTol = data->dynamic?0.5:0.0;
431: /* Setup the ksp */
432: if (data->ksp) {
433: /* Create the reference vector */
434: BVGetColumn(d->eps->V,0,&v0[0]);
435: v0[1] = v0[0];
436: VecCreateCompWithVecs(v0,data->ksp_max_size,NULL,&data->friends);
437: BVRestoreColumn(d->eps->V,0,&v0[0]);
439: /* Save the current pc and set a PCNONE */
440: KSPGetPC(data->ksp, &data->old_pc);
441: PetscObjectTypeCompare((PetscObject)data->old_pc,PCNONE,&t);
442: data->lastTol = 0.5;
443: if (t) data->old_pc = NULL;
444: else {
445: PetscObjectReference((PetscObject)data->old_pc);
446: PCCreate(PetscObjectComm((PetscObject)d->eps),&pc);
447: PCSetType(pc,PCSHELL);
448: PCSetOperators(pc,d->A,d->A);
449: PCSetReusePreconditioner(pc,PETSC_TRUE);
450: PCShellSetApply(pc,PCApply_dvd);
451: PCShellSetApplyBA(pc,PCApplyBA_dvd);
452: PCShellSetApplyTranspose(pc,PCApplyTranspose_dvd);
453: KSPSetPC(data->ksp,pc);
454: PCDestroy(&pc);
455: }
457: /* Create the (I-v*u')*K*(A-s*B) matrix */
458: MatGetSize(d->A,&rA,&cA);
459: MatGetLocalSize(d->A,&rlA,&clA);
460: MatCreateShell(PetscObjectComm((PetscObject)d->A),rlA*data->ksp_max_size,clA*data->ksp_max_size,rA*data->ksp_max_size,cA*data->ksp_max_size,data,&A);
461: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_dvd_jd);
462: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_dvd_jd);
463: MatShellSetOperation(A,MATOP_CREATE_VECS,(void(*)(void))MatCreateVecs_dvd_jd);
465: /* Try to avoid KSPReset */
466: KSPGetOperatorsSet(data->ksp,&t,NULL);
467: if (t) {
468: Mat M;
469: PetscInt rM;
470: KSPGetOperators(data->ksp,&M,NULL);
471: MatGetSize(M,&rM,NULL);
472: if (rM != rA*data->ksp_max_size) KSPReset(data->ksp);
473: }
474: EPS_KSPSetOperators(data->ksp,A,A);
475: KSPSetReusePreconditioner(data->ksp,PETSC_TRUE);
476: KSPSetUp(data->ksp);
477: MatDestroy(&A);
478: } else {
479: data->old_pc = NULL;
480: data->friends = NULL;
481: }
482: BVSetActiveColumns(data->KZ,0,0);
483: BVSetActiveColumns(data->U,0,0);
484: return 0;
485: }
487: /*
488: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
489: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
490: where
491: pX,pY, the right and left eigenvectors of the projected system
492: ld, the leading dimension of pX and pY
493: */
494: static PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
495: {
496: PetscInt n=i_e-i_s,size_KZ,V_new,rm,i,lv,kv,lKZ,kKZ;
497: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
498: const PetscScalar *array;
499: Mat M;
500: Vec u[2],v[2];
501: PetscBLASInt s,ldXKZ,info;
503: /* Check consistency */
504: BVGetActiveColumns(d->eps->V,&lv,&kv);
505: V_new = lv - data->size_cX;
506: PetscAssert(V_new<=data->old_size_X,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Consistency broken");
507: data->old_size_X = n;
508: data->size_cX = lv;
510: /* KZ <- KZ(rm:rm+max_cX-1) */
511: BVGetActiveColumns(data->KZ,&lKZ,&kKZ);
512: rm = PetscMax(V_new+lKZ,0);
513: if (rm > 0) {
514: for (i=0;i<lKZ;i++) {
515: BVCopyColumn(data->KZ,i+rm,i);
516: BVCopyColumn(data->U,i+rm,i);
517: }
518: }
520: /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
521: if (rm > 0) {
522: for (i=0;i<lKZ;i++) PetscArraycpy(&data->XKZ[i*data->ldXKZ+i],&data->XKZ[(i+rm)*data->ldXKZ+i+rm],lKZ);
523: }
524: lKZ = PetscMin(0,lKZ+V_new);
525: BVSetActiveColumns(data->KZ,lKZ,lKZ+n);
526: BVSetActiveColumns(data->U,lKZ,lKZ+n);
528: /* Compute X, KZ and KR */
529: BVGetColumn(data->U,lKZ,u);
530: if (n>1) BVGetColumn(data->U,lKZ+1,&u[1]);
531: BVGetColumn(data->KZ,lKZ,v);
532: if (n>1) BVGetColumn(data->KZ,lKZ+1,&v[1]);
533: d->improvex_jd_proj_uv(d,i_s,i_e,u,v,kr,theta,thetai,pX,pY,ld);
534: BVRestoreColumn(data->U,lKZ,u);
535: if (n>1) BVRestoreColumn(data->U,lKZ+1,&u[1]);
536: BVRestoreColumn(data->KZ,lKZ,v);
537: if (n>1) BVRestoreColumn(data->KZ,lKZ+1,&v[1]);
539: /* XKZ <- U'*KZ */
540: MatCreateSeqDense(PETSC_COMM_SELF,lKZ+n,lKZ+n,NULL,&M);
541: BVMatProject(data->KZ,NULL,data->U,M);
542: MatDenseGetArrayRead(M,&array);
543: for (i=lKZ;i<lKZ+n;i++) { /* upper part */
544: PetscArraycpy(&data->XKZ[data->ldXKZ*i],&array[i*(lKZ+n)],lKZ);
545: }
546: for (i=0;i<lKZ+n;i++) { /* lower part */
547: PetscArraycpy(&data->XKZ[data->ldXKZ*i+lKZ],&array[i*(lKZ+n)+lKZ],n);
548: }
549: MatDenseRestoreArrayRead(M,&array);
550: MatDestroy(&M);
552: /* iXKZ <- inv(XKZ) */
553: size_KZ = lKZ+n;
554: PetscBLASIntCast(lKZ+n,&s);
555: data->ldiXKZ = data->size_iXKZ = size_KZ;
556: for (i=0;i<size_KZ;i++) PetscArraycpy(&data->iXKZ[data->ldiXKZ*i],&data->XKZ[data->ldXKZ*i],size_KZ);
557: PetscBLASIntCast(data->ldiXKZ,&ldXKZ);
558: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
559: PetscCallBLAS("LAPACKgetrf",LAPACKgetrf_(&s,&s,data->iXKZ,&ldXKZ,data->iXKZPivots,&info));
560: PetscFPTrapPop();
561: SlepcCheckLapackInfo("getrf",info);
562: return 0;
563: }
565: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
566: {
567: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
568: PetscInt i,j,n,maxits,maxits0,lits,s,ld,k,max_size_D,lV,kV;
569: PetscScalar *pX,*pY;
570: PetscReal tol,tol0;
571: Vec *kr,kr_comp,D_comp,D[2],kr0[2];
572: PetscBool odd_situation = PETSC_FALSE;
574: BVGetActiveColumns(d->eps->V,&lV,&kV);
575: max_size_D = d->eps->ncv-kV;
576: /* Quick exit */
577: if ((max_size_D == 0) || r_e-r_s <= 0) {
578: *size_D = 0;
579: return 0;
580: }
582: n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
583: PetscAssert(n>0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"n == 0");
584: PetscAssert(data->size_X>=r_e-r_s,PETSC_COMM_SELF,PETSC_ERR_PLIB,"size_X < r_e-r_s");
586: DSGetLeadingDimension(d->eps->ds,&ld);
588: /* Restart lastTol if a new pair converged */
589: if (data->dynamic && data->size_cX < lV)
590: data->lastTol = 0.5;
592: for (i=0,s=0;i<n;i+=s) {
593: /* If the selected eigenvalue is complex, but the arithmetic is real... */
594: #if !defined(PETSC_USE_COMPLEX)
595: if (d->eigi[r_s+i] != 0.0) {
596: if (i+2 <= max_size_D) s=2;
597: else break;
598: } else
599: #endif
600: s=1;
602: data->r_s = r_s+i;
603: data->r_e = r_s+i+s;
604: SlepcVecPoolGetVecs(d->auxV,s,&kr);
606: /* Compute theta, maximum iterations and tolerance */
607: maxits = 0;
608: tol = 1;
609: for (j=0;j<s;j++) {
610: d->improvex_jd_lit(d,r_s+i+j,&data->theta[2*j],&data->thetai[j],&maxits0,&tol0);
611: maxits += maxits0;
612: tol *= tol0;
613: }
614: maxits/= s;
615: tol = data->dynamic?data->lastTol:PetscExpReal(PetscLogReal(tol)/s);
617: /* Compute u, v and kr */
618: k = r_s+i;
619: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
620: k = r_s+i;
621: DSVectors(d->eps->ds,DS_MAT_Y,&k,NULL);
622: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
623: DSGetArray(d->eps->ds,DS_MAT_Y,&pY);
624: dvd_improvex_jd_proj_cuv(d,r_s+i,r_s+i+s,kr,data->theta,data->thetai,pX,pY,ld);
625: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
626: DSRestoreArray(d->eps->ds,DS_MAT_Y,&pY);
628: /* Check if the first eigenpairs are converged */
629: if (i == 0) {
630: PetscInt oldnpreconv = d->npreconv;
631: d->preTestConv(d,0,r_s+s,r_s+s,&d->npreconv);
632: if (d->npreconv > oldnpreconv) break;
633: }
635: /* Test the odd situation of solving Ax=b with A=I */
636: #if !defined(PETSC_USE_COMPLEX)
637: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && data->thetai[0] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
638: #else
639: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
640: #endif
641: /* If JD */
642: if (data->ksp && !odd_situation) {
643: /* kr <- -kr */
644: for (j=0;j<s;j++) VecScale(kr[j],-1.0);
646: /* Compose kr and D */
647: kr0[0] = kr[0];
648: kr0[1] = (s==2 ? kr[1] : NULL);
649: VecCreateCompWithVecs(kr0,data->ksp_max_size,data->friends,&kr_comp);
650: BVGetColumn(d->eps->V,kV+i,&D[0]);
651: if (s==2) BVGetColumn(d->eps->V,kV+i+1,&D[1]);
652: else D[1] = NULL;
653: VecCreateCompWithVecs(D,data->ksp_max_size,data->friends,&D_comp);
654: VecCompSetSubVecs(data->friends,s,NULL);
656: /* Solve the correction equation */
657: KSPSetTolerances(data->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,maxits);
658: KSPSolve(data->ksp,kr_comp,D_comp);
659: KSPGetIterationNumber(data->ksp,&lits);
661: /* Destroy the composed ks and D */
662: VecDestroy(&kr_comp);
663: VecDestroy(&D_comp);
664: BVRestoreColumn(d->eps->V,kV+i,&D[0]);
665: if (s==2) BVRestoreColumn(d->eps->V,kV+i+1,&D[1]);
667: /* If GD */
668: } else {
669: BVGetColumn(d->eps->V,kV+i,&D[0]);
670: if (s==2) BVGetColumn(d->eps->V,kV+i+1,&D[1]);
671: for (j=0;j<s;j++) d->improvex_precond(d,r_s+i+j,kr[j],D[j]);
672: dvd_improvex_apply_proj(d,D,s);
673: BVRestoreColumn(d->eps->V,kV+i,&D[0]);
674: if (s==2) BVRestoreColumn(d->eps->V,kV+i+1,&D[1]);
675: }
676: /* Prevent that short vectors are discarded in the orthogonalization */
677: if (i == 0 && d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
678: for (j=0;j<s;j++) BVScaleColumn(d->eps->V,kV+i+j,1.0/d->eps->errest[d->nconv+r_s]);
679: }
680: SlepcVecPoolRestoreVecs(d->auxV,s,&kr);
681: }
682: *size_D = i;
683: if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
684: return 0;
685: }
687: PetscErrorCode dvd_improvex_jd(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs,PetscBool dynamic)
688: {
689: dvdImprovex_jd *data;
690: PetscBool useGD;
691: PC pc;
692: PetscInt size_P;
694: /* Setting configuration constrains */
695: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&useGD);
697: /* If the arithmetic is real and the problem is not Hermitian, then
698: the block size is incremented in one */
699: #if !defined(PETSC_USE_COMPLEX)
700: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) {
701: max_bs++;
702: b->max_size_P = PetscMax(b->max_size_P,2);
703: } else
704: #endif
705: {
706: b->max_size_P = PetscMax(b->max_size_P,1);
707: }
708: b->max_size_X = PetscMax(b->max_size_X,max_bs);
709: size_P = b->max_size_P;
711: /* Setup the preconditioner */
712: if (ksp) {
713: KSPGetPC(ksp,&pc);
714: dvd_static_precond_PC(d,b,pc);
715: } else dvd_static_precond_PC(d,b,NULL);
717: /* Setup the step */
718: if (b->state >= DVD_STATE_CONF) {
719: PetscNew(&data);
720: data->dynamic = dynamic;
721: PetscMalloc1(size_P*size_P,&data->XKZ);
722: PetscMalloc1(size_P*size_P,&data->iXKZ);
723: PetscMalloc1(size_P,&data->iXKZPivots);
724: data->ldXKZ = size_P;
725: data->size_X = b->max_size_X;
726: d->improveX_data = data;
727: data->ksp = useGD? NULL: ksp;
728: data->d = d;
729: d->improveX = dvd_improvex_jd_gen;
730: #if !defined(PETSC_USE_COMPLEX)
731: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) data->ksp_max_size = 2;
732: else
733: #endif
734: data->ksp_max_size = 1;
735: /* Create various vector basis */
736: BVDuplicateResize(d->eps->V,size_P,&data->KZ);
737: BVSetMatrix(data->KZ,NULL,PETSC_FALSE);
738: BVDuplicate(data->KZ,&data->U);
740: EPSDavidsonFLAdd(&d->startList,dvd_improvex_jd_start);
741: EPSDavidsonFLAdd(&d->endList,dvd_improvex_jd_end);
742: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_jd_d);
743: }
744: return 0;
745: }
747: #if !defined(PETSC_USE_COMPLEX)
748: static inline PetscErrorCode dvd_complex_rayleigh_quotient(Vec ur,Vec ui,Vec Axr,Vec Axi,Vec Bxr,Vec Bxi,PetscScalar *eigr,PetscScalar *eigi)
749: {
750: PetscScalar rAr,iAr,rAi,iAi,rBr,iBr,rBi,iBi,b0,b2,b4,b6,b7;
752: /* eigr = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k
753: eigi = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k
754: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */
755: VecDotBegin(Axr,ur,&rAr); /* r*A*r */
756: VecDotBegin(Axr,ui,&iAr); /* i*A*r */
757: VecDotBegin(Axi,ur,&rAi); /* r*A*i */
758: VecDotBegin(Axi,ui,&iAi); /* i*A*i */
759: VecDotBegin(Bxr,ur,&rBr); /* r*B*r */
760: VecDotBegin(Bxr,ui,&iBr); /* i*B*r */
761: VecDotBegin(Bxi,ur,&rBi); /* r*B*i */
762: VecDotBegin(Bxi,ui,&iBi); /* i*B*i */
763: VecDotEnd(Axr,ur,&rAr); /* r*A*r */
764: VecDotEnd(Axr,ui,&iAr); /* i*A*r */
765: VecDotEnd(Axi,ur,&rAi); /* r*A*i */
766: VecDotEnd(Axi,ui,&iAi); /* i*A*i */
767: VecDotEnd(Bxr,ur,&rBr); /* r*B*r */
768: VecDotEnd(Bxr,ui,&iBr); /* i*B*r */
769: VecDotEnd(Bxi,ur,&rBi); /* r*B*i */
770: VecDotEnd(Bxi,ui,&iBi); /* i*B*i */
771: b0 = rAr+iAi; /* rAr+iAi */
772: b2 = rAi-iAr; /* rAi-iAr */
773: b4 = rBr+iBi; /* rBr+iBi */
774: b6 = rBi-iBr; /* rBi-iBr */
775: b7 = b4*b4 + b6*b6; /* k */
776: *eigr = (b0*b4 + b2*b6) / b7; /* eig_r */
777: *eigi = (b2*b4 - b0*b6) / b7; /* eig_i */
778: return 0;
779: }
780: #endif
782: static inline PetscErrorCode dvd_compute_n_rr(PetscInt i_s,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,Vec *u,Vec *Ax,Vec *Bx)
783: {
784: PetscInt i;
785: PetscScalar b0,b1;
787: for (i=0; i<n; i++) {
788: #if !defined(PETSC_USE_COMPLEX)
789: if (eigi[i_s+i] != 0.0) {
790: PetscScalar eigr0=0.0,eigi0=0.0;
791: dvd_complex_rayleigh_quotient(u[i],u[i+1],Ax[i],Ax[i+1],Bx[i],Bx[i+1],&eigr0,&eigi0);
792: if (PetscAbsScalar(eigr[i_s+i]-eigr0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10 || PetscAbsScalar(eigi[i_s+i]-eigi0)/PetscAbsScalar(eigi[i_s+i]) > 1e-10) PetscInfo(u[0],"The eigenvalue %g%+gi is far from its Rayleigh quotient value %g%+gi\n",(double)eigr[i_s+i],(double)eigi[i_s+i],(double)eigr0,(double)eigi0);
793: i++;
794: } else
795: #endif
796: {
797: VecDotBegin(Ax[i],u[i],&b0);
798: VecDotBegin(Bx[i],u[i],&b1);
799: VecDotEnd(Ax[i],u[i],&b0);
800: VecDotEnd(Bx[i],u[i],&b1);
801: b0 = b0/b1;
802: if (PetscAbsScalar(eigr[i_s+i]-b0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10) PetscInfo(u[0],"The eigenvalue %g+%g is far from its Rayleigh quotient value %g+%g\n",(double)PetscRealPart(eigr[i_s+i]),(double)PetscImaginaryPart(eigr[i_s+i]),(double)PetscRealPart(b0),(double)PetscImaginaryPart(b0));
803: }
804: }
805: return 0;
806: }
808: /*
809: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
810: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
811: where
812: pX,pY, the right and left eigenvectors of the projected system
813: ld, the leading dimension of pX and pY
814: */
815: PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,Vec *v,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
816: {
817: PetscInt n = i_e-i_s,i;
818: PetscScalar *b;
819: Vec *Ax,*Bx,*r;
820: Mat M;
821: BV X;
823: BVDuplicateResize(d->eps->V,4,&X);
824: MatCreateSeqDense(PETSC_COMM_SELF,4,4,NULL,&M);
825: /* u <- X(i) */
826: dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);
828: /* v <- theta[0]A*u + theta[1]*B*u */
830: /* Bx <- B*X(i) */
831: Bx = kr;
832: if (d->BX) {
833: for (i=i_s; i<i_e; ++i) BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
834: } else {
835: for (i=0;i<n;i++) {
836: if (d->B) MatMult(d->B, u[i], Bx[i]);
837: else VecCopy(u[i], Bx[i]);
838: }
839: }
841: /* Ax <- A*X(i) */
842: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
843: Ax = r;
844: for (i=i_s; i<i_e; ++i) BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
846: /* v <- Y(i) */
847: for (i=i_s; i<i_e; ++i) BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,v[i-i_s],&pY[ld*i]);
849: /* Recompute the eigenvalue */
850: dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,v,Ax,Bx);
852: for (i=0;i<n;i++) {
853: #if !defined(PETSC_USE_COMPLEX)
854: if (d->eigi[i_s+i] != 0.0) {
855: /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i' 0 1 0
856: 0 theta_2i' 0 1
857: theta_2i+1 -thetai_i -eigr_i -eigi_i
858: thetai_i theta_2i+1 eigi_i -eigr_i ] */
859: MatDenseGetArrayWrite(M,&b);
860: b[0] = b[5] = PetscConj(theta[2*i]);
861: b[2] = b[7] = -theta[2*i+1];
862: b[6] = -(b[3] = thetai[i]);
863: b[1] = b[4] = 0.0;
864: b[8] = b[13] = 1.0/d->nX[i_s+i];
865: b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
866: b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
867: b[9] = b[12] = 0.0;
868: MatDenseRestoreArrayWrite(M,&b);
869: BVInsertVec(X,0,Ax[i]);
870: BVInsertVec(X,1,Ax[i+1]);
871: BVInsertVec(X,2,Bx[i]);
872: BVInsertVec(X,3,Bx[i+1]);
873: BVSetActiveColumns(X,0,4);
874: BVMultInPlace(X,M,0,4);
875: BVCopyVec(X,0,Ax[i]);
876: BVCopyVec(X,1,Ax[i+1]);
877: BVCopyVec(X,2,Bx[i]);
878: BVCopyVec(X,3,Bx[i+1]);
879: i++;
880: } else
881: #endif
882: {
883: /* [Ax_i Bx_i]*= [ theta_2i' 1/nX_i
884: theta_2i+1 -eig_i/nX_i ] */
885: MatDenseGetArrayWrite(M,&b);
886: b[0] = PetscConj(theta[i*2]);
887: b[1] = theta[i*2+1];
888: b[4] = 1.0/d->nX[i_s+i];
889: b[5] = -d->eigr[i_s+i]/d->nX[i_s+i];
890: MatDenseRestoreArrayWrite(M,&b);
891: BVInsertVec(X,0,Ax[i]);
892: BVInsertVec(X,1,Bx[i]);
893: BVSetActiveColumns(X,0,2);
894: BVMultInPlace(X,M,0,2);
895: BVCopyVec(X,0,Ax[i]);
896: BVCopyVec(X,1,Bx[i]);
897: }
898: }
899: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
901: /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
902: for (i=0;i<n;i++) d->improvex_precond(d,i_s+i,r[i],v[i]);
903: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);
905: /* kr <- P*(Ax - eig_i*Bx) */
906: d->calcpairs_proj_res(d,i_s,i_e,kr);
907: BVDestroy(&X);
908: MatDestroy(&M);
909: return 0;
910: }
912: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol)
913: {
914: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
915: PetscReal a;
917: a = SlepcAbsEigenvalue(d->eigr[i],d->eigi[i]);
919: if (d->nR[i] < data->fix*a) {
920: theta[0] = d->eigr[i];
921: theta[1] = 1.0;
922: #if !defined(PETSC_USE_COMPLEX)
923: *thetai = d->eigi[i];
924: #endif
925: } else {
926: theta[0] = d->target[0];
927: theta[1] = d->target[1];
928: #if !defined(PETSC_USE_COMPLEX)
929: *thetai = 0.0;
930: #endif
931: }
933: #if defined(PETSC_USE_COMPLEX)
934: if (thetai) *thetai = 0.0;
935: #endif
936: *maxits = data->maxits;
937: *tol = data->tol;
938: return 0;
939: }
941: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d,dvdBlackboard *b,PetscInt maxits,PetscReal tol,PetscReal fix)
942: {
943: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
945: /* Setup the step */
946: if (b->state >= DVD_STATE_CONF) {
947: data->maxits = maxits;
948: data->tol = tol;
949: data->fix = fix;
950: d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
951: }
952: return 0;
953: }
955: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d,dvdBlackboard *b)
956: {
957: /* Setup the step */
958: if (b->state >= DVD_STATE_CONF) {
959: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX;
960: }
961: return 0;
962: }
964: PetscErrorCode dvd_improvex_compute_X(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u_,PetscScalar *pX,PetscInt ld)
965: {
966: PetscInt n = i_e - i_s,i;
967: Vec *u;
969: if (u_) u = u_;
970: else if (d->correctXnorm) SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&u);
971: if (u_ || d->correctXnorm) {
972: for (i=0;i<n;i++) BVMultVec(d->eps->V,1.0,0.0,u[i],&pX[ld*(i+i_s)]);
973: }
974: /* nX(i) <- ||X(i)|| */
975: if (d->correctXnorm) {
976: for (i=0;i<n;i++) VecNormBegin(u[i],NORM_2,&d->nX[i_s+i]);
977: for (i=0;i<n;i++) VecNormEnd(u[i],NORM_2,&d->nX[i_s+i]);
978: #if !defined(PETSC_USE_COMPLEX)
979: for (i=0;i<n;i++) {
980: if (d->eigi[i_s+i] != 0.0) {
981: d->nX[i_s+i] = d->nX[i_s+i+1] = PetscSqrtScalar(d->nX[i_s+i]*d->nX[i_s+i]+d->nX[i_s+i+1]*d->nX[i_s+i+1]);
982: i++;
983: }
984: }
985: #endif
986: } else {
987: for (i=0;i<n;i++) d->nX[i_s+i] = 1.0;
988: }
989: if (d->correctXnorm && !u_) SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&u);
990: return 0;
991: }