Actual source code: dvdimprovex.c
slepc-3.7.0 2016-05-16
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: improve the eigenvectors X
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
12: SLEPc is free software: you can redistribute it and/or modify it under the
13: terms of version 3 of the GNU Lesser General Public License as published by
14: the Free Software Foundation.
16: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
17: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
18: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
19: more details.
21: You should have received a copy of the GNU Lesser General Public License
22: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
23: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
24: */
26: #include davidson.h
27: #include <slepcblaslapack.h>
29: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r ****/
31: typedef struct {
32: PetscInt size_X;
33: KSP ksp; /* correction equation solver */
34: Vec friends; /* reference vector for composite vectors */
35: PetscScalar theta[4],thetai[2]; /* the shifts used in the correction eq. */
36: PetscInt maxits; /* maximum number of iterations */
37: PetscInt r_s,r_e; /* the selected eigenpairs to improve */
38: PetscInt ksp_max_size; /* the ksp maximum subvectors size */
39: PetscReal tol; /* the maximum solution tolerance */
40: PetscReal lastTol; /* last tol for dynamic stopping criterion */
41: PetscReal fix; /* tolerance for using the approx. eigenvalue */
42: PetscBool dynamic; /* if the dynamic stopping criterion is applied */
43: dvdDashboard *d; /* the currect dvdDashboard reference */
44: PC old_pc; /* old pc in ksp */
45: BV KZ; /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
46: BV U; /* new X vectors */
47: PetscScalar *XKZ; /* X'*KZ */
48: PetscScalar *iXKZ; /* inverse of XKZ */
49: PetscInt ldXKZ; /* leading dimension of XKZ */
50: PetscInt size_iXKZ; /* size of iXKZ */
51: PetscInt ldiXKZ; /* leading dimension of iXKZ */
52: PetscInt size_cX; /* last value of d->size_cX */
53: PetscInt old_size_X; /* last number of improved vectors */
54: PetscBLASInt *iXKZPivots; /* array of pivots */
55: } dvdImprovex_jd;
59: /*
60: Compute (I - KZ*iXKZ*X')*V where,
61: V, the vectors to apply the projector,
62: cV, the number of vectors in V,
63: */
64: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV)
65: {
66: #if defined(PETSC_MISSING_LAPACK_GETRS)
68: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
69: #else
71: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
72: PetscInt i,ldh,k,l;
73: PetscScalar *h;
74: PetscBLASInt cV_,n,info,ld;
75: #if defined(PETSC_USE_COMPLEX)
76: PetscInt j;
77: #endif
80: if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
82: /* h <- X'*V */
83: PetscMalloc1(data->size_iXKZ*cV,&h);
84: ldh = data->size_iXKZ;
85: BVGetActiveColumns(data->U,&l,&k);
86: if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
87: BVSetActiveColumns(data->U,0,k);
88: for (i=0;i<cV;i++) {
89: BVDotVec(data->U,V[i],&h[ldh*i]);
90: #if defined(PETSC_USE_COMPLEX)
91: for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
92: #endif
93: }
94: BVSetActiveColumns(data->U,l,k);
96: /* h <- iXKZ\h */
97: PetscBLASIntCast(cV,&cV_);
98: PetscBLASIntCast(data->size_iXKZ,&n);
99: PetscBLASIntCast(data->ldiXKZ,&ld);
101: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
102: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
103: PetscFPTrapPop();
104: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRS %d",info);
106: /* V <- V - KZ*h */
107: BVSetActiveColumns(data->KZ,0,k);
108: for (i=0;i<cV;i++) {
109: BVMultVec(data->KZ,-1.0,1.0,V[i],&h[ldh*i]);
110: }
111: BVSetActiveColumns(data->KZ,l,k);
112: PetscFree(h);
113: return(0);
114: #endif
115: }
119: /*
120: Compute (I - X*iXKZ*KZ')*V where,
121: V, the vectors to apply the projector,
122: cV, the number of vectors in V,
123: */
124: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV)
125: {
126: #if defined(PETSC_MISSING_LAPACK_GETRS)
128: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
129: return(0);
130: #else
132: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
133: PetscInt i,ldh,k,l;
134: PetscScalar *h;
135: PetscBLASInt cV_, n, info, ld;
136: #if defined(PETSC_USE_COMPLEX)
137: PetscInt j;
138: #endif
141: if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
143: /* h <- KZ'*V */
144: PetscMalloc1(data->size_iXKZ*cV,&h);
145: ldh = data->size_iXKZ;
146: BVGetActiveColumns(data->U,&l,&k);
147: if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
148: BVSetActiveColumns(data->KZ,0,k);
149: for (i=0;i<cV;i++) {
150: BVDotVec(data->KZ,V[i],&h[ldh*i]);
151: #if defined(PETSC_USE_COMPLEX)
152: for (j=0;j<k;j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
153: #endif
154: }
155: BVSetActiveColumns(data->KZ,l,k);
157: /* h <- iXKZ\h */
158: PetscBLASIntCast(cV,&cV_);
159: PetscBLASIntCast(data->size_iXKZ,&n);
160: PetscBLASIntCast(data->ldiXKZ,&ld);
162: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
163: PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
164: PetscFPTrapPop();
165: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRS %d",info);
167: /* V <- V - U*h */
168: BVSetActiveColumns(data->U,0,k);
169: for (i=0;i<cV;i++) {
170: BVMultVec(data->U,-1.0,1.0,V[i],&h[ldh*i]);
171: }
172: BVSetActiveColumns(data->U,l,k);
173: PetscFree(h);
174: return(0);
175: #endif
176: }
180: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
181: {
183: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
186: if (data->friends) { VecDestroy(&data->friends); }
188: /* Restore the pc of ksp */
189: if (data->old_pc) {
190: KSPSetPC(data->ksp, data->old_pc);
191: PCDestroy(&data->old_pc);
192: }
193: return(0);
194: }
198: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
199: {
201: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
204: /* Free local data and objects */
205: PetscFree(data->XKZ);
206: PetscFree(data->iXKZ);
207: PetscFree(data->iXKZPivots);
208: BVDestroy(&data->KZ);
209: BVDestroy(&data->U);
210: PetscFree(data);
211: return(0);
212: }
216: /*
217: y <- theta[1]A*x - theta[0]*B*x
218: auxV, two auxiliary vectors
219: */
220: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y)
221: {
223: PetscInt n,i;
224: const Vec *Bx;
225: Vec *auxV;
228: n = data->r_e - data->r_s;
229: for (i=0;i<n;i++) {
230: MatMult(data->d->A,x[i],y[i]);
231: }
233: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
234: for (i=0;i<n;i++) {
235: #if !defined(PETSC_USE_COMPLEX)
236: if (data->d->eigi[data->r_s+i] != 0.0) {
237: if (data->d->B) {
238: MatMult(data->d->B,x[i],auxV[0]);
239: MatMult(data->d->B,x[i+1],auxV[1]);
240: Bx = auxV;
241: } else Bx = &x[i];
243: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i + ti_i*Bx_i+1;
244: y_i+1 t_2i+1*A*x_i+1 - ti_i*Bx_i - t_2i*Bx_i+1 ] */
245: VecAXPBYPCZ(y[i],-data->theta[2*i],data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
246: VecAXPBYPCZ(y[i+1],-data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
247: i++;
248: } else
249: #endif
250: {
251: if (data->d->B) {
252: MatMult(data->d->B,x[i],auxV[0]);
253: Bx = auxV;
254: } else Bx = &x[i];
255: VecAXPBY(y[i],-data->theta[i*2],data->theta[i*2+1],Bx[0]);
256: }
257: }
258: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
259: return(0);
260: }
264: /*
265: y <- theta[1]'*A'*x - theta[0]'*B'*x
266: */
267: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y)
268: {
270: PetscInt n,i;
271: const Vec *Bx;
272: Vec *auxV;
275: n = data->r_e - data->r_s;
276: for (i=0;i<n;i++) {
277: MatMultTranspose(data->d->A,x[i],y[i]);
278: }
280: SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
281: for (i=0;i<n;i++) {
282: #if !defined(PETSC_USE_COMPLEX)
283: if (data->d->eigi[data->r_s+i] != 0.0) {
284: if (data->d->B) {
285: MatMultTranspose(data->d->B,x[i],auxV[0]);
286: MatMultTranspose(data->d->B,x[i+1],auxV[1]);
287: Bx = auxV;
288: } else Bx = &x[i];
290: /* y_i <- [ t_2i+1*A*x_i - t_2i*Bx_i - ti_i*Bx_i+1;
291: y_i+1 t_2i+1*A*x_i+1 + ti_i*Bx_i - t_2i*Bx_i+1 ] */
292: VecAXPBYPCZ(y[i],-data->theta[2*i],-data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
293: VecAXPBYPCZ(y[i+1],data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
294: i++;
295: } else
296: #endif
297: {
298: if (data->d->B) {
299: MatMultTranspose(data->d->B,x[i],auxV[0]);
300: Bx = auxV;
301: } else Bx = &x[i];
302: VecAXPBY(y[i],PetscConj(-data->theta[i*2]),PetscConj(data->theta[i*2+1]),Bx[0]);
303: }
304: }
305: SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
306: return(0);
307: }
311: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w)
312: {
314: dvdImprovex_jd *data;
315: PetscInt n,i;
316: const Vec *inx,*outx,*wx;
317: Vec *auxV;
318: Mat A;
321: PCGetOperators(pc,&A,NULL);
322: MatShellGetContext(A,(void**)&data);
323: VecCompGetSubVecs(in,NULL,&inx);
324: VecCompGetSubVecs(out,NULL,&outx);
325: VecCompGetSubVecs(w,NULL,&wx);
326: n = data->r_e - data->r_s;
327: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
328: switch (side) {
329: case PC_LEFT:
330: /* aux <- theta[1]A*in - theta[0]*B*in */
331: dvd_aux_matmult(data,inx,auxV);
333: /* out <- K * aux */
334: for (i=0;i<n;i++) {
335: data->d->improvex_precond(data->d,data->r_s+i,auxV[i],outx[i]);
336: }
337: break;
338: case PC_RIGHT:
339: /* aux <- K * in */
340: for (i=0;i<n;i++) {
341: data->d->improvex_precond(data->d,data->r_s+i,inx[i],auxV[i]);
342: }
344: /* out <- theta[1]A*auxV - theta[0]*B*auxV */
345: dvd_aux_matmult(data,auxV,outx);
346: break;
347: case PC_SYMMETRIC:
348: /* aux <- K^{1/2} * in */
349: for (i=0;i<n;i++) {
350: PCApplySymmetricRight(data->old_pc,inx[i],auxV[i]);
351: }
353: /* wx <- theta[1]A*auxV - theta[0]*B*auxV */
354: dvd_aux_matmult(data,auxV,wx);
356: /* aux <- K^{1/2} * in */
357: for (i=0;i<n;i++) {
358: PCApplySymmetricLeft(data->old_pc,wx[i],outx[i]);
359: }
360: break;
361: default:
362: SETERRQ(PETSC_COMM_SELF,1,"Unsupported KSP side");
363: }
364: /* out <- out - v*(u'*out) */
365: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
366: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
367: return(0);
368: }
372: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out)
373: {
375: dvdImprovex_jd *data;
376: PetscInt n,i;
377: const Vec *inx, *outx;
378: Mat A;
381: PCGetOperators(pc,&A,NULL);
382: MatShellGetContext(A,(void**)&data);
383: VecCompGetSubVecs(in,NULL,&inx);
384: VecCompGetSubVecs(out,NULL,&outx);
385: n = data->r_e - data->r_s;
386: /* out <- K * in */
387: for (i=0;i<n;i++) {
388: data->d->improvex_precond(data->d,data->r_s+i,inx[i],outx[i]);
389: }
390: /* out <- out - v*(u'*out) */
391: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
392: return(0);
393: }
397: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out)
398: {
400: dvdImprovex_jd *data;
401: PetscInt n,i;
402: const Vec *inx, *outx;
403: Vec *auxV;
404: Mat A;
407: PCGetOperators(pc,&A,NULL);
408: MatShellGetContext(A,(void**)&data);
409: VecCompGetSubVecs(in,NULL,&inx);
410: VecCompGetSubVecs(out,NULL,&outx);
411: n = data->r_e - data->r_s;
412: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
413: /* auxV <- in */
414: for (i=0;i<n;i++) {
415: VecCopy(inx[i],auxV[i]);
416: }
417: /* auxV <- auxV - u*(v'*auxV) */
418: dvd_improvex_applytrans_proj(data->d,auxV,n);
419: /* out <- K' * aux */
420: for (i=0;i<n;i++) {
421: PCApplyTranspose(data->old_pc,auxV[i],outx[i]);
422: }
423: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
424: return(0);
425: }
429: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out)
430: {
432: dvdImprovex_jd *data;
433: PetscInt n;
434: const Vec *inx, *outx;
435: PCSide side;
438: MatShellGetContext(A,(void**)&data);
439: VecCompGetSubVecs(in,NULL,&inx);
440: VecCompGetSubVecs(out,NULL,&outx);
441: n = data->r_e - data->r_s;
442: /* out <- theta[1]A*in - theta[0]*B*in */
443: dvd_aux_matmult(data,inx,outx);
444: KSPGetPCSide(data->ksp,&side);
445: if (side == PC_RIGHT) {
446: /* out <- out - v*(u'*out) */
447: dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
448: }
449: return(0);
450: }
454: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out)
455: {
457: dvdImprovex_jd *data;
458: PetscInt n,i;
459: const Vec *inx,*outx,*r;
460: Vec *auxV;
461: PCSide side;
464: MatShellGetContext(A,(void**)&data);
465: VecCompGetSubVecs(in,NULL,&inx);
466: VecCompGetSubVecs(out,NULL,&outx);
467: n = data->r_e - data->r_s;
468: KSPGetPCSide(data->ksp,&side);
469: if (side == PC_RIGHT) {
470: /* auxV <- in */
471: SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
472: for (i=0;i<n;i++) {
473: VecCopy(inx[i],auxV[i]);
474: }
475: /* auxV <- auxV - v*(u'*auxV) */
476: dvd_improvex_applytrans_proj(data->d,auxV,n);
477: r = auxV;
478: } else r = inx;
479: /* out <- theta[1]A*r - theta[0]*B*r */
480: dvd_aux_matmulttrans(data,r,outx);
481: if (side == PC_RIGHT) {
482: SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
483: }
484: return(0);
485: }
489: static PetscErrorCode MatCreateVecs_dvd_jd(Mat A,Vec *right,Vec *left)
490: {
492: Vec *r,*l;
493: dvdImprovex_jd *data;
494: PetscInt n,i;
497: MatShellGetContext(A,(void**)&data);
498: n = data->ksp_max_size;
499: if (right) {
500: PetscMalloc1(n,&r);
501: }
502: if (left) {
503: PetscMalloc1(n,&l);
504: }
505: for (i=0;i<n;i++) {
506: MatCreateVecs(data->d->A,right?&r[i]:NULL,left?&l[i]:NULL);
507: }
508: if (right) {
509: VecCreateCompWithVecs(r,n,data->friends,right);
510: for (i=0;i<n;i++) {
511: VecDestroy(&r[i]);
512: }
513: }
514: if (left) {
515: VecCreateCompWithVecs(l,n,data->friends,left);
516: for (i=0;i<n;i++) {
517: VecDestroy(&l[i]);
518: }
519: }
521: if (right) {
522: PetscFree(r);
523: }
524: if (left) {
525: PetscFree(l);
526: }
527: return(0);
528: }
532: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
533: {
535: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
536: PetscInt rA, cA, rlA, clA;
537: Mat A;
538: PetscBool t;
539: PC pc;
540: Vec v0[2];
543: data->size_cX = data->old_size_X = 0;
544: data->lastTol = data->dynamic?0.5:0.0;
546: /* Setup the ksp */
547: if (data->ksp) {
548: /* Create the reference vector */
549: BVGetColumn(d->eps->V,0,&v0[0]);
550: v0[1] = v0[0];
551: VecCreateCompWithVecs(v0,data->ksp_max_size,NULL,&data->friends);
552: BVRestoreColumn(d->eps->V,0,&v0[0]);
553: PetscLogObjectParent((PetscObject)d->eps,(PetscObject)data->friends);
555: /* Save the current pc and set a PCNONE */
556: KSPGetPC(data->ksp, &data->old_pc);
557: PetscObjectTypeCompare((PetscObject)data->old_pc,PCNONE,&t);
558: data->lastTol = 0.5;
559: if (t) data->old_pc = 0;
560: else {
561: PetscObjectReference((PetscObject)data->old_pc);
562: PCCreate(PetscObjectComm((PetscObject)d->eps),&pc);
563: PCSetType(pc,PCSHELL);
564: PCSetOperators(pc,d->A,d->A);
565: PCSetReusePreconditioner(pc,PETSC_TRUE);
566: PCShellSetApply(pc,PCApply_dvd);
567: PCShellSetApplyBA(pc,PCApplyBA_dvd);
568: PCShellSetApplyTranspose(pc,PCApplyTranspose_dvd);
569: KSPSetPC(data->ksp,pc);
570: PCDestroy(&pc);
571: }
573: /* Create the (I-v*u')*K*(A-s*B) matrix */
574: MatGetSize(d->A,&rA,&cA);
575: MatGetLocalSize(d->A,&rlA,&clA);
576: 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);
577: MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_dvd_jd);
578: MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_dvd_jd);
579: MatShellSetOperation(A,MATOP_GET_VECS,(void(*)(void))MatCreateVecs_dvd_jd);
581: /* Try to avoid KSPReset */
582: KSPGetOperatorsSet(data->ksp,&t,NULL);
583: if (t) {
584: Mat M;
585: PetscInt rM;
586: KSPGetOperators(data->ksp,&M,NULL);
587: MatGetSize(M,&rM,NULL);
588: if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
589: }
590: KSPSetOperators(data->ksp,A,A);
591: KSPSetReusePreconditioner(data->ksp,PETSC_TRUE);
592: KSPSetUp(data->ksp);
593: MatDestroy(&A);
594: } else {
595: data->old_pc = 0;
596: data->friends = NULL;
597: }
598: BVSetActiveColumns(data->KZ,0,0);
599: BVSetActiveColumns(data->U,0,0);
600: return(0);
601: }
605: /*
606: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
607: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
608: where
609: pX,pY, the right and left eigenvectors of the projected system
610: ld, the leading dimension of pX and pY
611: */
612: 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)
613: {
614: #if defined(PETSC_MISSING_LAPACK_GETRF)
616: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable");
617: #else
619: PetscInt n=i_e-i_s,size_KZ,V_new,rm,i,lv,kv,lKZ,kKZ;
620: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
621: PetscScalar *array;
622: Mat M;
623: Vec u[2],v[2];
624: PetscBLASInt s,ldXKZ,info;
627: /* Check consistency */
628: BVGetActiveColumns(d->eps->V,&lv,&kv);
629: V_new = lv - data->size_cX;
630: if (V_new > data->old_size_X) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
631: data->old_size_X = n;
632: data->size_cX = lv;
634: /* KZ <- KZ(rm:rm+max_cX-1) */
635: BVGetActiveColumns(data->KZ,&lKZ,&kKZ);
636: rm = PetscMax(V_new+lKZ-d->max_cX_in_impr,0);
637: if (rm > 0) {
638: for (i=0;i<lKZ;i++) {
639: BVCopyColumn(data->KZ,i+rm,i);
640: BVCopyColumn(data->U,i+rm,i);
641: }
642: }
644: /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
645: if (rm > 0) {
646: for (i=0;i<lKZ;i++) {
647: PetscMemcpy(&data->XKZ[i*data->ldXKZ+i],&data->XKZ[(i+rm)*data->ldXKZ+i+rm],sizeof(PetscScalar)*lKZ);
648: }
649: }
650: lKZ = PetscMin(d->max_cX_in_impr,lKZ+V_new);
651: BVSetActiveColumns(data->KZ,lKZ,lKZ+n);
652: BVSetActiveColumns(data->U,lKZ,lKZ+n);
654: /* Compute X, KZ and KR */
655: BVGetColumn(data->U,lKZ,u);
656: if (n>1) { BVGetColumn(data->U,lKZ+1,&u[1]); }
657: BVGetColumn(data->KZ,lKZ,v);
658: if (n>1) { BVGetColumn(data->KZ,lKZ+1,&v[1]); }
659: d->improvex_jd_proj_uv(d,i_s,i_e,u,v,kr,theta,thetai,pX,pY,ld);
660: BVRestoreColumn(data->U,lKZ,u);
661: if (n>1) { BVRestoreColumn(data->U,lKZ+1,&u[1]); }
662: BVRestoreColumn(data->KZ,lKZ,v);
663: if (n>1) { BVRestoreColumn(data->KZ,lKZ+1,&v[1]); }
665: /* XKZ <- U'*KZ */
666: MatCreateSeqDense(PETSC_COMM_SELF,lKZ+n,lKZ+n,NULL,&M);
667: BVMatProject(data->KZ,NULL,data->U,M);
668: MatDenseGetArray(M,&array);
669: for (i=lKZ;i<lKZ+n;i++) { /* upper part */
670: PetscMemcpy(&data->XKZ[data->ldXKZ*i],&array[i*(lKZ+n)],lKZ*sizeof(PetscScalar));
671: }
672: for (i=0;i<lKZ+n;i++) { /* lower part */
673: PetscMemcpy(&data->XKZ[data->ldXKZ*i+lKZ],&array[i*(lKZ+n)+lKZ],n*sizeof(PetscScalar));
674: }
675: MatDenseRestoreArray(M,&array);
676: MatDestroy(&M);
678: /* iXKZ <- inv(XKZ) */
679: size_KZ = lKZ+n;
680: PetscBLASIntCast(lKZ+n,&s);
681: data->ldiXKZ = data->size_iXKZ = size_KZ;
682: for (i=0;i<size_KZ;i++) {
683: PetscMemcpy(&data->iXKZ[data->ldiXKZ*i],&data->XKZ[data->ldXKZ*i],sizeof(PetscScalar)*size_KZ);
684: }
685: PetscBLASIntCast(data->ldiXKZ,&ldXKZ);
686: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
687: PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&s,&s,data->iXKZ,&ldXKZ,data->iXKZPivots,&info));
688: PetscFPTrapPop();
689: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRF %d",info);
690: return(0);
691: #endif
692: }
696: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
697: {
698: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
700: PetscInt i,j,n,maxits,maxits0,lits,s,ld,k,max_size_D,lV,kV;
701: PetscScalar *pX,*pY;
702: PetscReal tol,tol0;
703: Vec *kr,kr_comp,D_comp,D[2],kr0[2];
704: PetscBool odd_situation = PETSC_FALSE;
707: BVGetActiveColumns(d->eps->V,&lV,&kV);
708: max_size_D = d->eps->ncv-kV;
709: /* Quick exit */
710: if ((max_size_D == 0) || r_e-r_s <= 0) {
711: *size_D = 0;
712: return(0);
713: }
715: n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
716: if (n == 0) SETERRQ(PETSC_COMM_SELF,1,"n == 0");
717: if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1,"size_X < r_e-r_s");
719: DSGetLeadingDimension(d->eps->ds,&ld);
721: /* Restart lastTol if a new pair converged */
722: if (data->dynamic && data->size_cX < lV)
723: data->lastTol = 0.5;
725: for (i=0,s=0;i<n;i+=s) {
726: /* If the selected eigenvalue is complex, but the arithmetic is real... */
727: #if !defined(PETSC_USE_COMPLEX)
728: if (d->eigi[i] != 0.0) {
729: if (i+2 <= max_size_D) s=2;
730: else break;
731: } else
732: #endif
733: s=1;
735: data->r_s = r_s+i;
736: data->r_e = r_s+i+s;
737: SlepcVecPoolGetVecs(d->auxV,s,&kr);
739: /* Compute theta, maximum iterations and tolerance */
740: maxits = 0;
741: tol = 1;
742: for (j=0;j<s;j++) {
743: d->improvex_jd_lit(d,r_s+i+j,&data->theta[2*j],&data->thetai[j],&maxits0,&tol0);
744: maxits += maxits0;
745: tol *= tol0;
746: }
747: maxits/= s;
748: tol = data->dynamic?data->lastTol:PetscExpReal(PetscLogReal(tol)/s);
750: /* Compute u, v and kr */
751: k = r_s+i;
752: DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
753: k = r_s+i;
754: DSVectors(d->eps->ds,DS_MAT_Y,&k,NULL);
755: DSGetArray(d->eps->ds,DS_MAT_X,&pX);
756: DSGetArray(d->eps->ds,DS_MAT_Y,&pY);
757: dvd_improvex_jd_proj_cuv(d,r_s+i,r_s+i+s,kr,data->theta,data->thetai,pX,pY,ld);
758: DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
759: DSRestoreArray(d->eps->ds,DS_MAT_Y,&pY);
761: /* Check if the first eigenpairs are converged */
762: if (i == 0) {
763: d->preTestConv(d,0,s,s,&d->npreconv);
764: if (d->npreconv > 0) break;
765: }
767: /* Test the odd situation of solving Ax=b with A=I */
768: #if !defined(PETSC_USE_COMPLEX)
769: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && data->thetai[0] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
770: #else
771: odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
772: #endif
773: /* If JD */
774: if (data->ksp && !odd_situation) {
775: /* kr <- -kr */
776: for (j=0;j<s;j++) {
777: VecScale(kr[j],-1.0);
778: }
780: /* Compose kr and D */
781: kr0[0] = kr[0];
782: kr0[1] = (s==2 ? kr[1] : NULL);
783: VecCreateCompWithVecs(kr0,data->ksp_max_size,data->friends,&kr_comp);
784: BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
785: if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
786: else D[1] = NULL;
787: VecCreateCompWithVecs(D,data->ksp_max_size,data->friends,&D_comp);
788: VecCompSetSubVecs(data->friends,s,NULL);
790: /* Solve the correction equation */
791: KSPSetTolerances(data->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,maxits);
792: KSPSolve(data->ksp,kr_comp,D_comp);
793: KSPGetIterationNumber(data->ksp,&lits);
795: /* Destroy the composed ks and D */
796: VecDestroy(&kr_comp);
797: VecDestroy(&D_comp);
798: BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
799: if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
801: /* If GD */
802: } else {
803: BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
804: if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
805: for (j=0;j<s;j++) {
806: d->improvex_precond(d,r_s+i+j,kr[j],D[j]);
807: }
808: dvd_improvex_apply_proj(d,D,s);
809: BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
810: if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
811: }
812: /* Prevent that short vectors are discarded in the orthogonalization */
813: if (i == 0 && d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
814: for (j=0;j<s;j++) {
815: BVScaleColumn(d->eps->V,kV+r_s+i+j,1.0/d->eps->errest[d->nconv+r_s]);
816: }
817: }
818: SlepcVecPoolRestoreVecs(d->auxV,s,&kr);
819: }
820: *size_D = i;
821: if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
822: return(0);
823: }
827: PetscErrorCode dvd_improvex_jd(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs,PetscInt cX_impr,PetscBool dynamic)
828: {
830: dvdImprovex_jd *data;
831: PetscBool useGD;
832: PC pc;
833: PetscInt size_P;
836: /* Setting configuration constrains */
837: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&useGD);
839: /* If the arithmetic is real and the problem is not Hermitian, then
840: the block size is incremented in one */
841: #if !defined(PETSC_USE_COMPLEX)
842: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) {
843: max_bs++;
844: b->max_size_P = PetscMax(b->max_size_P,2);
845: } else
846: #endif
847: {
848: b->max_size_P = PetscMax(b->max_size_P,1);
849: }
850: b->max_size_X = PetscMax(b->max_size_X,max_bs);
851: size_P = b->max_size_P+cX_impr;
853: /* Setup the preconditioner */
854: if (ksp) {
855: KSPGetPC(ksp,&pc);
856: dvd_static_precond_PC(d,b,pc);
857: } else {
858: dvd_static_precond_PC(d,b,0);
859: }
861: /* Setup the step */
862: if (b->state >= DVD_STATE_CONF) {
863: PetscNewLog(d->eps,&data);
864: data->dynamic = dynamic;
865: d->max_cX_in_impr = cX_impr;
866: PetscMalloc1(size_P*size_P,&data->XKZ);
867: PetscMalloc1(size_P*size_P,&data->iXKZ);
868: PetscMalloc1(size_P,&data->iXKZPivots);
869: data->ldXKZ = size_P;
870: data->size_X = b->max_size_X;
871: d->improveX_data = data;
872: data->ksp = useGD? NULL: ksp;
873: data->d = d;
874: d->improveX = dvd_improvex_jd_gen;
875: #if !defined(PETSC_USE_COMPLEX)
876: if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) data->ksp_max_size = 2;
877: else
878: #endif
879: data->ksp_max_size = 1;
880: /* Create various vector basis */
881: BVDuplicateResize(d->eps->V,size_P,&data->KZ);
882: BVSetMatrix(data->KZ,NULL,PETSC_FALSE);
883: BVDuplicate(data->KZ,&data->U);
885: EPSDavidsonFLAdd(&d->startList,dvd_improvex_jd_start);
886: EPSDavidsonFLAdd(&d->endList,dvd_improvex_jd_end);
887: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_jd_d);
888: }
889: return(0);
890: }
892: #if !defined(PETSC_USE_COMPLEX)
895: PETSC_STATIC_INLINE PetscErrorCode dvd_complex_rayleigh_quotient(Vec ur,Vec ui,Vec Axr,Vec Axi,Vec Bxr,Vec Bxi,PetscScalar *eigr,PetscScalar *eigi)
896: {
898: PetscScalar rAr,iAr,rAi,iAi,rBr,iBr,rBi,iBi,b0,b2,b4,b6,b7;
901: /* eigr = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k
902: eigi = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k
903: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */
904: VecDotBegin(Axr,ur,&rAr); /* r*A*r */
905: VecDotBegin(Axr,ui,&iAr); /* i*A*r */
906: VecDotBegin(Axi,ur,&rAi); /* r*A*i */
907: VecDotBegin(Axi,ui,&iAi); /* i*A*i */
908: VecDotBegin(Bxr,ur,&rBr); /* r*B*r */
909: VecDotBegin(Bxr,ui,&iBr); /* i*B*r */
910: VecDotBegin(Bxi,ur,&rBi); /* r*B*i */
911: VecDotBegin(Bxi,ui,&iBi); /* i*B*i */
912: VecDotEnd(Axr,ur,&rAr); /* r*A*r */
913: VecDotEnd(Axr,ui,&iAr); /* i*A*r */
914: VecDotEnd(Axi,ur,&rAi); /* r*A*i */
915: VecDotEnd(Axi,ui,&iAi); /* i*A*i */
916: VecDotEnd(Bxr,ur,&rBr); /* r*B*r */
917: VecDotEnd(Bxr,ui,&iBr); /* i*B*r */
918: VecDotEnd(Bxi,ur,&rBi); /* r*B*i */
919: VecDotEnd(Bxi,ui,&iBi); /* i*B*i */
920: b0 = rAr+iAi; /* rAr+iAi */
921: b2 = rAi-iAr; /* rAi-iAr */
922: b4 = rBr+iBi; /* rBr+iBi */
923: b6 = rBi-iBr; /* rBi-iBr */
924: b7 = b4*b4 + b6*b6; /* k */
925: *eigr = (b0*b4 + b2*b6) / b7; /* eig_r */
926: *eigi = (b2*b4 - b0*b6) / b7; /* eig_i */
927: return(0);
928: }
929: #endif
933: PETSC_STATIC_INLINE PetscErrorCode dvd_compute_n_rr(PetscInt i_s,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,Vec *u,Vec *Ax,Vec *Bx)
934: {
936: PetscInt i;
937: PetscScalar b0,b1;
940: for (i=0; i<n; i++) {
941: #if !defined(PETSC_USE_COMPLEX)
942: if (eigi[i_s+i] != 0.0) {
943: PetscScalar eigr0=0.0,eigi0=0.0;
944: dvd_complex_rayleigh_quotient(u[i],u[i+1],Ax[i],Ax[i+1],Bx[i],Bx[i+1],&eigr0,&eigi0);
945: 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) {
946: PetscInfo4(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);
947: }
948: i++;
949: } else
950: #endif
951: {
952: VecDotBegin(Ax[i],u[i],&b0);
953: VecDotBegin(Bx[i],u[i],&b1);
954: VecDotEnd(Ax[i],u[i],&b0);
955: VecDotEnd(Bx[i],u[i],&b1);
956: b0 = b0/b1;
957: if (PetscAbsScalar(eigr[i_s+i]-b0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10) {
958: PetscInfo4(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));
959: }
960: }
961: }
962: return(0);
963: }
967: /*
968: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
969: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
970: where
971: pX,pY, the right and left eigenvectors of the projected system
972: ld, the leading dimension of pX and pY
973: */
974: static 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)
975: {
977: PetscInt n = i_e-i_s,i;
978: PetscScalar *b;
979: Vec *Ax,*Bx,*r;
980: Mat M;
981: BV X;
984: BVDuplicateResize(d->eps->V,4,&X);
985: MatCreateSeqDense(PETSC_COMM_SELF,4,4,NULL,&M);
986: /* u <- X(i) */
987: dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);
989: /* v <- theta[0]A*u + theta[1]*B*u */
991: /* Bx <- B*X(i) */
992: Bx = kr;
993: if (d->BX) {
994: for (i=i_s; i<i_e; ++i) {
995: BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
996: }
997: } else {
998: for (i=0;i<n;i++) {
999: if (d->B) {
1000: MatMult(d->B, u[i], Bx[i]);
1001: } else {
1002: VecCopy(u[i], Bx[i]);
1003: }
1004: }
1005: }
1007: /* Ax <- A*X(i) */
1008: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
1009: Ax = r;
1010: for (i=i_s; i<i_e; ++i) {
1011: BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
1012: }
1014: /* v <- Y(i) */
1015: for (i=i_s; i<i_e; ++i) {
1016: BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,v[i-i_s],&pY[ld*i]);
1017: }
1019: /* Recompute the eigenvalue */
1020: dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,v,Ax,Bx);
1022: for (i=0;i<n;i++) {
1023: #if !defined(PETSC_USE_COMPLEX)
1024: if (d->eigi[i_s+i] != 0.0) {
1025: /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i' 0 1 0
1026: 0 theta_2i' 0 1
1027: theta_2i+1 -thetai_i -eigr_i -eigi_i
1028: thetai_i theta_2i+1 eigi_i -eigr_i ] */
1029: MatDenseGetArray(M,&b);
1030: b[0] = b[5] = PetscConj(theta[2*i]);
1031: b[2] = b[7] = -theta[2*i+1];
1032: b[6] = -(b[3] = thetai[i]);
1033: b[1] = b[4] = 0.0;
1034: b[8] = b[13] = 1.0/d->nX[i_s+i];
1035: b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
1036: b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
1037: b[9] = b[12] = 0.0;
1038: MatDenseRestoreArray(M,&b);
1039: BVInsertVec(X,0,Ax[i]);
1040: BVInsertVec(X,1,Ax[i+1]);
1041: BVInsertVec(X,2,Bx[i]);
1042: BVInsertVec(X,3,Bx[i+1]);
1043: BVSetActiveColumns(X,0,4);
1044: BVMultInPlace(X,M,0,4);
1045: BVCopyVec(X,0,Ax[i]);
1046: BVCopyVec(X,1,Ax[i+1]);
1047: BVCopyVec(X,2,Bx[i]);
1048: BVCopyVec(X,3,Bx[i+1]);
1049: i++;
1050: } else
1051: #endif
1052: {
1053: /* [Ax_i Bx_i]*= [ theta_2i' 1/nX_i
1054: theta_2i+1 -eig_i/nX_i ] */
1055: MatDenseGetArray(M,&b);
1056: b[0] = PetscConj(theta[i*2]);
1057: b[1] = theta[i*2+1];
1058: b[4] = 1.0/d->nX[i_s+i];
1059: b[5] = -d->eigr[i_s+i]/d->nX[i_s+i];
1060: MatDenseRestoreArray(M,&b);
1061: BVInsertVec(X,0,Ax[i]);
1062: BVInsertVec(X,1,Bx[i]);
1063: BVSetActiveColumns(X,0,2);
1064: BVMultInPlace(X,M,0,2);
1065: BVCopyVec(X,0,Ax[i]);
1066: BVCopyVec(X,1,Bx[i]);
1067: }
1068: }
1069: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
1071: /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
1072: for (i=0;i<n;i++) {
1073: d->improvex_precond(d,i_s+i,r[i],v[i]);
1074: }
1075: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);
1077: /* kr <- P*(Ax - eig_i*Bx) */
1078: d->calcpairs_proj_res(d,i_s,i_e,kr);
1079: BVDestroy(&X);
1080: MatDestroy(&M);
1081: return(0);
1082: }
1086: /*
1087: Compute: u <- K^{-1}*X, v <- X,
1088: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
1089: where
1090: pX,pY, the right and left eigenvectors of the projected system
1091: ld, the leading dimension of pX and pY
1092: */
1093: static PetscErrorCode dvd_improvex_jd_proj_uv_KXX(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,Vec *v,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
1094: {
1096: PetscInt n = i_e - i_s,i;
1097: PetscScalar *b;
1098: Vec *Ax,*Bx,*r;
1099: Mat M;
1100: BV X;
1103: BVDuplicateResize(d->eps->V,4,&X);
1104: MatCreateSeqDense(PETSC_COMM_SELF,4,2,NULL,&M);
1105: /* [v u] <- X(i) Y(i) */
1106: dvd_improvex_compute_X(d,i_s,i_e,v,pX,ld);
1107: for (i=i_s; i<i_e; ++i) {
1108: BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,u[i-i_s],&pY[ld*i]);
1109: }
1111: /* Bx <- B*X(i) */
1112: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
1113: Bx = r;
1114: if (d->BX) {
1115: for (i=i_s; i<i_e; ++i) {
1116: BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
1117: }
1118: } else {
1119: if (d->B) {
1120: for (i=0;i<n;i++) {
1121: MatMult(d->B,v[i],Bx[i]);
1122: }
1123: } else Bx = v;
1124: }
1126: /* Ax <- A*X(i) */
1127: Ax = kr;
1128: for (i=i_s; i<i_e; ++i) {
1129: BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
1130: }
1132: /* Recompute the eigenvalue */
1133: dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,u,Ax,Bx);
1135: for (i=0;i<n;i++) {
1136: if (d->eigi[i_s+i] == 0.0) {
1137: /* kr <- Ax -eig*Bx */
1138: VecAXPBY(kr[i],-d->eigr[i_s+i]/d->nX[i_s+i],1.0/d->nX[i_s+i],Bx[i]);
1139: } else {
1140: /* [kr_i kr_i+1 r_i r_i+1]*= [ 1 0
1141: 0 1
1142: -eigr_i -eigi_i
1143: eigi_i -eigr_i] */
1144: MatDenseGetArray(M,&b);
1145: b[0] = b[5] = 1.0/d->nX[i_s+i];
1146: b[2] = b[7] = -d->eigr[i_s+i]/d->nX[i_s+i];
1147: b[6] = -(b[3] = d->eigi[i_s+i]/d->nX[i_s+i]);
1148: b[1] = b[4] = 0.0;
1149: MatDenseRestoreArray(M,&b);
1150: BVInsertVec(X,0,kr[i]);
1151: BVInsertVec(X,1,kr[i+1]);
1152: BVInsertVec(X,2,r[i]);
1153: BVInsertVec(X,3,r[i+1]);
1154: BVSetActiveColumns(X,0,4);
1155: BVMultInPlace(X,M,0,2);
1156: BVCopyVec(X,0,kr[i]);
1157: BVCopyVec(X,1,kr[i+1]);
1158: i++;
1159: }
1160: }
1161: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
1163: /* kr <- P*kr */
1164: d->calcpairs_proj_res(d,i_s,i_e,r);
1165: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);
1167: /* u <- K^{-1} X(i) */
1168: for (i=0;i<n;i++) {
1169: d->improvex_precond(d,i_s+i,v[i],u[i]);
1170: }
1171: return(0);
1172: }
1176: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol)
1177: {
1178: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1179: PetscReal a;
1182: a = SlepcAbsEigenvalue(d->eigr[i],d->eigi[i]);
1184: if (d->nR[i]/a < data->fix) {
1185: theta[0] = d->eigr[i];
1186: theta[1] = 1.0;
1187: #if !defined(PETSC_USE_COMPLEX)
1188: *thetai = d->eigi[i];
1189: #endif
1190: } else {
1191: theta[0] = d->target[0];
1192: theta[1] = d->target[1];
1193: #if !defined(PETSC_USE_COMPLEX)
1194: *thetai = 0.0;
1195: #endif
1196: }
1198: #if defined(PETSC_USE_COMPLEX)
1199: if (thetai) *thetai = 0.0;
1200: #endif
1201: *maxits = data->maxits;
1202: *tol = data->tol;
1203: return(0);
1204: }
1208: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d,dvdBlackboard *b,PetscInt maxits,PetscReal tol,PetscReal fix)
1209: {
1210: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1213: /* Setup the step */
1214: if (b->state >= DVD_STATE_CONF) {
1215: data->maxits = maxits;
1216: data->tol = tol;
1217: data->fix = fix;
1218: d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
1219: }
1220: return(0);
1221: }
1225: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d,dvdBlackboard *b,ProjType_t p)
1226: {
1228: /* Setup the step */
1229: if (b->state >= DVD_STATE_CONF) {
1230: switch (p) {
1231: case DVD_PROJ_KXX:
1232: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KXX;
1233: break;
1234: case DVD_PROJ_KZX:
1235: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX;
1236: break;
1237: }
1238: }
1239: return(0);
1240: }
1244: PetscErrorCode dvd_improvex_compute_X(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u_,PetscScalar *pX,PetscInt ld)
1245: {
1247: PetscInt n = i_e - i_s,i;
1248: Vec *u;
1251: if (u_) u = u_;
1252: else if (d->correctXnorm) {
1253: SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&u);
1254: }
1255: if (u_ || d->correctXnorm) {
1256: for (i=0; i<n; i++) {
1257: BVMultVec(d->eps->V,1.0,0.0,u[i],&pX[ld*(i+i_s)]);
1258: }
1259: }
1260: /* nX(i) <- ||X(i)|| */
1261: if (d->correctXnorm) {
1262: for (i=0; i<n; i++) {
1263: VecNormBegin(u[i],NORM_2,&d->nX[i_s+i]);
1264: }
1265: for (i=0; i<n; i++) {
1266: VecNormEnd(u[i],NORM_2,&d->nX[i_s+i]);
1267: }
1268: #if !defined(PETSC_USE_COMPLEX)
1269: for (i=0;i<n;i++) {
1270: if (d->eigi[i_s+i] != 0.0) {
1271: 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]);
1272: i++;
1273: }
1274: }
1275: #endif
1276: } else {
1277: for (i=0;i<n;i++) d->nX[i_s+i] = 1.0;
1278: }
1279: if (d->correctXnorm && !u_) {
1280: SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&u);
1281: }
1282: return(0);
1283: }