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-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
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: PetscObjectReference((PetscObject)P);
65: PCSetOperators(pc, P, P, str);
66: MatDestroy(P);
67: PCSetUp(pc);
69: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
71: /* Else, use no preconditioner */
72: } else
73: d->improvex_precond = dvd_precond_none;
74: }
76: return(0);
77: }
81: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
82: {
83: PetscErrorCode ierr;
84: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
88: /* Free local data */
89: if (dvdpc->pc) { PCDestroy(dvdpc->pc); }
90: PetscFree(d->improvex_precond_data);
91: d->improvex_precond_data = PETSC_NULL;
93: return(0);
94: }
97: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
98: Vec Px)
99: {
100: PetscErrorCode ierr;
101: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
105: PCApply(dvdpc->pc, x, Px);
106:
107: return(0);
108: }
110: typedef struct {
111: Vec diagA, diagB;
112: } dvdJacobiPrecond;
113: /*
114: Create the Jacobi preconditioner for Generalized Eigenproblems
115: */
116: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d, dvdBlackboard *b)
117: {
118: PetscErrorCode ierr;
119: dvdJacobiPrecond
120: *dvdjp;
121: PetscTruth t;
125: /* Check if the problem matrices support GetDiagonal */
126: MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
127: if (t && d->B) {
128: MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
129: }
131: /* Setting configuration constrains */
132: b->own_vecs+= t?( (d->B == 0)?1:2 ) : 0;
134: /* Setup the step */
135: if (b->state >= DVD_STATE_CONF) {
136: if (t) {
137: PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
138: dvdjp->diagA = *b->free_vecs; b->free_vecs++;
139: MatGetDiagonal(d->A,dvdjp->diagA);
140: if (d->B) {
141: dvdjp->diagB = *b->free_vecs; b->free_vecs++;
142: MatGetDiagonal(d->B, dvdjp->diagB);
143: } else
144: dvdjp->diagB = 0;
145: d->improvex_precond_data = dvdjp;
146: d->improvex_precond = dvd_jacobi_precond_0;
148: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
150: /* Else, use no preconditioner */
151: } else
152: d->improvex_precond = dvd_precond_none;
153: }
155: return(0);
156: }
158: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
159: {
160: PetscErrorCode ierr;
161: dvdJacobiPrecond
162: *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
166: /* Compute inv(D - eig)*x */
167: if (dvdjp->diagB == 0) {
168: /* Px <- diagA - l */
169: VecCopy(dvdjp->diagA, Px);
170: VecShift(Px, -d->eigr[i]);
171: } else {
172: /* Px <- diagA - l*diagB */
173: VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
174:
175: }
177: /* Px(i) <- x/Px(i) */
178: VecPointwiseDivide(Px, x, Px);
180: return(0);
181: }
183: /*
184: Create a trivial preconditioner
185: */
186: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
187: {
188: PetscErrorCode ierr;
192: VecCopy(x, Px);
194: return(0);
195: }
198: /*
199: Use of PETSc profiler functions
200: */
202: /* Define stages */
203: #define DVD_STAGE_INITV 0
204: #define DVD_STAGE_NEWITER 1
205: #define DVD_STAGE_CALCPAIRS 2
206: #define DVD_STAGE_IMPROVEX 3
207: #define DVD_STAGE_UPDATEV 4
208: #define DVD_STAGE_ORTHV 5
210: PetscErrorCode dvd_profiler_d(dvdDashboard *d);
212: typedef struct {
213: PetscErrorCode (*old_initV)(struct _dvdDashboard*);
214: PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
215: PetscErrorCode (*old_improveX)(struct _dvdDashboard*, Vec *D,
216: PetscInt max_size_D, PetscInt r_s,
217: PetscInt r_e, PetscInt *size_D);
218: PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
219: PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
220: } DvdProfiler;
222: PetscLogStage stages[6] = {0,0,0,0,0,0};
224: /*** Other things ****/
228: PetscErrorCode dvd_prof_init() {
229: PetscErrorCode ierr;
233: if (!stages[0]) {
234: PetscLogStageRegister("Dvd_step_initV", &stages[DVD_STAGE_INITV]);
235:
236: PetscLogStageRegister("Dvd_step_calcPairs",
237: &stages[DVD_STAGE_CALCPAIRS]);
238: PetscLogStageRegister("Dvd_step_improveX",
239: &stages[DVD_STAGE_IMPROVEX]);
240: PetscLogStageRegister("Dvd_step_updateV",
241: &stages[DVD_STAGE_UPDATEV]);
242: PetscLogStageRegister("Dvd_step_orthV",
243: &stages[DVD_STAGE_ORTHV]);
244: }
245: dvd_blas_prof_init();
246:
247: return(0);
248: }
250: PetscErrorCode dvd_initV_prof(dvdDashboard* d) {
251: DvdProfiler *p = (DvdProfiler*)d->prof_data;
252: PetscErrorCode ierr;
256: PetscLogStagePush(stages[DVD_STAGE_INITV]);
257: p->old_initV(d);
258: PetscLogStagePop();
260: return(0);
261: }
263: PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d) {
264: DvdProfiler *p = (DvdProfiler*)d->prof_data;
265: PetscErrorCode ierr;
269: PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
270: p->old_calcPairs(d);
271: PetscLogStagePop();
273: return(0);
274: }
276: PetscErrorCode dvd_improveX_prof(dvdDashboard* d, Vec *D, PetscInt max_size_D,
277: PetscInt r_s, PetscInt r_e, PetscInt *size_D) {
278: DvdProfiler *p = (DvdProfiler*)d->prof_data;
279: PetscErrorCode ierr;
283: PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
284: p->old_improveX(d, D, max_size_D, r_s, r_e, size_D);
285: PetscLogStagePop();
287: return(0);
288: }
290: PetscErrorCode dvd_updateV_prof(dvdDashboard *d) {
291: DvdProfiler *p = (DvdProfiler*)d->prof_data;
292: PetscErrorCode ierr;
296: PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
297: p->old_updateV(d);
298: PetscLogStagePop();
300: return(0);
301: }
303: PetscErrorCode dvd_orthV_prof(dvdDashboard *d) {
304: DvdProfiler *p = (DvdProfiler*)d->prof_data;
305: PetscErrorCode ierr;
309: PetscLogStagePush(stages[DVD_STAGE_ORTHV]);
310: p->old_orthV(d);
311: PetscLogStagePop();
313: return(0);
314: }
318: PetscErrorCode dvd_profiler(dvdDashboard *d, dvdBlackboard *b)
319: {
320: PetscErrorCode ierr;
321: DvdProfiler *p;
325: /* Setup the step */
326: if (b->state >= DVD_STATE_CONF) {
327: if (d->prof_data) {
328: PetscFree(d->prof_data);
329: }
330: PetscMalloc(sizeof(DvdProfiler), &p);
331: d->prof_data = p;
332: p->old_initV = d->initV; d->initV = dvd_initV_prof;
333: p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
334: p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
335: p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;
337: DVD_FL_ADD(d->destroyList, dvd_profiler_d);
338: }
340: return(0);
341: }
345: PetscErrorCode dvd_profiler_d(dvdDashboard *d)
346: {
347: PetscErrorCode ierr;
348: DvdProfiler *p = (DvdProfiler*)d->prof_data;
352: /* Free local data */
353: PetscFree(p);
354: d->prof_data = PETSC_NULL;
356: return(0);
357: }
361: /*
362: Configure the harmonics.
363: switch(mode) {
364: DVD_HARM_RR: harmonic RR
365: DVD_HARM_RRR: relative harmonic RR
366: DVD_HARM_REIGS: rightmost eigenvalues
367: DVD_HARM_LEIGS: largest eigenvalues
368: }
369: fixedTarged, if true use the target instead of the best eigenvalue
370: target, the fixed target to be used
371: */
372: typedef struct {
373: PetscScalar
374: Wa, Wb, /* span{W} = span{Wa*AV - Wb*BV} */
375: Pa, Pb; /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
376: PetscTruth
377: withTarget;
378: HarmType_t
379: mode;
381: /* old values of eps */
382: EPSWhich
383: old_which;
384: PetscErrorCode
385: (*old_which_func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,
386: PetscInt*,void*);
387: void
388: *old_which_ctx;
389: } dvdHarmonic;
391: PetscErrorCode dvd_harm_start(dvdDashboard *d);
392: PetscErrorCode dvd_harm_end(dvdDashboard *d);
393: PetscErrorCode dvd_harm_d(dvdDashboard *d);
394: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t);
395: PetscErrorCode dvd_harm_updateW(dvdDashboard *d);
396: PetscErrorCode dvd_harm_proj(dvdDashboard *d);
397: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
398: PetscScalar br, PetscScalar bi, PetscInt *r,
399: void *ctx);
400: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d);
402: PetscErrorCode dvd_harm_conf(dvdDashboard *d, dvdBlackboard *b,
403: HarmType_t mode, PetscTruth fixedTarget,
404: PetscScalar t)
405: {
406: PetscErrorCode ierr;
407: dvdHarmonic *dvdh;
411: /* Set the problem to GNHEP */
412: // TODO: d->G maybe is upper triangular due to biorthogonality of V and W
413: d->sEP = d->sA = d->sB = 0;
415: /* Setup the step */
416: if (b->state >= DVD_STATE_CONF) {
417: PetscMalloc(sizeof(dvdHarmonic), &dvdh);
418: dvdh->withTarget = fixedTarget;
419: dvdh->mode = mode;
420: if (fixedTarget) dvd_harm_transf(dvdh, t);
421: d->calcpairs_W_data = dvdh;
422: d->calcpairs_W = dvd_harm_updateW;
423: d->calcpairs_proj_trans = dvd_harm_proj;
424: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
426: DVD_FL_ADD(d->startList, dvd_harm_start);
427: DVD_FL_ADD(d->endList, dvd_harm_end);
428: DVD_FL_ADD(d->destroyList, dvd_harm_d);
429: }
431: return(0);
432: }
437: PetscErrorCode dvd_harm_d(dvdDashboard *d)
438: {
439: PetscErrorCode ierr;
443: /* Free local data */
444: PetscFree(d->calcpairs_W_data);
445: d->calcpairs_W_data = PETSC_NULL;
447: return(0);
448: }
453: PetscErrorCode dvd_harm_start(dvdDashboard *d)
454: {
455: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
459: /* Overload the eigenpairs selection routine */
460: data->old_which = d->eps->which;
461: data->old_which_func = d->eps->which_func;
462: data->old_which_ctx = d->eps->which_ctx;
463: d->eps->which = EPS_WHICH_USER;
464: d->eps->which_func = dvd_harm_sort;
465: d->eps->which_ctx = data;
467: return(0);
468: }
473: PetscErrorCode dvd_harm_end(dvdDashboard *d)
474: {
475: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
479: /* Restore the eigenpairs selection routine */
480: d->eps->which = data->old_which;
481: d->eps->which_func = data->old_which_func;
482: d->eps->which_ctx = data->old_which_ctx;
484: return(0);
485: }
488: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t)
489: {
492: switch(dvdh->mode) {
493: case DVD_HARM_RR: /* harmonic RR */
494: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
495: case DVD_HARM_RRR: /* relative harmonic RR */
496: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
497: case DVD_HARM_REIGS: /* rightmost eigenvalues */
498: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
499: break;
500: case DVD_HARM_LEIGS: /* largest eigenvalues */
501: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
502: case DVD_HARM_NONE:
503: default:
504: SETERRQ(1, "The harmonic type is not supported!");
505: }
507: /* Check the transformation does not change the sign of the imaginary part */
508: #if !defined(PETSC_USE_COMPLEX)
509: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
510: dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
511: #endif
513: return(0);
514: }
516: PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
517: {
518: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
519: PetscErrorCode ierr;
520: PetscInt i;
524: /* Update the target if it is necessary */
525: if (!data->withTarget) dvd_harm_transf(data, d->eigr[0]);
526:
527: for(i=d->V_new_s; i<d->V_new_e; i++) {
528: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
529: VecCopy(d->AV[i], d->W[i]);
530: VecAXPBY(d->W[i], -data->Wb, data->Wa, (d->BV?d->BV:d->V)[i]);
531:
532: }
534: return(0);
535: }
537: PetscErrorCode dvd_harm_proj(dvdDashboard *d)
538: {
539: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
540: PetscInt i,j;
544: if (d->sH != d->sG) {
545: SETERRQ(1, "Error: Projected matrices H and G must have the same structure!");
546: PetscFunctionReturn(1);
547: }
549: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
550: if (DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)) /* Upper triangular part */
551: for(i=d->V_new_s; i<d->V_new_e; i++)
552: for(j=0; j<=i; j++) {
553: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
554: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
555: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
556: }
557: if (DVD_ISNOT(d->sH,DVD_MAT_UTRIANG)) /* Lower triangular part */
558: for(i=0; i<d->V_new_e; i++)
559: for(j=PetscMax(d->V_new_s,i+(DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)?1:0));
560: j<d->V_new_e; j++) {
561: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
562: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
563: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
564: }
566: return(0);
567: }
569: PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data, PetscScalar *ar,
570: PetscScalar *ai)
571: {
572: PetscScalar xr;
573: #if !defined(PETSC_USE_COMPLEX)
574: PetscScalar xi, k;
575: #endif
579: if(!ar) SETERRQ(1, "The real part has to be present!");
580: xr = *ar;
582: #if !defined(PETSC_USE_COMPLEX)
583: if(!ai) SETERRQ(1, "The imaginary part has to be present!");
584: xi = *ai;
586: if (xi != 0.0) {
587: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
588: data->Wa*data->Wa*xi*xi;
589: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
590: data->Wb*data->Wa*(xr*xr + xi*xi))/k;
591: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
592: } else
593: #endif
594: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
596: return(0);
597: }
600: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
601: PetscScalar br, PetscScalar bi, PetscInt *r,
602: void *ctx)
603: {
604: dvdHarmonic *data = (dvdHarmonic*)ctx;
605: PetscErrorCode ierr;
609: /* Back-transform the harmonic values */
610: dvd_harm_backtrans(data, &ar, &ai);
611: dvd_harm_backtrans(data, &br, &bi);
613: /* Compare values using the user options for the eigenpairs selection */
614: eps->which = data->old_which;
615: eps->which_func = data->old_which_func;
616: eps->which_ctx = data->old_which_ctx;
617: EPSCompareEigenvalues(eps, ar, ai, br, bi, r);
619: /* Restore the eps values */
620: eps->which = EPS_WHICH_USER;
621: eps->which_func = dvd_harm_sort;
622: eps->which_ctx = data;
624: return(0);
625: }
627: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
628: {
629: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
630: PetscInt i;
634: for(i=0; i<d->size_H; i++)
635: dvd_harm_backtrans(data, &d->eigr[i], &d->eigi[i]);
637: return(0);
638: }