Actual source code: dvd_improvex.c
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: improve the eigenvectors X
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2010, Universidad Politecnica de Valencia, Spain
10: This file is part of SLEPc.
11:
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
28: #include slepcvec.h
30: PetscErrorCode dvd_improvex_PfuncV(dvdDashboard *d, void *funcV, Vec *D,
31: PetscInt max_size_D, PetscInt r_s, PetscInt r_e, Vec *auxV,
32: PetscScalar *auxS);
33: PetscErrorCode dvd_matmult_jd(Mat A, Vec in, Vec out);
34: PetscErrorCode dvd_matgetvecs_jd(Mat A, Vec *right, Vec *left);
35: PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d);
36: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
37: PetscInt max_size_D, PetscInt r_s,
38: PetscInt r_e, PetscInt *size_D);
39: PetscErrorCode dvd_improvex_jd_proj_uv_KBXX(dvdDashboard *d, PetscInt i_s,
40: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
41: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
42: PetscErrorCode dvd_improvex_jd_proj_uv_KBXY(dvdDashboard *d, PetscInt i_s,
43: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
44: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
45: PetscErrorCode dvd_improvex_jd_proj_uv_KBXZ(dvdDashboard *d, PetscInt i_s,
46: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
47: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
48: PetscErrorCode dvd_improvex_jd_proj_uv_KBXZY(dvdDashboard *d, PetscInt i_s,
49: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
50: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
51: PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d, PetscInt i,
52: PetscScalar* theta, PetscScalar* thetai,
53: PetscInt *maxits, PetscReal *tol);
54: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
55: PetscScalar *pY, PetscInt ld_,
56: PetscScalar *auxS, PetscInt size_auxS);
58: #define size_Z (64*4)
60: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r ****/
62: typedef struct {
63: PetscInt size_X;
64: void
65: *old_improveX_data; /* old improveX_data */
66: improveX_type
67: old_improveX; /* old improveX */
68: KSP ksp; /* correction equation solver */
69: Vec *u,*v, /* the projector is (I-v*u') */
70: friends, /* reference vector for composite vectors */
71: *auxV; /* auxiliar vectors */
72: PetscScalar *theta,
73: *thetai; /* the shifts used in the correction eq. */
74: PetscInt maxits, /* maximum number of iterations */
75: r_s, r_e, /* the selected eigenpairs to improve */
76: n_uv, /* number of vectors u, v and kr */
77: ksp_max_size; /* the ksp maximum subvectors size */
78: PetscReal tol, /* the maximum solution tolerance */
79: fix; /* tolerance for using the approx. eigenvalue */
80: dvdDashboard
81: *d; /* the currect dvdDashboard reference */
82: } dvdImprovex_jd;
86: PetscErrorCode dvd_improvex_jd(dvdDashboard *d, dvdBlackboard *b, KSP ksp,
87: PetscInt max_bs)
88: {
89: PetscErrorCode ierr;
90: dvdImprovex_jd *data;
91: PetscInt rA, cA, rlA, clA;
92: Mat A;
93: PetscTruth t;
97: /* Setting configuration constrains */
98: /* If the arithmetic is real and the problem is not Hermitian, then
99: the block size is incremented in one */
100: #ifndef PETSC_USE_COMPLEX
101: if (DVD_ISNOT(d->sEP, DVD_EP_HERMITIAN)) {
102: max_bs++;
103: b->max_size_X = PetscMax(b->max_size_X, max_bs);
104: }
105: #endif
106: b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X*4); /* u,v,kr,
107: auxV */
108: b->max_size_auxS = PetscMax(b->max_size_auxS,
109: b->max_size_X*3 + /* theta, thetai */
110: 2*b->max_size_V*b->max_size_V + /* pX, pY */
111: 11*b->max_size_V+4*b->max_size_V*b->max_size_V
112: /* dvd_improvex_get_eigenvectors */
113: );
115: /* Setup the step */
116: if (b->state >= DVD_STATE_CONF) {
117: PetscMalloc(sizeof(dvdImprovex_jd), &data);
119: data->size_X = b->max_size_X;
120: data->old_improveX_data = d->improveX_data;
121: d->improveX_data = data;
122: data->old_improveX = d->improveX;
123: PetscTypeCompare((PetscObject)ksp, KSPPREONLY, &t);
124: data->ksp = t?0:ksp;
125: data->d = d;
126: d->improveX = dvd_improvex_jd_gen;
128: /* Create the reference vector */
129: VecCreateCompWithVecs(d->V, max_bs, PETSC_NULL, &data->friends);
130:
132: /* Setup the ksp */
133: if(data->ksp) {
134: /* Create the (I-v*u')*K*(A-s*B) matrix */
135: MatGetSize(d->A, &rA, &cA);
136: MatGetLocalSize(d->A, &rlA, &clA);
137: MatCreateShell(((PetscObject)d->A)->comm, rlA*max_bs, clA*max_bs,
138: rA*max_bs, cA*max_bs, data, &A);
139: MatShellSetOperation(A, MATOP_MULT,
140: (void(*)(void))dvd_matmult_jd);
141: MatShellSetOperation(A, MATOP_GET_VECS,
142: (void(*)(void))dvd_matgetvecs_jd);
143:
144:
145: data->ksp_max_size = max_bs;
146: KSPSetOperators(data->ksp, A, A, SAME_PRECONDITIONER);
147:
148: KSPSetUp(data->ksp);
149: MatDestroy(A);
150: }
152: DVD_FL_ADD(d->destroyList, dvd_improvex_jd_d);
153: }
155: return(0);
156: }
161: PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
162: {
163: PetscErrorCode ierr;
164: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
167:
168: /* Restore changes in dvdDashboard */
169: d->improveX_data = data->old_improveX_data;
171: /* Free local data and objects */
172: VecDestroy(data->friends);
173: PetscFree(data);
175: return(0);
176: }
181: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
182: PetscInt max_size_D, PetscInt r_s,
183: PetscInt r_e, PetscInt *size_D)
184: {
185: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
186: PetscErrorCode ierr;
187: PetscInt i, j, n, maxits, maxits0, lits, s;
188: PetscScalar a, *pX, *pY, *auxS = d->auxS;
189: PetscReal tol, tol0;
190: Vec *u, *v, *kr, kr_comp, D_comp;
194: /* Quick exit */
195: if ((max_size_D == 0) || r_e-r_s <= 0) {
196: /* Callback old improveX */
197: if (data->old_improveX) {
198: d->improveX_data = data->old_improveX_data;
199: data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
200: d->improveX_data = data;
201: }
202: return(0);
203: }
204:
205: n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
206: if (n == 0) SETERRQ(1, "n == 0!\n");
207: if (data->size_X < r_e-r_s) SETERRQ(1, "size_X < r_e-r_s!\n");
209: /* Compute the eigenvectors of the selected pairs */
210: pX = auxS; auxS+= d->size_H*d->size_H;
211: pY = auxS; auxS+= d->size_H*d->size_H;
212: dvd_improvex_get_eigenvectors(d, pX, pY, d->size_H, auxS,
213: d->size_auxS-(auxS-d->auxS));
214:
216: for(i=0, s=0; i<n; i+=s) {
217: /* If the selected eigenvalue is complex, but the arithmetic is real... */
218: #ifndef PETSC_USE_COMPLEX
219: if (PetscAbsScalar(d->eigi[i] != 0.0)) {
220: if (i+2 <= max_size_D) s=2; else break;
221: } else
222: #endif
223: s=1;
225: data->auxV = d->auxV;
226: data->r_s = r_s+i; data->r_e = r_s+i+s;
227: data->n_uv = s;
228: data->theta = auxS; data->thetai = auxS+2*s;
230: /* Compute theta, maximum iterations and tolerance */
231: maxits = 0; tol = 1;
232: for(j=0; j<s; j++) {
233: d->improvex_jd_lit(d, r_s+i+j, &data->theta[2*j],
234: &data->thetai[j], &maxits0, &tol0);
235:
236: maxits+= maxits0; tol*= tol0;
237: }
238: maxits/= s; tol = exp(log(tol)/s);
240: /* Compute u, v and kr */
241: d->improvex_jd_proj_uv(d, r_s+i, r_s+i+s, &u, &v, &kr,
242: &data->auxV, data->theta, data->thetai,
243: &pX[d->size_H*(r_s+i)],
244: &pY[d->size_H*(r_s+i)], d->size_H);
245:
246: data->u = u; data->v = v;
248: /* Check if the first eigenpairs are converged */
249: if (i == 0) {
250: d->preTestConv(d, 0, s, s, PETSC_NULL, PETSC_NULL, &d->npreconv);
251:
252: if (d->npreconv > 0) break;
253: }
255: /* Compute kr <- kr - u*(v'*kr) */
256: for(j=0; j<s; j++) {
257: VecDot(kr[j], v[j], &a);
258: VecAXPY(kr[j], -a, u[j]);
259: VecScale(kr[j], -1.0);
260: }
262: /* If KSP==0, D <- kr */
263: if (!data->ksp) {
264: for(j=0; j<s; j++) {
265: VecCopy(kr[j], D[j+i]);
266: }
267: } else {
268: /* Compouse kr and D */
269: VecCreateCompWithVecs(kr, data->ksp_max_size, data->friends,
270: &kr_comp);
271: VecCreateCompWithVecs(&D[i], data->ksp_max_size, data->friends,
272: &D_comp);
273: VecCompSetVecs(data->friends, PETSC_NULL, s);
274:
275: /* Solve the correction equation */
276: KSPSetTolerances(data->ksp, tol, PETSC_DEFAULT, PETSC_DEFAULT,
277: maxits);
278: KSPSolve(data->ksp, kr_comp, D_comp);
279: KSPGetIterationNumber(data->ksp, &lits);
280: d->eps->OP->lineariterations+= lits;
281:
282: /* Destroy the composed ks and D */
283: VecDestroy(kr_comp);
284: VecDestroy(D_comp);
285: }
286: }
287: *size_D = i;
288:
289: /* Callback old improveX */
290: if (data->old_improveX) {
291: d->improveX_data = data->old_improveX_data;
292: data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
293: d->improveX_data = data;
294: }
296: return(0);
297: }
302: PetscErrorCode dvd_matmult_jd(Mat A, Vec in, Vec out)
303: {
304: PetscErrorCode ierr;
305: dvdImprovex_jd *data;
306: PetscScalar a;
307: PetscInt n, i;
308: const Vec *inx, *outx, *Bx;
312: MatShellGetContext(A, (void**)&data);
313: VecCompGetVecs(in, &inx, PETSC_NULL);
314: VecCompGetVecs(out, &outx, PETSC_NULL);
315: n = data->r_e - data->r_s;
317: /* aux <- theta[1]A*in - theta[0]*B*in */
318: for(i=0; i<n; i++) {
319: MatMult(data->d->A, inx[i], data->auxV[i]);
320: }
321: if (data->d->B) {
322: for(i=0; i<n; i++) {
323: MatMult(data->d->B, inx[i], outx[i]);
324: }
325: Bx = outx;
326: } else
327: Bx = inx;
329: for(i=0; i<n; i++) {
330: #ifndef PETSC_USE_COMPLEX
331: if(data->d->eigi[data->r_s+i] != 0.0) {
332: /* aux_i <- [ t_2i+1*A*inx_i - t_2i*Bx_i + ti_i*Bx_i+1;
333: aux_i+1 t_2i+1*A*inx_i+1 - ti_i*Bx_i - t_2i*Bx_i+1 ] */
334: VecAXPBYPCZ(data->auxV[i], -data->theta[2*i], data->thetai[i],
335: data->theta[2*i+1], Bx[i], Bx[i+1]);
336:
337: VecAXPBYPCZ(data->auxV[i+1], -data->thetai[i],
338: -data->theta[2*i], data->theta[2*i+1], Bx[i],
339: Bx[i+1]);
340: i++;
341: } else
342: #endif
343: {
344: VecAXPBY(data->auxV[i], -data->theta[i*2], data->theta[i*2+1],
345: Bx[i]);
346: }
347: }
349: /* out <- K * aux */
350: for(i=0; i<n; i++) {
351: data->d->improvex_precond(data->d, data->r_s+i, data->auxV[i],
352: outx[i]);
353: }
355: /* out <- out - u*(v'*out) */
356: for(i=0; i<n; i++) {
357: VecDot(outx[i], data->v[i], &a);
358: VecAXPY(outx[i], -a, data->u[i]);
359: }
361: return(0);
362: }
367: PetscErrorCode dvd_matgetvecs_jd(Mat A, Vec *right, Vec *left)
368: {
369: PetscErrorCode ierr;
370: Vec *r, *l;
371: dvdImprovex_jd *data;
372: PetscInt n, i;
376: MatShellGetContext(A, (void**)&data);
377: n = data->ksp_max_size;
378: if (right) {
379: PetscMalloc(sizeof(Vec)*n, &r);
380: }
381: if (left) {
382: PetscMalloc(sizeof(Vec)*n, &l);
383: }
384: for (i=0; i<n; i++) {
385: MatGetVecs(data->d->A, right?&r[i]:PETSC_NULL,
386: left?&l[i]:PETSC_NULL);
387: }
388: if(right) {
389: VecCreateCompWithVecs(r, n, data->friends, right);
390: for (i=0; i<n; i++) {
391: VecDestroy(r[i]);
392: }
393: }
394: if(left) {
395: VecCreateCompWithVecs(l, n, data->friends, left);
396: for (i=0; i<n; i++) {
397: VecDestroy(l[i]);
398: }
399: }
401: if (right) {
402: PetscFree(r);
403: }
404: if (left) {
405: PetscFree(l);
406: }
408: return(0);
409: }
414: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d, dvdBlackboard *b,
415: ProjType_t p)
416: {
419: /* Setup the step */
420: if (b->state >= DVD_STATE_CONF) {
421: switch(p) {
422: case DVD_PROJ_KBXX:
423: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KBXX; break;
424: case DVD_PROJ_KBXY:
425: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KBXY; break;
426: case DVD_PROJ_KBXZ:
427: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KBXZ; break;
428: case DVD_PROJ_KBXZY:
429: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KBXZY; break;
430: }
431: }
433: return(0);
434: }
436: #define DVD_COMPLEX_RAYLEIGH_QUOTIENT(ur,ui,Axr,Axi,Bxr,Bxi,eigr,eigi,b,ierr)\
437: { \
438: VecDot((Axr), (ur), &(b)[0]); /* r*A*r */ \
439: VecDot((Axr), (ui), &(b)[1]); /* i*A*r */ \
440: VecDot((Axi), (ur), &(b)[2]); /* r*A*i */ \
441: VecDot((Axi), (ui), &(b)[3]); /* i*A*i */ \
442: VecDot((Bxr), (ur), &(b)[4]); /* r*B*r */ \
443: VecDot((Bxr), (ui), &(b)[5]); /* i*B*r */ \
444: VecDot((Bxi), (ur), &(b)[6]); /* r*B*i */ \
445: VecDot((Bxi), (ui), &(b)[7]); /* i*B*i */ \
446: (b)[0] = (b)[0]+(b)[3]; /* rAr+iAi */ \
447: (b)[2] = (b)[2]-(b)[1]; /* rAi-iAr */ \
448: (b)[4] = (b)[4]+(b)[7]; /* rBr+iBi */ \
449: (b)[6] = (b)[6]-(b)[5]; /* rBi-iBr */ \
450: (b)[7] = (b)[4]*(b)[4] + (b)[6]*(b)[6]; /* k */ \
451: *(eigr) = ((b)[0]*(b)[4] + (b)[2]*(b)[6]) / (b)[7]; /* eig_r */ \
452: *(eigi) = ((b)[2]*(b)[4] - (b)[0]*(b)[6]) / (b)[7]; /* eig_i */ \
453: }
455: #if !defined(PETSC_USE_COMPLEX)
456: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
457: for((i)=0; (i)<(n); (i)++) { \
458: if ((eigi)[(i_s)+(i)] != 0.0) { \
459: /* eig_r = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k \
460: eig_i = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k \
461: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */ \
462: DVD_COMPLEX_RAYLEIGH_QUOTIENT((u)[(i)], (u)[(i)+1], (Ax)[(i)], \
463: (Ax)[(i)+1], (Bx)[(i)], (Bx)[(i)+1], &(b)[8], &(b)[9], (b), (ierr)); \
464: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[8])/ \
465: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10 || \
466: PetscAbsScalar((eigi)[(i_s)+(i)] - (b)[9])/ \
467: PetscAbsScalar((eigi)[(i_s)+(i)]) > 1e-10 ) { \
468: (ierr) = PetscInfo4((eps), "The eigenvalue %g+%g is far from its "\
469: "Rayleigh quotient value %g+%g\n", \
470: (eigr)[(i_s)+(i)], \
471: (eigi)[(i_s)+1], (b)[8], (b)[9]); \
472: } \
473: (i)++; \
474: } \
475: }
476: #else
477: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
478: for((i)=0; (i)<(n); (i)++) { \
479: (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]); \
480: (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]); \
481: (b)[0] = (b)[0]/(b)[1]; \
482: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
483: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10 ) { \
484: (ierr) = PetscInfo4((eps), "The eigenvalue %g+%g is far from its " \
485: "Rayleigh quotient value %g+%g\n", \
486: PetscRealPart((eigr)[(i_s)+(i)]), \
487: PetscImaginaryPart((eigr)[(i_s)+(i)]), PetscRealPart((b)[0]), \
488: PetscImaginaryPart((b)[0])); \
489: } \
490: }
491: #endif
493: #if !defined(PETSC_USE_COMPLEX)
494: #define DVD_NORM_FOR_UV(x) PetscAbsScalar(x)
495: #else
496: #define DVD_NORM_FOR_UV(x) (x)
497: #endif
499: #define DVD_NORMALIZE_UV(u,v,ierr,a) { \
500: (ierr) = VecDot((u), (v), &(a)); \
501: if ((a) == 0.0) { \
502: SETERRQ(1, "Error: inappropriate approximate eigenvector norm!"); \
503: } \
504: if ((u) == (v)) { \
505: VecScale((u), 1.0/PetscSqrtScalar(DVD_NORM_FOR_UV(a))); \
506: \
507: } else { \
508: VecScale((u), 1.0/(a)); \
509: } \
510: }
513: /*
514: Compute: u <- K^{-1}*B*X, v <- (theta[0]*A+theta[1]*B)*X,
515: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1],
516: where
517: auxV, 4*(i_e-i_s) auxiliar global vectors
518: pX,pY, the right and left eigenvectors of the projected system
519: ld, the leading dimension of pX and pY
520: */
523: PetscErrorCode dvd_improvex_jd_proj_uv_KBXZ(dvdDashboard *d, PetscInt i_s,
524: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
525: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
526: {
527: PetscErrorCode ierr;
528: PetscInt n = i_e - i_s, i;
529: PetscScalar a, b[16];
530: Vec *Ax, *Bx, *r, *auxV = *auxV_, X[4];
531: /* The memory manager doen't allow to call a subroutines */
532: PetscScalar Z[size_Z];
536: /* Book space for u, v and kr */
537: *u = auxV; auxV+= n;
538: *kr = auxV; auxV+= n;
539: *v = auxV; auxV+= n;
541: /* Ax = v, Bx = r = auxV */
542: Ax = *v;
543: r = Bx = auxV; auxV+= n;
545: /* Ax <- A*X(i) */
546: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV, d->size_AV, pX, ld,
547: d->size_H, n);
549: /* Bx <- B*X(i) */
550: for(i=i_s; i<i_e; i++) d->nX[i] = 1.0;
551: if (d->BV) {
552: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV, d->size_BV, pX, ld,
553: d->size_H, n);
554: } else {
555: SlepcUpdateVectorsZ(d->B?*u:Bx, 0.0, 1.0, d->V, d->size_V, pX, ld,
556: d->size_H, n);
557: if (d->B) {
558: for(i=0; i<n; i++) {
559: MatMult(d->B, (*u)[i], Bx[i]);
560: }
561: }
562: }
564: /* Recompute the eigenvalue */
565: SlepcUpdateVectorsZ(*u, 0.0, 1.0, d->W?d->W:d->V, d->size_V, pY, ld,
566: d->size_H, n);
567: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, *u, Ax, Bx, b, ierr);
569: /* u <- K^{-1} Bx */
570: for(i=0; i<n; i++) {
571: d->improvex_precond(d, i_s+i, Bx[i], (*u)[i]);
572: }
574: for(i=0; i<n; i++) {
575: if (d->eigi[i_s+i] == 0.0) {
576: /* [v r] <- [Ax Bx][theta_2i' 1 ]
577: [theta_2i+1 -eig] */
578: b[0] = PetscConj(theta[i*2]); b[1] = theta[2*i+1];
579: b[2] = 1.0; b[3] = -d->eigr[i_s+i];
580: SlepcUpdateVectorsS(&(*v)[i], n, 0.0, 1.0, &(*v)[i], 2*n, n,
581: b, 2, 2, 2);
582: } else {
583: /* [v_i v_i+1 r_i r_i+1]*= [tau_0' 0 1 0
584: 0 tau_0' 0 1
585: tau_1 0 -eigr_i -eigi_i
586: 0 tau_1 eigi_i -eigr_i ] */
587: b[0] = b[5] = PetscConj(theta[2*i]);
588: b[2] = b[7] = theta[2*i+1];
589: b[8] = b[13] = 1.0;
590: b[10] = b[15] = -d->eigr[i_s+i];
591: b[14] = -(b[11] = d->eigi[i_s+i]);
592: b[1] = b[3] = b[4] = b[6] = b[9] = b[12] = 0.0;
593: X[0] = (*v)[i]; X[1] = (*v)[i+1]; X[2] = r[i]; X[3] = r[i+1];
594: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 4, Z, size_Z);
595:
596: i++;
597: }
598: }
600: /* kr <- K^{-1}*r */
601: d->calcpairs_proj_res(d, i_s, i_e, r);
602: for(i=0; i<n; i++) {
603: d->improvex_precond(d, i_s+i, r[i], (*kr)[i]);
604: }
606: /* Free r */
607: auxV-= n;
609: /* Normalize the projector */
610: for(i=0; i<n; i++) DVD_NORMALIZE_UV((*u)[i],(*v)[i],ierr,a);
612: /* Return the next free vector */
613: *auxV_ = auxV;
615: return(0);
616: }
618: /*
619: Compute: u <- K^{-1}*B*X, v <- (theta[0]*A+theta[1]*B)*Y,
620: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
621: where
622: auxV, 4*(i_e-i_s) auxiliar global vectors
623: pX,pY, the right and left eigenvectors of the projected system
624: ld, the leading dimension of pX and pY
625: */
628: PetscErrorCode dvd_improvex_jd_proj_uv_KBXZY(dvdDashboard *d, PetscInt i_s,
629: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
630: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
631: {
632: PetscErrorCode ierr;
633: PetscInt n = i_e - i_s, i;
634: PetscScalar a, b[16];
635: Vec *Ax, *Bx, *r, *auxV = *auxV_, X[4];
636: /* The memory manager doen't allow to call a subroutines */
637: PetscScalar Z[size_Z];
641: /* Book space for u, v and kr */
642: *u = auxV; auxV+= n;
643: *kr = auxV; auxV+= n;
644: *v = auxV; auxV+= n;
645: r = auxV; auxV+= n;
647: /* u <- Y(i) */
648: SlepcUpdateVectorsZ(*u, 0.0, 1.0, d->W?d->W:d->V, d->size_V, pY, ld,
649: d->size_H, n);
651: /* v <- theta[0]A*u + theta[1]*B*u */
652: for(i=0; i<n; i++) {
653: MatMult(d->A, (*u)[i], (*v)[i]);
654: }
655: if (d->B) {
656: for(i=0; i<n; i++) {
657: MatMult(d->B, (*u)[i], (*kr)[i]);
658: }
659: Bx = *kr;
660: } else
661: Bx = *u;
663: for(i=0; i<n; i++) {
664: #ifndef PETSC_USE_COMPLEX
665: if(d->eigi[i_s+i] != 0.0) {
666: /* [v_i v_i+1 Bx_i Bx_i+1]*= [ theta_2i' 0
667: 0 theta_2i'
668: theta_2i+1 -thetai_i
669: thetai_i theta_2i+1 ] */
670: b[0] = b[5] = PetscConj(theta[2*i]);
671: b[2] = b[7] = -theta[2*i+1];
672: b[6] = -(b[3] = thetai[i]);
673: b[1] = b[4] = 0.0;
674: X[0] = (*v)[i]; X[1] = (*v)[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
675: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
676:
677: i++;
678: } else
679: #endif
680: {
681: /* v_i <- v_i*theta_2i' + Bx_i*theta_2i+1 */
682: VecAXPBY((*v)[i], theta[i*2+1], PetscConj(theta[i*2]), Bx[i]);
683:
684: }
685: }
687: /* Bx <- B*X(i) */
688: Bx = *kr;
689: for(i=i_s; i<i_e; i++) d->nX[i] = 1.0;
690: if (d->BV) {
691: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV, d->size_BV, pX, ld,
692: d->size_H, n);
693: } else {
694: SlepcUpdateVectorsZ(d->B?r:Bx, 0.0, 1.0, d->V, d->size_V, pX, ld,
695: d->size_H, n);
696: if (d->B) {
697: for(i=0; i<n; i++) {
698: MatMult(d->B, r[i], Bx[i]);
699: }
700: }
701: }
703: /* Ax <- A*X(i) */
704: Ax = r;
705: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV, d->size_AV, pX, ld,
706: d->size_H, n);
708: /* Recompute the eigenvalue */
709: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, *u, Ax, Bx, b, ierr);
711: /* u <- K^{-1} Bx */
712: for(i=0; i<n; i++) {
713: d->improvex_precond(d, i_s+i, Bx[i], (*u)[i]);
714: }
716: for(i=0; i<n; i++) {
717: if (d->eigi[i_s+i] == 0.0) {
718: /* r <- Ax -eig*Bx */
719: VecAXPBY(r[i], -d->eigr[i_s+i], 1.0, Bx[i]);
720: } else {
721: /* [r_i r_i+1 kr_i kr_i+1]*= [ 1 0
722: 0 1
723: -eigr_i -eigi_i
724: eigi_i -eigr_i] */
725: b[0] = b[5] = 1.0;
726: b[2] = b[7] = -d->eigr[i_s+i];
727: b[6] = -(b[3] = d->eigi[i_s+i]);
728: b[1] = b[4] = 0.0;
729: X[0] = r[i]; X[1] = r[i+1]; X[2] = (*kr)[i]; X[3] = (*kr)[i+1];
730: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
731:
732: i++;
733: }
734: }
736: /* kr <- K^{-1}*r */
737: d->calcpairs_proj_res(d, i_s, i_e, r);
738: for(i=0; i<n; i++) {
739: d->improvex_precond(d, i_s+i, r[i], (*kr)[i]);
740: }
742: /* Free r */
743: auxV-= n;
745: /* Normalize the projector */
746: for(i=0; i<n; i++) DVD_NORMALIZE_UV((*u)[i],(*v)[i],ierr,a);
748: /* Return the next free vector */
749: *auxV_ = auxV;
751: return(0);
752: }
754: /*
755: Compute: u <- K^{-1}*B*X, v <- X,
756: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
757: where
758: auxV, 4*(i_e-i_s) auxiliar global vectors
759: pX,pY, the right and left eigenvectors of the projected system
760: ld, the leading dimension of pX and pY
761: */
764: PetscErrorCode dvd_improvex_jd_proj_uv_KBXX(dvdDashboard *d, PetscInt i_s,
765: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
766: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
767: {
768: PetscErrorCode ierr;
769: PetscInt n = i_e - i_s, i;
770: PetscScalar a, b[16];
771: Vec *Ax, *Bx, *r, *auxV = *auxV_, X[4];
772: /* The memory manager doen't allow to call a subroutines */
773: PetscScalar Z[size_Z];
777: /* Book space for u, v and kr */
778: *u = auxV; auxV+= n;
779: *kr = auxV; auxV+= n;
780: *v = auxV; auxV+= n;
781: r = auxV; auxV+= n;
784: /* [v u] <- X(i) Y(i) */
785: SlepcUpdateVectorsZ(*v, 0.0, 1.0, d->V, d->size_V, pX, ld,
786: d->size_H, n);
787: SlepcUpdateVectorsZ(*u, 0.0, 1.0, d->W?d->W:d->V, d->size_V, pY, ld,
788: d->size_H, n);
790: /* Bx <- B*X(i) */
791: Bx = *kr;
792: for(i=i_s; i<i_e; i++) d->nX[i] = 1.0;
793: if (d->BV) {
794: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV, d->size_BV, pX, ld,
795: d->size_H, n);
796: } else {
797: if (d->B) {
798: for(i=0; i<n; i++) {
799: MatMult(d->B, (*v)[i], Bx[i]);
800: }
801: } else
802: Bx = *v;
803: }
805: /* Ax <- A*X(i) */
806: Ax = r;
807: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV, d->size_AV, pX, ld,
808: d->size_H, n);
810: /* Recompute the eigenvalue */
811: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, *u, Ax, Bx, b, ierr);
813: /* u <- K^{-1} Bx */
814: for(i=0; i<n; i++) {
815: d->improvex_precond(d, i_s+i, Bx[i], (*u)[i]);
816: }
818: for(i=0; i<n; i++) {
819: if (d->eigi[i_s+i] == 0.0) {
820: /* r <- Ax -eig*Bx */
821: VecAXPBY(r[i], -d->eigr[i_s+i], 1.0, Bx[i]);
822: } else {
823: /* [r_i r_i+1 kr_i kr_i+1]*= [ 1 0
824: 0 1
825: -eigr_i -eigi_i
826: eigi_i -eigr_i] */
827: b[0] = b[5] = 1.0;
828: b[2] = b[7] = -d->eigr[i_s+i];
829: b[6] = -(b[3] = d->eigi[i_s+i]);
830: b[1] = b[4] = 0.0;
831: X[0] = r[i]; X[1] = r[i+1]; X[2] = (*kr)[i]; X[3] = (*kr)[i+1];
832: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
833:
834: i++;
835: }
836: }
838: /* kr <- K^{-1}*r */
839: d->calcpairs_proj_res(d, i_s, i_e, r);
840: for(i=0; i<n; i++) {
841: d->improvex_precond(d, i_s+i, r[i], (*kr)[i]);
842: }
844: /* Free r */
845: auxV-= n;
847: /* Normalize the projector */
848: for(i=0; i<n; i++) DVD_NORMALIZE_UV((*u)[i],(*v)[i],ierr,a);
850: /* Return the next free vector */
851: *auxV_ = auxV;
853: return(0);
854: }
856: /*
857: Compute: u <- K^{-1}*B*X, v <- Y,
858: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- V*pY[i_s..i_e-1]
859: where
860: auxV, 4*(i_e-i_s) auxiliar global vectors
861: pX,pY, the right and left eigenvectors of the projected system
862: ld, the leading dimension of pX and pY
863: */
866: PetscErrorCode dvd_improvex_jd_proj_uv_KBXY(dvdDashboard *d, PetscInt i_s,
867: PetscInt i_e, Vec **u, Vec **v, Vec **kr, Vec **auxV_, PetscScalar *theta,
868: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
869: {
870: PetscErrorCode ierr;
871: PetscInt n = i_e - i_s, i;
872: PetscScalar a, b[16];
873: Vec *Ax, *Bx, *r, *auxV = *auxV_, X[4];
874: /* The memory manager doen't allow to call a subroutines */
875: PetscScalar Z[size_Z];
879: /* Book space for u, v and kr */
880: *u = auxV; auxV+= n;
881: *kr = auxV; auxV+= n;
882: *v = auxV; auxV+= n;
883: r = auxV; auxV+= n;
886: /* v <- Y(i) */
887: SlepcUpdateVectorsZ(*v, 0.0, 1.0, d->W?d->W:d->V, d->size_V, pY, ld,
888: d->size_H, n);
890: /* Bx <- B*X(i) */
891: Bx = *kr;
892: for(i=i_s; i<i_e; i++) d->nX[i] = 1.0;
893: if (d->BV) {
894: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV, d->size_BV, pX, ld,
895: d->size_H, n);
896: } else {
897: SlepcUpdateVectorsZ(d->B?*u:Bx, 0.0, 1.0, d->V, d->size_V, pX, ld,
898: d->size_H, n);
899: if (d->B) {
900: for(i=0; i<n; i++) {
901: MatMult(d->B, (*u)[i], Bx[i]);
902: }
903: }
904: }
906: /* Ax <- A*X(i) */
907: Ax = r;
908: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV, d->size_AV, pX, ld,
909: d->size_H, n);
911: /* Recompute the eigenvalue */
912: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, *v, Ax, Bx, b, ierr);
914: /* u <- K^{-1} Bx */
915: for(i=0; i<n; i++) {
916: d->improvex_precond(d, i_s+i, Bx[i], (*u)[i]);
917: }
919: for(i=0; i<n; i++) {
920: if (d->eigi[i_s+i] == 0.0) {
921: /* r <- Ax -eig*Bx */
922: VecAXPBY(r[i], -d->eigr[i_s+i], 1.0, Bx[i]);
923: } else {
924: /* [r_i r_i+1 kr_i kr_i+1]*= [ 1 0
925: 0 1
926: -eigr_i -eigi_i
927: eigi_i -eigr_i] */
928: b[0] = b[5] = 1.0;
929: b[2] = b[7] = -d->eigr[i_s+i];
930: b[6] = -(b[3] = d->eigi[i_s+i]);
931: b[1] = b[4] = 0.0;
932: X[0] = r[i]; X[1] = r[i+1]; X[2] = (*kr)[i]; X[3] = (*kr)[i+1];
933: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
934:
935: i++;
936: }
937: }
939: /* kr <- K^{-1}*r */
940: d->calcpairs_proj_res(d, i_s, i_e, r);
941: for(i=0; i<n; i++) {
942: d->improvex_precond(d, i_s+i, r[i], (*kr)[i]);
943: }
945: /* Free r */
946: auxV-= n;
948: /* Normalize the projector */
949: for(i=0; i<n; i++) DVD_NORMALIZE_UV((*u)[i],(*v)[i],ierr,a);
951: /* Return the next free vector */
952: *auxV_ = auxV;
954: return(0);
955: }
962: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d, dvdBlackboard *b,
963: PetscInt maxits, PetscReal tol,
964: PetscReal fix)
965: {
966: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
970: /* Setup the step */
971: if (b->state >= DVD_STATE_CONF) {
972: data->maxits = maxits;
973: data->tol = tol;
974: data->fix = fix;
975: d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
976: }
978: return(0);
979: }
984: PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d, PetscInt i,
985: PetscScalar* theta, PetscScalar* thetai, PetscInt *maxits, PetscReal *tol)
986: {
987: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
988: PetscReal a;
992: #ifndef PETSC_USE_COMPLEX
993: a = sqrt(d->eigr[i]*d->eigr[i]+d->eigi[i]*d->eigi[i]);
994: #else
995: a = PetscAbsScalar(d->eigr[i]);
996: #endif
998: if (d->nR[i]/a < data->fix) {
999: theta[0] = d->eigr[i];
1000: theta[1] = 1.0;
1001: #ifndef PETSC_USE_COMPLEX
1002: *thetai = d->eigi[i];
1003: #endif
1004: } else {
1005: theta[0] = d->target[0];
1006: theta[1] = d->target[1];
1007: #ifndef PETSC_USE_COMPLEX
1008: *thetai = 0.0;
1009: #endif
1010: }
1012: #ifdef PETSC_USE_COMPLEX
1013: if(thetai) *thetai = 0.0;
1014: #endif
1015: *maxits = data->maxits;
1016: *tol = data->tol;
1018: return(0);
1019: }
1022: /**** Patterns implementation *************************************************/
1024: typedef PetscInt (*funcV0_t)(dvdDashboard*, PetscInt, PetscInt, Vec*);
1025: typedef PetscInt (*funcV1_t)(dvdDashboard*, PetscInt, PetscInt, Vec*,
1026: PetscScalar*, Vec);
1028: /* Compute D <- K^{-1} * funcV[r_s..r_e] */
1031: PetscErrorCode dvd_improvex_PfuncV(dvdDashboard *d, void *funcV, Vec *D,
1032: PetscInt max_size_D, PetscInt r_s, PetscInt r_e, Vec *auxV,
1033: PetscScalar *auxS)
1034: {
1035: PetscErrorCode ierr;
1036: PetscInt i;
1040: if (max_size_D >= r_e-r_s+1) {
1041: /* The optimized version needs one vector extra of D */
1042: /* D(1:r.size) = R(r_s:r_e-1) */
1043: if (auxS) ((funcV1_t)funcV)(d, r_s, r_e, D+1, auxS, auxV[0]);
1044: else ((funcV0_t)funcV)(d, r_s, r_e, D+1);
1046: /* D = K^{-1} * R */
1047: for (i=0; i<r_e-r_s; i++) {
1048: d->improvex_precond(d, i+r_s, D[i+1], D[i]);
1049: }
1050: } else if (max_size_D == r_e-r_s) {
1051: /* Non-optimized version */
1052: /* auxV <- R[r_e-1] */
1053: if (auxS) ((funcV1_t)funcV)(d, r_e-1, r_e, auxV, auxS, auxV[1]);
1054: else ((funcV0_t)funcV)(d, r_e-1, r_e, auxV);
1056: /* D(1:r.size-1) = R(r_s:r_e-2) */
1057: if (auxS) ((funcV1_t)funcV)(d, r_s, r_e-1, D+1, auxS, auxV[1]);
1058: else ((funcV0_t)funcV)(d, r_s, r_e-1, D+1);
1060: /* D = K^{-1} * R */
1061: for (i=0; i<r_e-r_s-1; i++) {
1062: d->improvex_precond(d, i+r_s, D[i+1], D[i]);
1063: }
1064: d->improvex_precond(d, r_e-1, auxV[0], D[r_e-r_s-1]);
1065: } else {
1066: SETERRQ(1, "Problem: r_e-r_s > max_size_D!");
1067: }
1069: return(0);
1070: }
1073: /* Compute the left and right projected eigenvectors where,
1074: pX, the returned right eigenvectors
1075: pY, the returned left eigenvectors,
1076: ld_, the leading dimension of pX and pY,
1077: auxS, auxiliar vector of size length 6*d->size_H
1078: */
1081: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
1082: PetscScalar *pY, PetscInt ld, PetscScalar *auxS, PetscInt size_auxS)
1083: {
1084: PetscErrorCode ierr;
1085:
1088: SlepcDenseCopy(pY, ld, d->T?d->pY:d->pX, d->ldpX, d->size_H,
1089: d->size_H);
1090: SlepcDenseCopy(pX, ld, d->pX, d->ldpX, d->size_H, d->size_H);
1091:
1092:
1093: /* [qX, qY] <- eig(S, T); pX <- d->pX * qX; pY <- d->pY * qY */
1094: dvd_compute_eigenvectors(d->size_H, d->S, d->ldS, d->T, d->ldT, pX,
1095: ld, pY, ld, auxS, size_auxS, PETSC_TRUE);
1096:
1098: /* 2-Normalize the columns of pX an pY */
1099: SlepcDenseNorm(pX, ld, d->size_H, d->size_H, d->eigi);
1100: SlepcDenseNorm(pY, ld, d->size_H, d->size_H, d->eigi);
1102: return(0);
1103: }