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: PCSetOperators(pc, P, P, str);
65: PCSetUp(pc);
67: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
69: /* Else, use no preconditioner */
70: } else
71: d->improvex_precond = dvd_precond_none;
72: }
74: return(0);
75: }
79: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
80: {
81: PetscErrorCode ierr;
82: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
86: /* Free local data */
87: if (dvdpc->pc) { PCDestroy(dvdpc->pc); }
88: PetscFree(d->improvex_precond_data);
89: d->improvex_precond_data = PETSC_NULL;
91: return(0);
92: }
95: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
96: Vec Px)
97: {
98: PetscErrorCode ierr;
99: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
103: PCApply(dvdpc->pc, x, Px);
104:
105: return(0);
106: }
108: typedef struct {
109: Vec diagA, diagB;
110: } dvdJacobiPrecond;
111: /*
112: Create the Jacobi preconditioner for Generalized Eigenproblems
113: */
114: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d, dvdBlackboard *b)
115: {
116: PetscErrorCode ierr;
117: dvdJacobiPrecond
118: *dvdjp;
119: PetscTruth t;
123: /* Check if the problem matrices support GetDiagonal */
124: MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
125: if ((t == PETSC_TRUE) && d->B) {
126: MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
127: }
129: /* Setting configuration constrains */
130: b->own_vecs+= t==PETSC_TRUE?( (d->B == 0)?1:2 ) : 0;
132: /* Setup the step */
133: if (b->state >= DVD_STATE_CONF) {
134: if (t == PETSC_TRUE) {
135: PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
136: dvdjp->diagA = *b->free_vecs; b->free_vecs++;
137: MatGetDiagonal(d->A,dvdjp->diagA);
138: if (d->B) {
139: dvdjp->diagB = *b->free_vecs; b->free_vecs++;
140: MatGetDiagonal(d->B, dvdjp->diagB);
141: } else
142: dvdjp->diagB = 0;
143: d->improvex_precond_data = dvdjp;
144: d->improvex_precond = dvd_jacobi_precond_0;
146: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
148: /* Else, use no preconditioner */
149: } else
150: d->improvex_precond = dvd_precond_none;
151: }
153: return(0);
154: }
156: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
157: {
158: PetscErrorCode ierr;
159: dvdJacobiPrecond
160: *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
164: /* Compute inv(D - eig)*x */
165: if (dvdjp->diagB == 0) {
166: /* Px <- diagA - l */
167: VecCopy(dvdjp->diagA, Px);
168: VecShift(Px, -d->eigr[i]);
169: } else {
170: /* Px <- diagA - l*diagB */
171: VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
172:
173: }
175: /* Px(i) <- x/Px(i) */
176: VecPointwiseDivide(Px, x, Px);
178: return(0);
179: }
181: /*
182: Create a trivial preconditioner
183: */
184: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
185: {
186: PetscErrorCode ierr;
190: VecCopy(x, Px);
192: return(0);
193: }
196: /*
197: Use of PETSc profiler functions
198: */
200: /* Define stages */
201: #define DVD_STAGE_INITV 0
202: #define DVD_STAGE_NEWITER 1
203: #define DVD_STAGE_CALCPAIRS 2
204: #define DVD_STAGE_IMPROVEX 3
205: #define DVD_STAGE_UPDATEV 4
206: #define DVD_STAGE_ORTHV 5
208: PetscErrorCode dvd_profiler_d(dvdDashboard *d);
210: typedef struct {
211: PetscErrorCode (*old_initV)(struct _dvdDashboard*);
212: PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
213: PetscErrorCode (*old_improveX)(struct _dvdDashboard*, Vec *D,
214: PetscInt max_size_D, PetscInt r_s,
215: PetscInt r_e, PetscInt *size_D);
216: PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
217: PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
218: } DvdProfiler;
220: PetscLogStage stages[6] = {0,0,0,0,0,0};
222: /*** Other things ****/
226: PetscErrorCode dvd_prof_init() {
227: PetscErrorCode ierr;
231: if (!stages[0]) {
232: PetscLogStageRegister("Dvd_step_initV", &stages[DVD_STAGE_INITV]);
233:
234: PetscLogStageRegister("Dvd_step_calcPairs",
235: &stages[DVD_STAGE_CALCPAIRS]);
236: PetscLogStageRegister("Dvd_step_improveX",
237: &stages[DVD_STAGE_IMPROVEX]);
238: PetscLogStageRegister("Dvd_step_updateV",
239: &stages[DVD_STAGE_UPDATEV]);
240: PetscLogStageRegister("Dvd_step_orthV",
241: &stages[DVD_STAGE_ORTHV]);
242: }
243: dvd_blas_prof_init();
244:
245: return(0);
246: }
248: PetscErrorCode dvd_initV_prof(dvdDashboard* d) {
249: DvdProfiler *p = (DvdProfiler*)d->prof_data;
250: PetscErrorCode ierr;
254: PetscLogStagePush(stages[DVD_STAGE_INITV]);
255: p->old_initV(d);
256: PetscLogStagePop();
258: return(0);
259: }
261: PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d) {
262: DvdProfiler *p = (DvdProfiler*)d->prof_data;
263: PetscErrorCode ierr;
267: PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
268: p->old_calcPairs(d);
269: PetscLogStagePop();
271: return(0);
272: }
274: PetscErrorCode dvd_improveX_prof(dvdDashboard* d, Vec *D, PetscInt max_size_D,
275: PetscInt r_s, PetscInt r_e, PetscInt *size_D) {
276: DvdProfiler *p = (DvdProfiler*)d->prof_data;
277: PetscErrorCode ierr;
281: PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
282: p->old_improveX(d, D, max_size_D, r_s, r_e, size_D);
283: PetscLogStagePop();
285: return(0);
286: }
288: PetscErrorCode dvd_updateV_prof(dvdDashboard *d) {
289: DvdProfiler *p = (DvdProfiler*)d->prof_data;
290: PetscErrorCode ierr;
294: PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
295: p->old_updateV(d);
296: PetscLogStagePop();
298: return(0);
299: }
301: PetscErrorCode dvd_orthV_prof(dvdDashboard *d) {
302: DvdProfiler *p = (DvdProfiler*)d->prof_data;
303: PetscErrorCode ierr;
307: PetscLogStagePush(stages[DVD_STAGE_ORTHV]);
308: p->old_orthV(d);
309: PetscLogStagePop();
311: return(0);
312: }
316: PetscErrorCode dvd_profiler(dvdDashboard *d, dvdBlackboard *b)
317: {
318: PetscErrorCode ierr;
319: DvdProfiler *p;
323: /* Setup the step */
324: if (b->state >= DVD_STATE_CONF) {
325: if (d->prof_data) {
326: PetscFree(d->prof_data);
327: }
328: PetscMalloc(sizeof(DvdProfiler), &p);
329: d->prof_data = p;
330: p->old_initV = d->initV; d->initV = dvd_initV_prof;
331: p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
332: p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
333: p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;
335: DVD_FL_ADD(d->destroyList, dvd_profiler_d);
336: }
338: return(0);
339: }
343: PetscErrorCode dvd_profiler_d(dvdDashboard *d)
344: {
345: PetscErrorCode ierr;
346: DvdProfiler *p = (DvdProfiler*)d->prof_data;
350: /* Free local data */
351: PetscFree(p);
352: d->prof_data = PETSC_NULL;
354: return(0);
355: }
359: /*
360: Configure the harmonics.
361: switch(mode) {
362: DVD_HARM_RR: harmonic RR
363: DVD_HARM_RRR: relative harmonic RR
364: DVD_HARM_REIGS: rightmost eigenvalues
365: DVD_HARM_LEIGS: largest eigenvalues
366: }
367: fixedTarged, if true use the target instead of the best eigenvalue
368: target, the fixed target to be used
369: */
370: typedef struct {
371: PetscScalar
372: Wa, Wb, /* span{W} = span{Wa*AV - Wb*BV} */
373: Pa, Pb; /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
374: PetscTruth
375: withTarget;
376: HarmType_t
377: mode;
379: /* old values of eps */
380: EPSWhich
381: old_which;
382: PetscErrorCode
383: (*old_which_func)(EPS,PetscScalar,PetscScalar,PetscScalar,PetscScalar,
384: PetscInt*,void*);
385: void
386: *old_which_ctx;
387: } dvdHarmonic;
389: PetscErrorCode dvd_harm_start(dvdDashboard *d);
390: PetscErrorCode dvd_harm_end(dvdDashboard *d);
391: PetscErrorCode dvd_harm_d(dvdDashboard *d);
392: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t);
393: PetscErrorCode dvd_harm_updateW(dvdDashboard *d);
394: PetscErrorCode dvd_harm_proj(dvdDashboard *d);
395: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
396: PetscScalar br, PetscScalar bi, PetscInt *r,
397: void *ctx);
398: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d);
400: PetscErrorCode dvd_harm_conf(dvdDashboard *d, dvdBlackboard *b,
401: HarmType_t mode, PetscTruth fixedTarget,
402: PetscScalar t)
403: {
404: PetscErrorCode ierr;
405: dvdHarmonic *dvdh;
409: /* Set the problem to GNHEP */
410: // TODO: d->G maybe is upper triangular due to biorthogonality of V and W
411: d->sEP = d->sA = d->sB = 0;
413: /* Setup the step */
414: if (b->state >= DVD_STATE_CONF) {
415: PetscMalloc(sizeof(dvdHarmonic), &dvdh);
416: dvdh->withTarget = fixedTarget;
417: dvdh->mode = mode;
418: if (fixedTarget == PETSC_TRUE) dvd_harm_transf(dvdh, t);
419: d->calcpairs_W_data = dvdh;
420: d->calcpairs_W = dvd_harm_updateW;
421: d->calcpairs_proj_trans = dvd_harm_proj;
422: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
424: DVD_FL_ADD(d->startList, dvd_harm_start);
425: DVD_FL_ADD(d->endList, dvd_harm_end);
426: DVD_FL_ADD(d->destroyList, dvd_harm_d);
427: }
429: return(0);
430: }
435: PetscErrorCode dvd_harm_d(dvdDashboard *d)
436: {
437: PetscErrorCode ierr;
441: /* Free local data */
442: PetscFree(d->calcpairs_W_data);
443: d->calcpairs_W_data = PETSC_NULL;
445: return(0);
446: }
451: PetscErrorCode dvd_harm_start(dvdDashboard *d)
452: {
453: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
457: /* Overload the eigenpairs selection routine */
458: data->old_which = d->eps->which;
459: data->old_which_func = d->eps->which_func;
460: data->old_which_ctx = d->eps->which_ctx;
461: d->eps->which = EPS_WHICH_USER;
462: d->eps->which_func = dvd_harm_sort;
463: d->eps->which_ctx = data;
465: return(0);
466: }
471: PetscErrorCode dvd_harm_end(dvdDashboard *d)
472: {
473: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
477: /* Restore the eigenpairs selection routine */
478: d->eps->which = data->old_which;
479: d->eps->which_func = data->old_which_func;
480: d->eps->which_ctx = data->old_which_ctx;
482: return(0);
483: }
486: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t)
487: {
490: switch(dvdh->mode) {
491: case DVD_HARM_RR: /* harmonic RR */
492: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
493: case DVD_HARM_RRR: /* relative harmonic RR */
494: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
495: case DVD_HARM_REIGS: /* rightmost eigenvalues */
496: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
497: break;
498: case DVD_HARM_LEIGS: /* largest eigenvalues */
499: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
500: case DVD_HARM_NONE:
501: default:
502: SETERRQ(1, "The harmonic type is not supported!");
503: }
505: /* Check the transformation does not change the sign of the imaginary part */
506: #if !defined(PETSC_USE_COMPLEX)
507: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
508: dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
509: #endif
511: return(0);
512: }
514: PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
515: {
516: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
517: PetscErrorCode ierr;
518: PetscInt i;
522: /* Update the target if it is necessary */
523: if (data->withTarget == PETSC_FALSE) dvd_harm_transf(data, d->eigr[0]);
524:
525: for(i=d->V_new_s; i<d->V_new_e; i++) {
526: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
527: VecCopy(d->AV[i], d->W[i]);
528: VecAXPBY(d->W[i], -data->Wb, data->Wa, (d->BV?d->BV:d->V)[i]);
529:
530: }
532: return(0);
533: }
535: PetscErrorCode dvd_harm_proj(dvdDashboard *d)
536: {
537: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
538: PetscInt i,j;
542: if (d->sH != d->sG) {
543: SETERRQ(1, "Error: Projected matrices H and G must have the same structure!");
544: PetscFunctionReturn(1);
545: }
547: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
548: if (DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)) /* Upper triangular part */
549: for(i=d->V_new_s; i<d->V_new_e; i++)
550: for(j=0; j<=i; j++) {
551: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
552: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
553: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
554: }
555: if (DVD_ISNOT(d->sH,DVD_MAT_UTRIANG)) /* Lower triangular part */
556: for(i=0; i<d->V_new_e; i++)
557: for(j=PetscMax(d->V_new_s,i+(DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)?1:0));
558: j<d->V_new_e; j++) {
559: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
560: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
561: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
562: }
564: return(0);
565: }
567: PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data, PetscScalar *ar,
568: PetscScalar *ai)
569: {
570: PetscScalar xr;
571: #if !defined(PETSC_USE_COMPLEX)
572: PetscScalar xi, k;
573: #endif
577: if(!ar) SETERRQ(1, "The real part has to be present!");
578: xr = *ar;
580: #if !defined(PETSC_USE_COMPLEX)
581: if(!ai) SETERRQ(1, "The imaginary part has to be present!");
582: xi = *ai;
584: if (xi != 0.0) {
585: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
586: data->Wa*data->Wa*xi*xi;
587: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
588: data->Wb*data->Wa*(xr*xr + xi*xi))/k;
589: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
590: } else
591: #endif
592: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
594: return(0);
595: }
598: PetscErrorCode dvd_harm_sort(EPS eps, PetscScalar ar, PetscScalar ai,
599: PetscScalar br, PetscScalar bi, PetscInt *r,
600: void *ctx)
601: {
602: dvdHarmonic *data = (dvdHarmonic*)ctx;
603: PetscErrorCode ierr;
607: /* Back-transform the harmonic values */
608: dvd_harm_backtrans(data, &ar, &ai);
609: dvd_harm_backtrans(data, &br, &bi);
611: /* Compare values using the user options for the eigenpairs selection */
612: eps->which = data->old_which;
613: eps->which_func = data->old_which_func;
614: eps->which_ctx = data->old_which_ctx;
615: EPSCompareEigenvalues(eps, ar, ai, br, bi, r);
617: /* Restore the eps values */
618: eps->which = EPS_WHICH_USER;
619: eps->which_func = dvd_harm_sort;
620: eps->which_ctx = data;
622: return(0);
623: }
625: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
626: {
627: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
628: PetscInt i;
632: for(i=0; i<d->size_H; i++)
633: dvd_harm_backtrans(data, &d->eigr[i], &d->eigi[i]);
635: return(0);
636: }