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-2011, Universitat 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 <slepcvec.h>
28: #include <slepcblaslapack.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_start(dvdDashboard *d);
37: PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d);
38: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
39: PetscInt max_size_D, PetscInt r_s,
40: PetscInt r_e, PetscInt *size_D);
41: PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d, PetscInt i_s,
42: PetscInt i_e, Vec **u, Vec **v, Vec *kr, Vec **auxV, PetscScalar **auxS,
43: PetscScalar *theta, PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY,
44: PetscInt ld);
45: PetscErrorCode dvd_improvex_jd_proj_uv_KXX(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_KZX(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);
57: PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d, Vec *V, PetscInt cV, PetscScalar *auxS);
59: #define size_Z (64*4)
61: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r ****/
63: typedef struct {
64: PetscInt size_X;
65: void
66: *old_improveX_data; /* old improveX_data */
67: improveX_type
68: old_improveX; /* old improveX */
69: KSP ksp; /* correction equation solver */
70: Vec
71: friends, /* reference vector for composite vectors */
72: *auxV; /* auxiliar vectors */
73: PetscScalar *auxS, /* auxiliar scalars */
74: *theta,
75: *thetai; /* the shifts used in the correction eq. */
76: PetscInt maxits, /* maximum number of iterations */
77: r_s, r_e, /* the selected eigenpairs to improve */
78: ksp_max_size; /* the ksp maximum subvectors size */
79: PetscReal tol, /* the maximum solution tolerance */
80: lastTol, /* last tol for dynamic stopping criterion */
81: fix; /* tolerance for using the approx. eigenvalue */
82: PetscBool
83: dynamic; /* if the dynamic stopping criterion is applied */
84: dvdDashboard
85: *d; /* the currect dvdDashboard reference */
86: PC old_pc; /* old pc in ksp */
87: Vec
88: *u, /* new X vectors */
89: *real_KZ, /* original KZ */
90: *KZ; /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
91: PetscScalar
92: *XKZ, /* X'*KZ */
93: *iXKZ; /* inverse of XKZ */
94: PetscInt
95: ldXKZ, /* leading dimension of XKZ */
96: size_iXKZ, /* size of iXKZ */
97: ldiXKZ, /* leading dimension of iXKZ */
98: size_KZ, /* size of converged KZ */
99: size_real_KZ, /* original size of KZ */
100: size_cX, /* last value of d->size_cX */
101: old_size_X; /* last number of improved vectors */
102: PetscBLASInt
103: *iXKZPivots; /* array of pivots */
104: } dvdImprovex_jd;
106: #define _Ceil(A,B) ((A)/(B)+((A)%(B)==0?0:1))
107: #define FromIntToScalar(S) ((PetscInt)_Ceil((S)*sizeof(PetscBLASInt),sizeof(PetscScalar)))
111: PetscErrorCode dvd_improvex_jd(dvdDashboard *d, dvdBlackboard *b, KSP ksp,
112: PetscInt max_bs, PetscInt cX_impr, PetscBool dynamic)
113: {
114: PetscErrorCode ierr;
115: dvdImprovex_jd *data;
116: PetscBool t, herm = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
117: PC pc;
118: PetscInt size_P;
122: /* Setting configuration constrains */
123: PetscTypeCompare((PetscObject)ksp, KSPPREONLY, &t);
124: if (t) ksp = PETSC_NULL;
126: /* If the arithmetic is real and the problem is not Hermitian, then
127: the block size is incremented in one */
128: #if !defined(PETSC_USE_COMPLEX)
129: if (!herm) {
130: max_bs++;
131: b->max_size_P = PetscMax(b->max_size_P, 2);
132: } else
133: #endif
134: b->max_size_P = PetscMax(b->max_size_P, 1);
135: b->max_size_X = PetscMax(b->max_size_X, max_bs);
136: size_P = b->max_size_P+cX_impr;
137: b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X*(ksp?3:2));
138: /* u, kr?, auxV */
139: b->own_scalars+= size_P*size_P; /* XKZ */
140: b->max_size_auxS = PetscMax(b->max_size_auxS,
141: (herm?0:1)*2*b->max_size_proj*b->max_size_proj + /* pX, pY */
142: b->max_size_X*3 + /* theta, thetai */
143: size_P*size_P + /* iXKZ */
144: FromIntToScalar(size_P) + /* iXkZPivots */
145: PetscMax(PetscMax(
146: 3*b->max_size_proj*b->max_size_X, /* dvd_improvex_apply_proj */
147: 8*cX_impr*b->max_size_X), /* dvd_improvex_jd_proj_cuv_KZX */
148: (herm?0:1)*6*b->max_size_proj)); /* dvd_improvex_get_eigenvectors */
149: b->own_vecs+= size_P; /* KZ */
151: /* Setup the preconditioner */
152: if (ksp) {
153: KSPGetPC(ksp, &pc);
154: dvd_static_precond_PC(d, b, pc);
155: } else {
156: dvd_static_precond_PC(d, b, 0);
157: }
159: /* Setup the step */
160: if (b->state >= DVD_STATE_CONF) {
161: PetscMalloc(sizeof(dvdImprovex_jd), &data);
162: data->dynamic = dynamic;
163: data->size_real_KZ = size_P;
164: data->real_KZ = b->free_vecs; b->free_vecs+= data->size_real_KZ;
165: d->max_cX_in_impr = cX_impr;
166: data->XKZ = b->free_scalars; b->free_scalars+= size_P*size_P;
167: data->ldXKZ = size_P;
168: data->size_X = b->max_size_X;
169: data->old_improveX_data = d->improveX_data;
170: d->improveX_data = data;
171: data->old_improveX = d->improveX;
172: data->ksp = ksp;
173: data->d = d;
174: d->improveX = dvd_improvex_jd_gen;
175: data->ksp_max_size = max_bs;
177: DVD_FL_ADD(d->startList, dvd_improvex_jd_start);
178: DVD_FL_ADD(d->endList, dvd_improvex_jd_end);
179: DVD_FL_ADD(d->destroyList, dvd_improvex_jd_d);
180: }
182: return(0);
183: }
188: PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
189: {
190: PetscErrorCode ierr;
191: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
192: PetscInt rA, cA, rlA, clA;
193: Mat A;
194: PetscBool t;
195: PC pc;
199: data->KZ = data->real_KZ;
200: data->size_KZ = data->size_cX = data->old_size_X = 0;
201: data->lastTol = data->dynamic?0.5:0.0;
203: /* Setup the ksp */
204: if(data->ksp) {
205: /* Create the reference vector */
206: VecCreateCompWithVecs(d->V, data->ksp_max_size, PETSC_NULL,
207: &data->friends);
209: /* Save the current pc and set a PCNONE */
210: KSPGetPC(data->ksp, &data->old_pc);
211: PetscTypeCompare((PetscObject)data->old_pc, PCNONE, &t);
212:
213: data->lastTol = 0.5;
214: if (t) {
215: data->old_pc = 0;
216: } else {
217: PetscObjectReference((PetscObject)data->old_pc);
218: PCCreate(((PetscObject)d->eps)->comm, &pc);
219: PCSetType(pc, PCNONE);
220: PCSetOperators(pc, d->A, d->A, SAME_PRECONDITIONER);
221: KSPSetPC(data->ksp, pc);
222: PCDestroy(&pc);
223: }
225: /* Create the (I-v*u')*K*(A-s*B) matrix */
226: MatGetSize(d->A, &rA, &cA);
227: MatGetLocalSize(d->A, &rlA, &clA);
228: MatCreateShell(((PetscObject)d->A)->comm, rlA*data->ksp_max_size,
229: clA*data->ksp_max_size, rA*data->ksp_max_size,
230: cA*data->ksp_max_size, data, &A);
231: MatShellSetOperation(A, MATOP_MULT,
232: (void(*)(void))dvd_matmult_jd);
233: MatShellSetOperation(A, MATOP_GET_VECS,
234: (void(*)(void))dvd_matgetvecs_jd);
235:
237: /* Try to avoid KSPReset */
238: KSPGetOperatorsSet(data->ksp,&t,PETSC_NULL);
239: if (t) {
240: Mat M;
241: PetscInt rM;
242: KSPGetOperators(data->ksp,&M,PETSC_NULL,PETSC_NULL);
243: MatGetSize(M,&rM,PETSC_NULL);
244: if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
245: }
246: KSPSetOperators(data->ksp, A, A, SAME_PRECONDITIONER);
247:
248: KSPSetUp(data->ksp);
249: MatDestroy(&A);
250: } else {
251: data->old_pc = 0;
252: data->friends = PETSC_NULL;
253: }
254:
255: return(0);
256: }
261: PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
262: {
263: PetscErrorCode ierr;
264: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
268: if (data->friends) { VecDestroy(&data->friends); }
270: /* Restore the pc of ksp */
271: if (data->old_pc) {
272: KSPSetPC(data->ksp, data->old_pc);
273: PCDestroy(&data->old_pc);
274: }
276: return(0);
277: }
282: PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
283: {
284: PetscErrorCode ierr;
285: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
288:
289: /* Restore changes in dvdDashboard */
290: d->improveX_data = data->old_improveX_data;
292: /* Free local data and objects */
293: PetscFree(data);
295: return(0);
296: }
301: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
302: PetscInt max_size_D, PetscInt r_s,
303: PetscInt r_e, PetscInt *size_D)
304: {
305: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
306: PetscErrorCode ierr;
307: PetscInt i, j, n, maxits, maxits0, lits, s, ldpX;
308: PetscScalar *pX, *pY, *auxS = d->auxS, *auxS0;
309: PetscReal tol, tol0;
310: Vec *u, *v, *kr, kr_comp, D_comp;
314: /* Quick exit */
315: if ((max_size_D == 0) || r_e-r_s <= 0) {
316: /* Callback old improveX */
317: if (data->old_improveX) {
318: d->improveX_data = data->old_improveX_data;
319: data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
320: d->improveX_data = data;
321: }
322: return(0);
323: }
324:
325: n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
326: if (n == 0) SETERRQ(PETSC_COMM_SELF,1, "n == 0!\n");
327: if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1, "size_X < r_e-r_s!\n");
329: /* Compute the eigenvectors of the selected pairs */
330: if (DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
331: pX = d->pX;
332: pY = d->pY?d->pY:d->pX;
333: ldpX = d->ldpX;
334: } else {
335: pX = auxS; auxS+= d->size_H*d->size_H;
336: pY = auxS; auxS+= d->size_H*d->size_H;
337: dvd_improvex_get_eigenvectors(d, pX, pY, d->size_H, auxS,
338: d->size_auxS-(auxS-d->auxS));
339:
340: ldpX = d->size_H;
341: }
343: /* Restart lastTol if a new pair converged */
344: if (data->dynamic && data->size_cX < d->size_cX)
345: data->lastTol = 0.5;
347: for(i=0, s=0, auxS0=auxS; i<n; i+=s) {
348: /* If the selected eigenvalue is complex, but the arithmetic is real... */
349: #if !defined(PETSC_USE_COMPLEX)
350: if (PetscAbsScalar(d->eigi[i] != 0.0)) {
351: if (i+2 <= max_size_D) s=2; else break;
352: } else
353: #endif
354: s=1;
356: data->auxV = d->auxV;
357: data->r_s = r_s+i; data->r_e = r_s+i+s;
358: auxS = auxS0;
359: data->theta = auxS; auxS+= 2*s;
360: data->thetai = auxS; auxS+= s;
361: /* If GD, kr = D */
362: if (!data->ksp) {
363: kr = &D[i];
365: /* If JD, kr = auxV */
366: } else {
367: kr = data->auxV; data->auxV+= s;
368: }
370: /* Compute theta, maximum iterations and tolerance */
371: maxits = 0; tol = 1;
372: for(j=0; j<s; j++) {
373: d->improvex_jd_lit(d, r_s+i+j, &data->theta[2*j],
374: &data->thetai[j], &maxits0, &tol0);
375:
376: maxits+= maxits0; tol*= tol0;
377: }
378: maxits/= s; tol = data->dynamic?data->lastTol:exp(log(tol)/s);
380: /* Compute u, v and kr */
381: dvd_improvex_jd_proj_cuv(d, r_s+i, r_s+i+s, &u, &v, kr,
382: &data->auxV, &auxS, data->theta, data->thetai,
383: &pX[d->size_H*(r_s+i+d->cX_in_H)], &pY[d->size_H*(r_s+i+d->cX_in_H)], ldpX);
384:
385: data->u = u;
387: /* Check if the first eigenpairs are converged */
388: if (i == 0) {
389: d->preTestConv(d, 0, s, s, PETSC_NULL, PETSC_NULL, &d->npreconv);
390:
391: if (d->npreconv > 0) break;
392: }
394: /* Compute kr <- kr - v*(u'*kr) */
395: dvd_improvex_apply_proj(d, kr, s, auxS);
397: /* If JD */
398: if (data->ksp) {
399: data->auxS = auxS;
401: /* kr <- -kr */
402: for(j=0; j<s; j++) {
403: VecScale(kr[j], -1.0);
404: }
406: /* Compouse kr and D */
407: VecCreateCompWithVecs(kr, data->ksp_max_size, data->friends,
408: &kr_comp);
409: VecCreateCompWithVecs(&D[i], data->ksp_max_size, data->friends,
410: &D_comp);
411: VecCompSetSubVecs(data->friends,s,PETSC_NULL);
412:
413: /* Solve the correction equation */
414: KSPSetTolerances(data->ksp, tol, PETSC_DEFAULT, PETSC_DEFAULT,
415: maxits);
416: KSPSolve(data->ksp, kr_comp, D_comp);
417: KSPGetIterationNumber(data->ksp, &lits);
418: d->eps->OP->lineariterations+= lits;
419:
420: /* Destroy the composed ks and D */
421: VecDestroy(&kr_comp);
422: VecDestroy(&D_comp);
423: }
424: }
425: *size_D = i;
426: if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
427:
428: /* Callback old improveX */
429: if (data->old_improveX) {
430: d->improveX_data = data->old_improveX_data;
431: data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
432: d->improveX_data = data;
433: }
435: return(0);
436: }
441: PetscErrorCode dvd_matmult_jd(Mat A, Vec in, Vec out)
442: {
443: PetscErrorCode ierr;
444: dvdImprovex_jd *data;
445: PetscInt n, i;
446: const Vec *inx, *outx, *Bx;
450: MatShellGetContext(A, (void**)&data);
451: VecCompGetSubVecs(in,PETSC_NULL,&inx);
452: VecCompGetSubVecs(out,PETSC_NULL,&outx);
453: n = data->r_e - data->r_s;
455: /* aux <- theta[1]A*in - theta[0]*B*in */
456: for(i=0; i<n; i++) {
457: MatMult(data->d->A, inx[i], data->auxV[i]);
458: }
459: if (data->d->B) {
460: for(i=0; i<n; i++) {
461: MatMult(data->d->B, inx[i], outx[i]);
462: }
463: Bx = outx;
464: } else
465: Bx = inx;
467: for(i=0; i<n; i++) {
468: #if !defined(PETSC_USE_COMPLEX)
469: if(data->d->eigi[data->r_s+i] != 0.0) {
470: /* aux_i <- [ t_2i+1*A*inx_i - t_2i*Bx_i + ti_i*Bx_i+1;
471: aux_i+1 t_2i+1*A*inx_i+1 - ti_i*Bx_i - t_2i*Bx_i+1 ] */
472: VecAXPBYPCZ(data->auxV[i], -data->theta[2*i], data->thetai[i],
473: data->theta[2*i+1], Bx[i], Bx[i+1]);
474:
475: VecAXPBYPCZ(data->auxV[i+1], -data->thetai[i],
476: -data->theta[2*i], data->theta[2*i+1], Bx[i],
477: Bx[i+1]);
478: i++;
479: } else
480: #endif
481: {
482: VecAXPBY(data->auxV[i], -data->theta[i*2], data->theta[i*2+1],
483: Bx[i]);
484: }
485: }
487: /* out <- K * aux */
488: for(i=0; i<n; i++) {
489: data->d->improvex_precond(data->d, data->r_s+i, data->auxV[i],
490: outx[i]);
491: }
493: /* out <- out - v*(u'*out) */
494: dvd_improvex_apply_proj(data->d, (Vec*)outx, n, data->auxS);
496: return(0);
497: }
502: PetscErrorCode dvd_matgetvecs_jd(Mat A, Vec *right, Vec *left)
503: {
504: PetscErrorCode ierr;
505: Vec *r, *l;
506: dvdImprovex_jd *data;
507: PetscInt n, i;
511: MatShellGetContext(A, (void**)&data);
512: n = data->ksp_max_size;
513: if (right) {
514: PetscMalloc(sizeof(Vec)*n, &r);
515: }
516: if (left) {
517: PetscMalloc(sizeof(Vec)*n, &l);
518: }
519: for (i=0; i<n; i++) {
520: MatGetVecs(data->d->A, right?&r[i]:PETSC_NULL,
521: left?&l[i]:PETSC_NULL);
522: }
523: if(right) {
524: VecCreateCompWithVecs(r, n, data->friends, right);
525: for (i=0; i<n; i++) {
526: VecDestroy(&r[i]);
527: }
528: }
529: if(left) {
530: VecCreateCompWithVecs(l, n, data->friends, left);
531: for (i=0; i<n; i++) {
532: VecDestroy(&l[i]);
533: }
534: }
536: if (right) {
537: PetscFree(r);
538: }
539: if (left) {
540: PetscFree(l);
541: }
543: return(0);
544: }
549: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d, dvdBlackboard *b,
550: ProjType_t p)
551: {
554: /* Setup the step */
555: if (b->state >= DVD_STATE_CONF) {
556: switch(p) {
557: case DVD_PROJ_KXX:
558: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KXX; break;
559: case DVD_PROJ_KZX:
560: d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX; break;
561: }
562: }
564: return(0);
565: }
569: /*
570: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
571: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
572: where
573: auxV, 4*(i_e-i_s) auxiliar global vectors
574: pX,pY, the right and left eigenvectors of the projected system
575: ld, the leading dimension of pX and pY
576: auxS: max(8*bs*max_cX_in_proj,size_V*size_V)
577: */
578: PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d, PetscInt i_s,
579: PetscInt i_e, Vec **u, Vec **v, Vec *kr, Vec **auxV, PetscScalar **auxS,
580: PetscScalar *theta, PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY,
581: PetscInt ld)
582: {
583: PetscErrorCode ierr;
584: PetscInt n = i_e - i_s, size_KZ, V_new, rm, i;
585: PetscScalar *wS0, *wS1;
586: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
587: PetscBLASInt s, ldXKZ, info;
591: /* Check consistency */
592: V_new = d->size_cX - data->size_cX;
593: if (V_new > data->old_size_X) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
594: data->old_size_X = n;
596: /* KZ <- KZ(rm:rm+max_cX-1) */
597: rm = PetscMax(V_new+data->size_KZ-d->max_cX_in_impr, 0);
598: for (i=0; i<d->max_cX_in_impr; i++) {
599: VecCopy(data->KZ[i+rm], data->KZ[i]);
600: }
601: data->size_cX = d->size_cX;
603: /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
604: for (i=0; i<d->max_cX_in_impr; i++) {
605: SlepcDenseCopy(&data->XKZ[i*data->ldXKZ+i], data->ldXKZ, &data->XKZ[(i+rm)*data->ldXKZ+i+rm], data->ldXKZ, data->size_KZ, 1);
606: }
607: data->size_KZ = PetscMin(d->max_cX_in_impr, data->size_KZ+V_new);
609: /* Compute X, KZ and KR */
610: *u = *auxV; *auxV+= n;
611: *v = &data->KZ[data->size_KZ];
612: d->improvex_jd_proj_uv(d, i_s, i_e, *u, *v, kr, *auxV, theta, thetai,
613: pX, pY, ld);
615: /* XKZ <- X'*KZ */
616: size_KZ = data->size_KZ+n;
617: wS0 = *auxS; wS1 = wS0+2*n*data->size_KZ+n*n;
618: VecsMult(data->XKZ, 0, data->ldXKZ, d->V-data->size_KZ, 0, data->size_KZ, data->KZ, data->size_KZ, size_KZ, wS0, wS1);
619: VecsMult(&data->XKZ[data->size_KZ], 0, data->ldXKZ, *u, 0, n, data->KZ, 0, size_KZ, wS0, wS1);
621: /* iXKZ <- inv(XKZ) */
622: s = PetscBLASIntCast(size_KZ);
623: data->ldiXKZ = data->size_iXKZ = size_KZ;
624: data->iXKZ = *auxS; *auxS+= size_KZ*size_KZ;
625: data->iXKZPivots = (PetscBLASInt*)*auxS;
626: *auxS += FromIntToScalar(size_KZ);
627: SlepcDenseCopy(data->iXKZ,data->ldiXKZ,data->XKZ,data->ldXKZ,size_KZ,size_KZ);
628: ldXKZ = PetscBLASIntCast(data->ldiXKZ);
629: LAPACKgetrf_(&s, &s, data->iXKZ, &ldXKZ, data->iXKZPivots, &info);
630: if (info) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_LIB, "Error in Lapack XGETRF %d", info);
632: return(0);
633: }
635: #define DVD_COMPLEX_RAYLEIGH_QUOTIENT(ur,ui,Axr,Axi,Bxr,Bxi,eigr,eigi,b,ierr)\
636: { \
637: VecDot((Axr), (ur), &(b)[0]); /* r*A*r */ \
638: VecDot((Axr), (ui), &(b)[1]); /* i*A*r */ \
639: VecDot((Axi), (ur), &(b)[2]); /* r*A*i */ \
640: VecDot((Axi), (ui), &(b)[3]); /* i*A*i */ \
641: VecDot((Bxr), (ur), &(b)[4]); /* r*B*r */ \
642: VecDot((Bxr), (ui), &(b)[5]); /* i*B*r */ \
643: VecDot((Bxi), (ur), &(b)[6]); /* r*B*i */ \
644: VecDot((Bxi), (ui), &(b)[7]); /* i*B*i */ \
645: (b)[0] = (b)[0]+(b)[3]; /* rAr+iAi */ \
646: (b)[2] = (b)[2]-(b)[1]; /* rAi-iAr */ \
647: (b)[4] = (b)[4]+(b)[7]; /* rBr+iBi */ \
648: (b)[6] = (b)[6]-(b)[5]; /* rBi-iBr */ \
649: (b)[7] = (b)[4]*(b)[4] + (b)[6]*(b)[6]; /* k */ \
650: *(eigr) = ((b)[0]*(b)[4] + (b)[2]*(b)[6]) / (b)[7]; /* eig_r */ \
651: *(eigi) = ((b)[2]*(b)[4] - (b)[0]*(b)[6]) / (b)[7]; /* eig_i */ \
652: }
654: #if !defined(PETSC_USE_COMPLEX)
655: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
656: for((i)=0; (i)<(n); (i)++) { \
657: if ((eigi)[(i_s)+(i)] != 0.0) { \
658: /* eig_r = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k \
659: eig_i = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k \
660: k = (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr) */ \
661: DVD_COMPLEX_RAYLEIGH_QUOTIENT((u)[(i)], (u)[(i)+1], (Ax)[(i)], \
662: (Ax)[(i)+1], (Bx)[(i)], (Bx)[(i)+1], &(b)[8], &(b)[9], (b), (ierr)); \
663: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[8])/ \
664: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10 || \
665: PetscAbsScalar((eigi)[(i_s)+(i)] - (b)[9])/ \
666: PetscAbsScalar((eigi)[(i_s)+(i)]) > 1e-10 ) { \
667: (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its "\
668: "Rayleigh quotient value %G+%G\n", \
669: (eigr)[(i_s)+(i)], \
670: (eigi)[(i_s)+1], (b)[8], (b)[9]); \
671: } \
672: (i)++; \
673: } \
674: }
675: #else
676: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
677: for((i)=0; (i)<(n); (i)++) { \
678: (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]); \
679: (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]); \
680: (b)[0] = (b)[0]/(b)[1]; \
681: if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
682: PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10 ) { \
683: (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its " \
684: "Rayleigh quotient value %G+%G\n", \
685: PetscRealPart((eigr)[(i_s)+(i)]), \
686: PetscImaginaryPart((eigr)[(i_s)+(i)]), PetscRealPart((b)[0]), \
687: PetscImaginaryPart((b)[0])); \
688: } \
689: }
690: #endif
694: /*
695: Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
696: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
697: where
698: auxV, 4*(i_e-i_s) auxiliar global vectors
699: pX,pY, the right and left eigenvectors of the projected system
700: ld, the leading dimension of pX and pY
701: */
702: PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d, PetscInt i_s,
703: PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
704: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
705: {
706: PetscErrorCode ierr;
707: PetscInt n = i_e - i_s, i;
708: PetscScalar b[16];
709: Vec *Ax, *Bx, *r=auxV, X[4];
710: /* The memory manager doen't allow to call a subroutines */
711: PetscScalar Z[size_Z];
715: /* u <- X(i) */
716: SlepcUpdateVectorsZ(u, 0.0, 1.0, d->V-d->cX_in_H, d->size_V+d->cX_in_H, pX, ld, d->size_H, n);
717: /* nX(i) <- ||X(i)|| */
718: if (d->correctXnorm) {
719: for (i=0; i<n; i++) {
720: VecNormBegin(u[i], NORM_2, &d->nX[i_s+i]);
721: }
722: for (i=0; i<n; i++) {
723: VecNormEnd(u[i], NORM_2, &d->nX[i_s+i]);
724: }
725: #if !defined(PETSC_USE_COMPLEX)
726: for(i=0; i<n; i++)
727: if(d->eigi[i_s+i] != 0.0)
728: 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]);
729: #endif
730: } else {
731: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
732: }
734: /* v <- theta[0]A*u + theta[1]*B*u */
736: /* Bx <- B*X(i) */
737: Bx = kr;
738: if (d->BV) {
739: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV-d->cX_in_H, d->size_BV+d->cX_in_H, pX, ld, d->size_H, n);
740: } else {
741: for(i=0; i<n; i++) {
742: if (d->B) {
743: MatMult(d->B, u[i], Bx[i]);
744: } else {
745: VecCopy(u[i], Bx[i]);
746: }
747: }
748: }
750: /* Ax <- A*X(i) */
751: Ax = r;
752: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV-d->cX_in_H, d->size_AV+d->cX_in_H, pX, ld, d->size_H, n);
754: /* v <- Y(i) */
755: SlepcUpdateVectorsZ(v, 0.0, 1.0, (d->W?d->W:d->V)-d->cX_in_H, d->size_V+d->cX_in_H, pY, ld, d->size_H, n);
757: /* Recompute the eigenvalue */
758: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, v, Ax, Bx, b, ierr);
760: for(i=0; i<n; i++) {
761: #if !defined(PETSC_USE_COMPLEX)
762: if(d->eigi[i_s+i] != 0.0) {
763: /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i' 0 1 0
764: 0 theta_2i' 0 1
765: theta_2i+1 -thetai_i -eigr_i -eigi_i
766: thetai_i theta_2i+1 eigi_i -eigr_i ] */
767: b[0] = b[5] = PetscConj(theta[2*i]);
768: b[2] = b[7] = -theta[2*i+1];
769: b[6] = -(b[3] = thetai[i]);
770: b[1] = b[4] = 0.0;
771: b[8] = b[13] = 1.0/d->nX[i_s+i];
772: b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
773: b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
774: b[9] = b[12] = 0.0;
775: X[0] = Ax[i]; X[1] = Ax[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
776: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 4, Z, size_Z);
777:
778: i++;
779: } else
780: #endif
781: {
782: /* [Ax_i Bx_i]*= [ theta_2i' 1/nX_i
783: theta_2i+1 -eig_i/nX_i ] */
784: b[0] = PetscConj(theta[i*2]);
785: b[1] = theta[i*2+1];
786: b[2] = 1.0/d->nX[i_s+i];
787: b[3] = -d->eigr[i_s+i]/d->nX[i_s+i];
788: X[0] = Ax[i]; X[1] = Bx[i];
789: SlepcUpdateVectorsD(X, 2, 1.0, b, 2, 2, 2, Z, size_Z);
790:
791: }
792: }
793: for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
795: /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
796: for(i=0; i<n; i++) {
797: d->improvex_precond(d, i_s+i, r[i], v[i]);
798: }
800: /* kr <- K^{-1}*kr = K^{-1}*(Ax - eig_i*Bx) */
801: d->calcpairs_proj_res(d, i_s, i_e, Bx);
802: for(i=0; i<n; i++) {
803: VecCopy(Bx[i], r[i]);
804: d->improvex_precond(d, i_s+i, r[i], kr[i]);
805: }
807: return(0);
808: }
812: /*
813: Compute: u <- K^{-1}*X, v <- X,
814: kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
815: where
816: auxV, 4*(i_e-i_s) auxiliar global vectors
817: pX,pY, the right and left eigenvectors of the projected system
818: ld, the leading dimension of pX and pY
819: */
820: PetscErrorCode dvd_improvex_jd_proj_uv_KXX(dvdDashboard *d, PetscInt i_s,
821: PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
822: PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
823: {
824: PetscErrorCode ierr;
825: PetscInt n = i_e - i_s, i;
826: PetscScalar b[16];
827: Vec *Ax, *Bx, *r = auxV, X[4];
828: /* The memory manager doen't allow to call a subroutines */
829: PetscScalar Z[size_Z];
833: /* [v u] <- X(i) Y(i) */
834: SlepcUpdateVectorsZ(v, 0.0, 1.0, d->V-d->cX_in_H, d->size_V+d->cX_in_H, pX, ld, d->size_H, n);
835: SlepcUpdateVectorsZ(u, 0.0, 1.0, (d->W?d->W:d->V)-d->cX_in_H, d->size_V+d->cX_in_H, pY, ld, d->size_H, n);
837: /* Bx <- B*X(i) */
838: Bx = kr;
839: for(i=i_s; i<i_e; i++) d->nX[i] = 1.0;
840: if (d->BV) {
841: SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV-d->cX_in_H, d->size_BV+d->cX_in_H, pX, ld, d->size_H, n);
842: } else {
843: if (d->B) {
844: for(i=0; i<n; i++) {
845: MatMult(d->B, v[i], Bx[i]);
846: }
847: } else
848: Bx = v;
849: }
851: /* Ax <- A*X(i) */
852: Ax = r;
853: SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV-d->cX_in_H, d->size_AV+d->cX_in_H, pX, ld, d->size_H, n);
855: /* Recompute the eigenvalue */
856: DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, u, Ax, Bx, b, ierr);
858: for(i=0; i<n; i++) {
859: if (d->eigi[i_s+i] == 0.0) {
860: /* r <- Ax -eig*Bx */
861: VecAXPBY(r[i], -d->eigr[i_s+i], 1.0, Bx[i]);
862: } else {
863: /* [r_i r_i+1 kr_i kr_i+1]*= [ 1 0
864: 0 1
865: -eigr_i -eigi_i
866: eigi_i -eigr_i] */
867: b[0] = b[5] = 1.0;
868: b[2] = b[7] = -d->eigr[i_s+i];
869: b[6] = -(b[3] = d->eigi[i_s+i]);
870: b[1] = b[4] = 0.0;
871: X[0] = r[i]; X[1] = r[i+1]; X[2] = kr[i]; X[3] = kr[i+1];
872: SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
873:
874: i++;
875: }
876: }
878: /* kr <- K^{-1}*r */
879: d->calcpairs_proj_res(d, i_s, i_e, r);
880: for(i=0; i<n; i++) {
881: d->improvex_precond(d, i_s+i, r[i], kr[i]);
882: }
884: /* u <- K^{-1} X(i) */
885: for(i=0; i<n; i++) {
886: d->improvex_precond(d, i_s+i, v[i], u[i]);
887: }
889: return(0);
890: }
895: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d, dvdBlackboard *b,
896: PetscInt maxits, PetscReal tol,
897: PetscReal fix)
898: {
899: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
903: /* Setup the step */
904: if (b->state >= DVD_STATE_CONF) {
905: data->maxits = maxits;
906: data->tol = tol;
907: data->fix = fix;
908: d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
909: }
911: return(0);
912: }
917: PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d, PetscInt i,
918: PetscScalar* theta, PetscScalar* thetai, PetscInt *maxits, PetscReal *tol)
919: {
920: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
921: PetscReal a;
924: a = SlepcAbsEigenvalue(d->eigr[i],d->eigi[i]);
926: if (d->nR[i]/a < data->fix) {
927: theta[0] = d->eigr[i];
928: theta[1] = 1.0;
929: #if !defined(PETSC_USE_COMPLEX)
930: *thetai = d->eigi[i];
931: #endif
932: } else {
933: theta[0] = d->target[0];
934: theta[1] = d->target[1];
935: #if !defined(PETSC_USE_COMPLEX)
936: *thetai = 0.0;
937: #endif
938: }
940: #if defined(PETSC_USE_COMPLEX)
941: if(thetai) *thetai = 0.0;
942: #endif
943: *maxits = data->maxits;
944: *tol = data->tol;
946: return(0);
947: }
950: /**** Patterns implementation *************************************************/
952: typedef PetscInt (*funcV0_t)(dvdDashboard*, PetscInt, PetscInt, Vec*);
953: typedef PetscInt (*funcV1_t)(dvdDashboard*, PetscInt, PetscInt, Vec*,
954: PetscScalar*, Vec);
958: /* Compute D <- K^{-1} * funcV[r_s..r_e] */
959: PetscErrorCode dvd_improvex_PfuncV(dvdDashboard *d, void *funcV, Vec *D,
960: PetscInt max_size_D, PetscInt r_s, PetscInt r_e, Vec *auxV,
961: PetscScalar *auxS)
962: {
963: PetscErrorCode ierr;
964: PetscInt i;
968: if (max_size_D >= r_e-r_s+1) {
969: /* The optimized version needs one vector extra of D */
970: /* D(1:r.size) = R(r_s:r_e-1) */
971: if (auxS) ((funcV1_t)funcV)(d, r_s, r_e, D+1, auxS, auxV[0]);
972: else ((funcV0_t)funcV)(d, r_s, r_e, D+1);
974: /* D = K^{-1} * R */
975: for (i=0; i<r_e-r_s; i++) {
976: d->improvex_precond(d, i+r_s, D[i+1], D[i]);
977: }
978: } else if (max_size_D == r_e-r_s) {
979: /* Non-optimized version */
980: /* auxV <- R[r_e-1] */
981: if (auxS) ((funcV1_t)funcV)(d, r_e-1, r_e, auxV, auxS, auxV[1]);
982: else ((funcV0_t)funcV)(d, r_e-1, r_e, auxV);
984: /* D(1:r.size-1) = R(r_s:r_e-2) */
985: if (auxS) ((funcV1_t)funcV)(d, r_s, r_e-1, D+1, auxS, auxV[1]);
986: else ((funcV0_t)funcV)(d, r_s, r_e-1, D+1);
988: /* D = K^{-1} * R */
989: for (i=0; i<r_e-r_s-1; i++) {
990: d->improvex_precond(d, i+r_s, D[i+1], D[i]);
991: }
992: d->improvex_precond(d, r_e-1, auxV[0], D[r_e-r_s-1]);
993: } else {
994: SETERRQ(PETSC_COMM_SELF,1, "Problem: r_e-r_s > max_size_D!");
995: }
997: return(0);
998: }
1003: /* Compute the left and right projected eigenvectors where,
1004: pX, the returned right eigenvectors
1005: pY, the returned left eigenvectors,
1006: ld_, the leading dimension of pX and pY,
1007: auxS, auxiliar vector of size length 6*d->size_H
1008: */
1009: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
1010: PetscScalar *pY, PetscInt ld, PetscScalar *auxS, PetscInt size_auxS)
1011: {
1012: PetscErrorCode ierr;
1013:
1016: SlepcDenseCopy(pY, ld, d->T?d->pY:d->pX, d->ldpX, d->size_H,
1017: d->size_H);
1018: SlepcDenseCopy(pX, ld, d->pX, d->ldpX, d->size_H, d->size_H);
1019:
1020:
1021: /* [qX, qY] <- eig(S, T); pX <- d->pX * qX; pY <- d->pY * qY */
1022: dvd_compute_eigenvectors(d->size_H, d->S, d->ldS, d->T, d->ldT, pX,
1023: ld, pY, ld, auxS, size_auxS, PETSC_TRUE);
1024:
1026: /* 2-Normalize the columns of pX an pY */
1027: SlepcDenseNorm(pX, ld, d->size_H, d->size_H, d->eigi-d->cX_in_H);
1028: SlepcDenseNorm(pY, ld, d->size_H, d->size_H, d->eigi-d->cX_in_H);
1030: return(0);
1031: }
1035: /* Compute (I - KZ*iXKZ*X')*V where,
1036: V, the vectors to apply the projector,
1037: cV, the number of vectors in V,
1038: auxS, auxiliar vector of size length 3*size_iXKZ*cV
1039: */
1040: PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d, Vec *V, PetscInt cV, PetscScalar *auxS)
1041: {
1042: PetscErrorCode ierr;
1043: dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1044: PetscInt size_in = data->size_iXKZ*cV, i, ldh;
1045: PetscScalar *h, *in, *out;
1046: PetscBLASInt cV_, n, info, ld;
1047: DvdReduction r;
1048: DvdReductionChunk
1049: ops[4];
1050: DvdMult_copy_func
1051: sr[4];
1052:
1055: /* Check consistency */
1056: if (cV > 2) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
1058: /* h <- X'*V */
1059: h = auxS; in = h+size_in; out = in+size_in; ldh = data->size_iXKZ;
1060: SlepcAllReduceSumBegin(ops, 4, in, out, size_in, &r,
1061: ((PetscObject)d->V[0])->comm);
1062: for (i=0; i<cV; i++) {
1063: VecsMultS(&h[i*ldh],0,ldh,d->V-data->size_KZ,0,data->size_KZ,V+i,0,1,&r,&sr[i*2]);
1064: VecsMultS(&h[i*ldh+data->size_KZ],0,ldh,data->u,0,data->size_iXKZ-data->size_KZ,V+i,0,1,&r,&sr[i*2+1]);
1065: }
1066: SlepcAllReduceSumEnd(&r);
1068: /* h <- iXKZ\h */
1069: cV_ = PetscBLASIntCast(cV);
1070: n = PetscBLASIntCast(data->size_iXKZ);
1071: ld = PetscBLASIntCast(data->ldiXKZ);
1073: LAPACKgetrs_("N", &n, &cV_, data->iXKZ, &ld, data->iXKZPivots, h, &n, &info);
1074: if (info) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_LIB, "Error in Lapack XGETRS %d", info);
1076: /* V <- V - KZ*h */
1077: for (i=0; i<cV; i++) {
1078: SlepcUpdateVectorsZ(V+i,1.0,-1.0,data->KZ,data->size_iXKZ,&h[ldh*i],ldh,data->size_iXKZ,1);
1079: }
1081: return(0);
1082: }