Actual source code: dvd_utils.c
1: /*
2: SLEPc eigensolver: "davidson"
4: Some utils
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
28: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
29: Vec Px);
30: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px);
31: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px);
32: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d);
34: typedef struct {
35: PC pc;
36: } dvdPCWrapper;
37: /*
38: Create a static preconditioner from a PC
39: */
42: PetscErrorCode dvd_static_precond_PC(dvdDashboard *d, dvdBlackboard *b, PC pc)
43: {
44: PetscErrorCode ierr;
45: dvdPCWrapper *dvdpc;
46: Mat P;
47: MatStructure str;
51: /* Setup the step */
52: if (b->state >= DVD_STATE_CONF) {
53: /* If the preconditioner is valid */
54: if (pc) {
55: PetscMalloc(sizeof(dvdPCWrapper), &dvdpc);
56: dvdpc->pc = pc;
57: PetscObjectReference((PetscObject)pc);
58: d->improvex_precond_data = dvdpc;
59: d->improvex_precond = dvd_static_precond_PC_0;
61: /* PC saves the matrix associated with the linear system, and it has to
62: be initialize to a valid matrix */
63: PCGetOperators(pc, PETSC_NULL, &P, &str);
64: if (P) {
65: PetscObjectReference((PetscObject)P);
66: PCSetOperators(pc, P, P, str);
67: MatDestroy(&P);
68: } else
69: d->improvex_precond = dvd_precond_none;
71: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
73: /* Else, use no preconditioner */
74: } else
75: d->improvex_precond = dvd_precond_none;
76: }
78: return(0);
79: }
83: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
84: {
85: PetscErrorCode ierr;
86: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
90: /* Free local data */
91: if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
92: PetscFree(d->improvex_precond_data);
94: return(0);
95: }
100: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
101: Vec Px)
102: {
103: PetscErrorCode ierr;
104: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
108: PCApply(dvdpc->pc, x, Px);
109:
110: return(0);
111: }
113: typedef struct {
114: Vec diagA, diagB;
115: } dvdJacobiPrecond;
119: /*
120: Create the Jacobi preconditioner for Generalized Eigenproblems
121: */
122: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d, dvdBlackboard *b)
123: {
124: PetscErrorCode ierr;
125: dvdJacobiPrecond
126: *dvdjp;
127: PetscBool t;
131: /* Check if the problem matrices support GetDiagonal */
132: MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
133: if (t && d->B) {
134: MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
135: }
137: /* Setting configuration constrains */
138: b->own_vecs+= t?( (d->B == 0)?1:2 ) : 0;
140: /* Setup the step */
141: if (b->state >= DVD_STATE_CONF) {
142: if (t) {
143: PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
144: dvdjp->diagA = *b->free_vecs; b->free_vecs++;
145: MatGetDiagonal(d->A,dvdjp->diagA);
146: if (d->B) {
147: dvdjp->diagB = *b->free_vecs; b->free_vecs++;
148: MatGetDiagonal(d->B, dvdjp->diagB);
149: } else
150: dvdjp->diagB = 0;
151: d->improvex_precond_data = dvdjp;
152: d->improvex_precond = dvd_jacobi_precond_0;
154: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
156: /* Else, use no preconditioner */
157: } else
158: d->improvex_precond = dvd_precond_none;
159: }
161: return(0);
162: }
166: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
167: {
168: PetscErrorCode ierr;
169: dvdJacobiPrecond
170: *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
174: /* Compute inv(D - eig)*x */
175: if (dvdjp->diagB == 0) {
176: /* Px <- diagA - l */
177: VecCopy(dvdjp->diagA, Px);
178: VecShift(Px, -d->eigr[i]);
179: } else {
180: /* Px <- diagA - l*diagB */
181: VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
182:
183: }
185: /* Px(i) <- x/Px(i) */
186: VecPointwiseDivide(Px, x, Px);
188: return(0);
189: }
193: /*
194: Create a trivial preconditioner
195: */
196: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
197: {
198: PetscErrorCode ierr;
202: VecCopy(x, Px);
204: return(0);
205: }
208: /*
209: Use of PETSc profiler functions
210: */
212: /* Define stages */
213: #define DVD_STAGE_INITV 0
214: #define DVD_STAGE_NEWITER 1
215: #define DVD_STAGE_CALCPAIRS 2
216: #define DVD_STAGE_IMPROVEX 3
217: #define DVD_STAGE_UPDATEV 4
218: #define DVD_STAGE_ORTHV 5
220: PetscErrorCode dvd_profiler_d(dvdDashboard *d);
222: typedef struct {
223: PetscErrorCode (*old_initV)(struct _dvdDashboard*);
224: PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
225: PetscErrorCode (*old_improveX)(struct _dvdDashboard*, Vec *D,
226: PetscInt max_size_D, PetscInt r_s,
227: PetscInt r_e, PetscInt *size_D);
228: PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
229: PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
230: } DvdProfiler;
232: PetscLogStage stages[6] = {0,0,0,0,0,0};
234: /*** Other things ****/
238: PetscErrorCode dvd_prof_init()
239: {
240: PetscErrorCode ierr;
243: if (!stages[0]) {
244: PetscLogStageRegister("Dvd_step_initV", &stages[DVD_STAGE_INITV]);
245:
246: PetscLogStageRegister("Dvd_step_calcPairs",&stages[DVD_STAGE_CALCPAIRS]);
247: PetscLogStageRegister("Dvd_step_improveX",&stages[DVD_STAGE_IMPROVEX]);
248: PetscLogStageRegister("Dvd_step_updateV",&stages[DVD_STAGE_UPDATEV]);
249: PetscLogStageRegister("Dvd_step_orthV",&stages[DVD_STAGE_ORTHV]);
250: }
251: return(0);
252: }
256: PetscErrorCode dvd_initV_prof(dvdDashboard* d)
257: {
258: DvdProfiler *p = (DvdProfiler*)d->prof_data;
259: PetscErrorCode ierr;
263: PetscLogStagePush(stages[DVD_STAGE_INITV]);
264: p->old_initV(d);
265: PetscLogStagePop();
267: return(0);
268: }
272: PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d)
273: {
274: DvdProfiler *p = (DvdProfiler*)d->prof_data;
275: PetscErrorCode ierr;
279: PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
280: p->old_calcPairs(d);
281: PetscLogStagePop();
283: return(0);
284: }
288: PetscErrorCode dvd_improveX_prof(dvdDashboard* d, Vec *D, PetscInt max_size_D,
289: PetscInt r_s, PetscInt r_e, PetscInt *size_D)
290: {
291: DvdProfiler *p = (DvdProfiler*)d->prof_data;
292: PetscErrorCode ierr;
296: PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
297: p->old_improveX(d, D, max_size_D, r_s, r_e, size_D);
298: PetscLogStagePop();
300: return(0);
301: }
305: PetscErrorCode dvd_updateV_prof(dvdDashboard *d)
306: {
307: DvdProfiler *p = (DvdProfiler*)d->prof_data;
308: PetscErrorCode ierr;
312: PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
313: p->old_updateV(d);
314: PetscLogStagePop();
316: return(0);
317: }
321: PetscErrorCode dvd_orthV_prof(dvdDashboard *d)
322: {
323: DvdProfiler *p = (DvdProfiler*)d->prof_data;
324: PetscErrorCode ierr;
328: PetscLogStagePush(stages[DVD_STAGE_ORTHV]);
329: p->old_orthV(d);
330: PetscLogStagePop();
332: return(0);
333: }
337: PetscErrorCode dvd_profiler(dvdDashboard *d, dvdBlackboard *b)
338: {
339: PetscErrorCode ierr;
340: DvdProfiler *p;
344: /* Setup the step */
345: if (b->state >= DVD_STATE_CONF) {
346: PetscFree(d->prof_data);
347: PetscMalloc(sizeof(DvdProfiler), &p);
348: d->prof_data = p;
349: p->old_initV = d->initV; d->initV = dvd_initV_prof;
350: p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
351: p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
352: p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;
354: DVD_FL_ADD(d->destroyList, dvd_profiler_d);
355: }
357: return(0);
358: }
362: PetscErrorCode dvd_profiler_d(dvdDashboard *d)
363: {
364: PetscErrorCode ierr;
365: DvdProfiler *p = (DvdProfiler*)d->prof_data;
369: /* Free local data */
370: PetscFree(p);
372: return(0);
373: }
377: /*
378: Configure the harmonics.
379: switch(mode) {
380: DVD_HARM_RR: harmonic RR
381: DVD_HARM_RRR: relative harmonic RR
382: DVD_HARM_REIGS: rightmost eigenvalues
383: DVD_HARM_LEIGS: largest eigenvalues
384: }
385: fixedTarged, if true use the target instead of the best eigenvalue
386: target, the fixed target to be used
387: */
388: typedef struct {
389: PetscScalar
390: Wa, Wb, /* span{W} = span{Wa*AV - Wb*BV} */
391: Pa, Pb; /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
392: PetscBool
393: withTarget;
394: HarmType_t
395: mode;
397: /* old values of eps */
398: EPSWhich
399: old_which;
400: PetscErrorCode
401: (*old_which_func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,
402: PetscInt*,void*);
403: void
404: *old_which_ctx;
405: } dvdHarmonic;
407: PetscErrorCode dvd_harm_start(dvdDashboard *d);
408: PetscErrorCode dvd_harm_end(dvdDashboard *d);
409: PetscErrorCode dvd_harm_d(dvdDashboard *d);
410: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t);
411: PetscErrorCode dvd_harm_updateW(dvdDashboard *d);
412: PetscErrorCode dvd_harm_proj(dvdDashboard *d);
413: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
414: PetscScalar br, PetscScalar bi, PetscInt *r,
415: void *ctx);
416: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d);
420: PetscErrorCode dvd_harm_conf(dvdDashboard *d, dvdBlackboard *b,
421: HarmType_t mode, PetscBool fixedTarget,
422: PetscScalar t)
423: {
424: PetscErrorCode ierr;
425: dvdHarmonic *dvdh;
429: /* Set the problem to GNHEP:
430: d->G maybe is upper triangular due to biorthogonality of V and W */
431: d->sEP = d->sA = d->sB = 0;
433: /* Setup the step */
434: if (b->state >= DVD_STATE_CONF) {
435: PetscMalloc(sizeof(dvdHarmonic), &dvdh);
436: dvdh->withTarget = fixedTarget;
437: dvdh->mode = mode;
438: if (fixedTarget) dvd_harm_transf(dvdh, t);
439: d->calcpairs_W_data = dvdh;
440: d->calcpairs_W = dvd_harm_updateW;
441: d->calcpairs_proj_trans = dvd_harm_proj;
442: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
444: DVD_FL_ADD(d->startList, dvd_harm_start);
445: DVD_FL_ADD(d->endList, dvd_harm_end);
446: DVD_FL_ADD(d->destroyList, dvd_harm_d);
447: }
449: return(0);
450: }
455: PetscErrorCode dvd_harm_d(dvdDashboard *d)
456: {
457: PetscErrorCode ierr;
461: /* Free local data */
462: PetscFree(d->calcpairs_W_data);
464: return(0);
465: }
470: PetscErrorCode dvd_harm_start(dvdDashboard *d)
471: {
472: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
476: /* Overload the eigenpairs selection routine */
477: data->old_which = d->eps->which;
478: data->old_which_func = d->eps->which_func;
479: data->old_which_ctx = d->eps->which_ctx;
480: d->eps->which = EPS_WHICH_USER;
481: d->eps->which_func = dvd_harm_sort;
482: d->eps->which_ctx = data;
484: return(0);
485: }
490: PetscErrorCode dvd_harm_end(dvdDashboard *d)
491: {
492: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
496: /* Restore the eigenpairs selection routine */
497: d->eps->which = data->old_which;
498: d->eps->which_func = data->old_which_func;
499: d->eps->which_ctx = data->old_which_ctx;
501: return(0);
502: }
507: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t)
508: {
511: switch(dvdh->mode) {
512: case DVD_HARM_RR: /* harmonic RR */
513: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
514: case DVD_HARM_RRR: /* relative harmonic RR */
515: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
516: case DVD_HARM_REIGS: /* rightmost eigenvalues */
517: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
518: break;
519: case DVD_HARM_LEIGS: /* largest eigenvalues */
520: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
521: case DVD_HARM_NONE:
522: default:
523: SETERRQ(PETSC_COMM_WORLD,1, "Harmonic type not supported");
524: }
526: /* Check the transformation does not change the sign of the imaginary part */
527: #if !defined(PETSC_USE_COMPLEX)
528: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
529: dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
530: #endif
532: return(0);
533: }
537: PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
538: {
539: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
540: PetscErrorCode ierr;
541: PetscInt i;
545: /* Update the target if it is necessary */
546: if (!data->withTarget) dvd_harm_transf(data, d->eigr[0]);
547:
548: for(i=d->V_new_s; i<d->V_new_e; i++) {
549: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
550: VecCopy(d->AV[i], d->W[i]);
551: VecAXPBY(d->W[i], -data->Wb, data->Wa, (d->BV?d->BV:d->V)[i]);
552:
553: }
555: return(0);
556: }
560: PetscErrorCode dvd_harm_proj(dvdDashboard *d)
561: {
562: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
563: PetscInt i,j;
567: if (d->sH != d->sG) {
568: SETERRQ(PETSC_COMM_WORLD,1,"Projected matrices H and G must have the same structure");
569: PetscFunctionReturn(1);
570: }
572: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
573: if (DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)) /* Upper triangular part */
574: for(i=d->V_new_s+d->cX_in_H; i<d->V_new_e+d->cX_in_H; i++)
575: for(j=0; j<=i; j++) {
576: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
577: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
578: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
579: }
580: if (DVD_ISNOT(d->sH,DVD_MAT_UTRIANG)) /* Lower triangular part */
581: for(i=0; i<d->V_new_e+d->cX_in_H; i++)
582: for(j=PetscMax(d->V_new_s+d->cX_in_H,i+(DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)?1:0));
583: j<d->V_new_e+d->cX_in_H; j++) {
584: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
585: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
586: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
587: }
589: return(0);
590: }
594: PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data, PetscScalar *ar,
595: PetscScalar *ai)
596: {
597: PetscScalar xr;
598: #if !defined(PETSC_USE_COMPLEX)
599: PetscScalar xi, k;
600: #endif
604: xr = *ar;
605: #if !defined(PETSC_USE_COMPLEX)
607: xi = *ai;
609: if (xi != 0.0) {
610: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
611: data->Wa*data->Wa*xi*xi;
612: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
613: data->Wb*data->Wa*(xr*xr + xi*xi))/k;
614: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
615: } else
616: #endif
617: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
619: return(0);
620: }
625: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
626: PetscScalar br, PetscScalar bi, PetscInt *r,
627: void *ctx)
628: {
629: dvdHarmonic *data = (dvdHarmonic*)ctx;
630: PetscErrorCode ierr;
634: /* Back-transform the harmonic values */
635: dvd_harm_backtrans(data, &ar, &ai);
636: dvd_harm_backtrans(data, &br, &bi);
638: /* Compare values using the user options for the eigenpairs selection */
639: eps->which = data->old_which;
640: eps->which_func = data->old_which_func;
641: eps->which_ctx = data->old_which_ctx;
642: EPSCompareEigenvalues(eps, ar, ai, br, bi, r);
644: /* Restore the eps values */
645: eps->which = EPS_WHICH_USER;
646: eps->which_func = dvd_harm_sort;
647: eps->which_ctx = data;
649: return(0);
650: }
654: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
655: {
656: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
657: PetscInt i;
661: for(i=0; i<d->size_H; i++)
662: dvd_harm_backtrans(data, &d->eigr[i-d->cX_in_H], &d->eigi[i-d->cX_in_H]);
664: return(0);
665: }