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-2013, Universitat Politecnica de Valencia, Spain
10: This file is part of SLEPc.
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,Vec Px);
29: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px);
30: PetscErrorCode dvd_precond_none(dvdDashboard *d,PetscInt i,Vec x,Vec Px);
31: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d);
33: typedef struct {
34: PC pc;
35: } 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: PetscBool t0,t1,t2;
50: /* Setup the step */
51: if (b->state >= DVD_STATE_CONF) {
52: /* If the preconditioner is valid */
53: if (pc) {
54: PetscMalloc(sizeof(dvdPCWrapper),&dvdpc);
55: PetscLogObjectMemory(d->eps,sizeof(dvdPCWrapper));
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: PCGetOperatorsSet(pc,NULL,&t0);
64: PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t1);
65: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t2);
66: if (t0 && !t1) {
67: PCGetOperators(pc,NULL,&P,NULL);
68: PetscObjectReference((PetscObject)P);
69: PCSetOperators(pc,P,P,SAME_PRECONDITIONER);
70: MatDestroy(&P);
71: } else if (t2) {
72: PCSetOperators(pc,d->A,d->A,SAME_PRECONDITIONER);
73: } else {
74: d->improvex_precond = dvd_precond_none;
75: }
77: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
79: /* Else, use no preconditioner */
80: } else d->improvex_precond = dvd_precond_none;
81: }
82: return(0);
83: }
87: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
88: {
89: PetscErrorCode ierr;
90: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
93: /* Free local data */
94: if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
95: PetscFree(d->improvex_precond_data);
96: return(0);
97: }
101: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
102: {
103: PetscErrorCode ierr;
104: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
107: PCApply(dvdpc->pc, x, Px);
108: return(0);
109: }
111: typedef struct {
112: Vec diagA, diagB;
113: } dvdJacobiPrecond;
117: /*
118: Create the Jacobi preconditioner for Generalized Eigenproblems
119: */
120: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d,dvdBlackboard *b)
121: {
122: PetscErrorCode ierr;
123: dvdJacobiPrecond *dvdjp;
124: PetscBool t;
127: /* Check if the problem matrices support GetDiagonal */
128: MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
129: if (t && d->B) {
130: MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
131: }
133: /* Setting configuration constrains */
134: b->own_vecs += t?((d->B == 0)?1:2) : 0;
136: /* Setup the step */
137: if (b->state >= DVD_STATE_CONF) {
138: if (t) {
139: PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
140: PetscLogObjectMemory(d->eps,sizeof(dvdJacobiPrecond));
141: dvdjp->diagA = *b->free_vecs;
142: b->free_vecs++;
143: MatGetDiagonal(d->A,dvdjp->diagA);
144: if (d->B) {
145: dvdjp->diagB = *b->free_vecs;
146: b->free_vecs++;
147: MatGetDiagonal(d->B, dvdjp->diagB);
148: } else dvdjp->diagB = 0;
149: d->improvex_precond_data = dvdjp;
150: d->improvex_precond = dvd_jacobi_precond_0;
152: DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);
154: /* Else, use no preconditioner */
155: } else d->improvex_precond = dvd_precond_none;
156: }
157: return(0);
158: }
162: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
163: {
164: PetscErrorCode ierr;
165: dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
168: /* Compute inv(D - eig)*x */
169: if (dvdjp->diagB == 0) {
170: /* Px <- diagA - l */
171: VecCopy(dvdjp->diagA, Px);
172: VecShift(Px, -d->eigr[i]);
173: } else {
174: /* Px <- diagA - l*diagB */
175: VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
176: }
178: /* Px(i) <- x/Px(i) */
179: VecPointwiseDivide(Px, x, Px);
180: return(0);
181: }
185: /*
186: Create a trivial preconditioner
187: */
188: PetscErrorCode dvd_precond_none(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
189: {
190: PetscErrorCode ierr;
193: 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,PetscInt max_size_D,PetscInt r_s,PetscInt r_e,PetscInt *size_D);
216: PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
217: PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
218: } DvdProfiler;
220: static PetscLogStage stages[6] = {0,0,0,0,0,0};
222: /*** Other things ****/
226: PetscErrorCode dvd_prof_init()
227: {
228: PetscErrorCode ierr;
231: if (!stages[0]) {
232: PetscLogStageRegister("Dvd_step_initV",&stages[DVD_STAGE_INITV]);
233: PetscLogStageRegister("Dvd_step_calcPairs",&stages[DVD_STAGE_CALCPAIRS]);
234: PetscLogStageRegister("Dvd_step_improveX",&stages[DVD_STAGE_IMPROVEX]);
235: PetscLogStageRegister("Dvd_step_updateV",&stages[DVD_STAGE_UPDATEV]);
236: PetscLogStageRegister("Dvd_step_orthV",&stages[DVD_STAGE_ORTHV]);
237: }
238: return(0);
239: }
243: PetscErrorCode dvd_initV_prof(dvdDashboard* d)
244: {
245: DvdProfiler *p = (DvdProfiler*)d->prof_data;
246: PetscErrorCode ierr;
249: PetscLogStagePush(stages[DVD_STAGE_INITV]);
250: p->old_initV(d);
251: PetscLogStagePop();
252: return(0);
253: }
257: PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d)
258: {
259: DvdProfiler *p = (DvdProfiler*)d->prof_data;
260: PetscErrorCode ierr;
263: PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
264: p->old_calcPairs(d);
265: PetscLogStagePop();
266: return(0);
267: }
271: PetscErrorCode dvd_improveX_prof(dvdDashboard* d,Vec *D,PetscInt max_size_D,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
272: {
273: DvdProfiler *p = (DvdProfiler*)d->prof_data;
274: PetscErrorCode ierr;
277: PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
278: p->old_improveX(d, D, max_size_D, r_s, r_e, size_D);
279: PetscLogStagePop();
280: return(0);
281: }
285: PetscErrorCode dvd_updateV_prof(dvdDashboard *d)
286: {
287: DvdProfiler *p = (DvdProfiler*)d->prof_data;
288: PetscErrorCode ierr;
291: PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
292: p->old_updateV(d);
293: PetscLogStagePop();
294: return(0);
295: }
299: PetscErrorCode dvd_orthV_prof(dvdDashboard *d)
300: {
301: DvdProfiler *p = (DvdProfiler*)d->prof_data;
302: PetscErrorCode ierr;
305: PetscLogStagePush(stages[DVD_STAGE_ORTHV]);
306: p->old_orthV(d);
307: PetscLogStagePop();
308: return(0);
309: }
313: PetscErrorCode dvd_profiler(dvdDashboard *d,dvdBlackboard *b)
314: {
315: PetscErrorCode ierr;
316: DvdProfiler *p;
319: /* Setup the step */
320: if (b->state >= DVD_STATE_CONF) {
321: PetscFree(d->prof_data);
322: PetscMalloc(sizeof(DvdProfiler),&p);
323: PetscLogObjectMemory(d->eps,sizeof(DvdProfiler));
324: d->prof_data = p;
325: p->old_initV = d->initV; d->initV = dvd_initV_prof;
326: p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
327: p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
328: p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;
330: DVD_FL_ADD(d->destroyList, dvd_profiler_d);
331: }
332: return(0);
333: }
337: PetscErrorCode dvd_profiler_d(dvdDashboard *d)
338: {
339: PetscErrorCode ierr;
340: DvdProfiler *p = (DvdProfiler*)d->prof_data;
343: /* Free local data */
344: PetscFree(p);
345: return(0);
346: }
349: /*
350: Configure the harmonics.
351: switch (mode) {
352: DVD_HARM_RR: harmonic RR
353: DVD_HARM_RRR: relative harmonic RR
354: DVD_HARM_REIGS: rightmost eigenvalues
355: DVD_HARM_LEIGS: largest eigenvalues
356: }
357: fixedTarged, if true use the target instead of the best eigenvalue
358: target, the fixed target to be used
359: */
360: typedef struct {
361: PetscScalar
362: Wa, Wb, /* span{W} = span{Wa*AV - Wb*BV} */
363: Pa, Pb; /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
364: PetscBool
365: withTarget;
366: HarmType_t
367: mode;
368: } dvdHarmonic;
370: PetscErrorCode dvd_harm_start(dvdDashboard *d);
371: PetscErrorCode dvd_harm_end(dvdDashboard *d);
372: PetscErrorCode dvd_harm_d(dvdDashboard *d);
373: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh,PetscScalar t);
374: PetscErrorCode dvd_harm_updateW(dvdDashboard *d);
375: PetscErrorCode dvd_harm_proj(dvdDashboard *d);
376: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d);
377: PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi);
381: PetscErrorCode dvd_harm_conf(dvdDashboard *d,dvdBlackboard *b,HarmType_t mode,PetscBool fixedTarget,PetscScalar t)
382: {
383: PetscErrorCode ierr;
384: dvdHarmonic *dvdh;
387: /* Set the problem to GNHEP:
388: d->G maybe is upper triangular due to biorthogonality of V and W */
389: d->sEP = d->sA = d->sB = 0;
391: /* Setup the step */
392: if (b->state >= DVD_STATE_CONF) {
393: PetscMalloc(sizeof(dvdHarmonic),&dvdh);
394: PetscLogObjectMemory(d->eps,sizeof(dvdHarmonic));
395: dvdh->withTarget = fixedTarget;
396: dvdh->mode = mode;
397: if (fixedTarget) dvd_harm_transf(dvdh, t);
398: d->calcpairs_W_data = dvdh;
399: d->calcpairs_W = dvd_harm_updateW;
400: d->calcpairs_proj_trans = dvd_harm_proj;
401: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
402: d->calcpairs_eig_backtrans = dvd_harm_eig_backtrans;
404: DVD_FL_ADD(d->destroyList, dvd_harm_d);
405: }
406: return(0);
407: }
411: PetscErrorCode dvd_harm_d(dvdDashboard *d)
412: {
413: PetscErrorCode ierr;
416: /* Free local data */
417: PetscFree(d->calcpairs_W_data);
418: return(0);
419: }
423: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh,PetscScalar t)
424: {
426: switch (dvdh->mode) {
427: case DVD_HARM_RR: /* harmonic RR */
428: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
429: case DVD_HARM_RRR: /* relative harmonic RR */
430: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
431: case DVD_HARM_REIGS: /* rightmost eigenvalues */
432: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
433: break;
434: case DVD_HARM_LEIGS: /* largest eigenvalues */
435: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
436: case DVD_HARM_NONE:
437: default:
438: SETERRQ(PETSC_COMM_SELF,1, "Harmonic type not supported");
439: }
441: /* Check the transformation does not change the sign of the imaginary part */
442: #if !defined(PETSC_USE_COMPLEX)
443: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
444: dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
445: #endif
446: return(0);
447: }
451: PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
452: {
453: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
454: PetscErrorCode ierr;
455: PetscInt i;
458: /* Update the target if it is necessary */
459: if (!data->withTarget) {
460: dvd_harm_transf(data,d->eigr[0]);
461: }
463: for (i=d->V_new_s;i<d->V_new_e;i++) {
464: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
465: VecAXPBYPCZ(d->W[i],data->Wa,-data->Wb,0.0,d->AV[i],(d->BV?d->BV:d->V)[i]);
466: }
467: return(0);
468: }
472: PetscErrorCode dvd_harm_proj(dvdDashboard *d)
473: {
474: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
475: PetscInt i,j;
478: if (d->sH != d->sG) SETERRQ(PETSC_COMM_SELF,1,"Projected matrices H and G must have the same structure");
480: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
481: if (DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)) /* Upper triangular part */
482: for (i=d->V_new_s+d->cX_in_H;i<d->V_new_e+d->cX_in_H;i++)
483: for (j=0;j<=i;j++) {
484: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
485: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
486: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
487: }
488: if (DVD_ISNOT(d->sH,DVD_MAT_UTRIANG)) /* Lower triangular part */
489: for (i=0;i<d->V_new_e+d->cX_in_H;i++)
490: for (j=PetscMax(d->V_new_s+d->cX_in_H,i+(DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)?1:0));
491: j<d->V_new_e+d->cX_in_H; j++) {
492: PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
493: d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
494: d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
495: }
496: return(0);
497: }
501: PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data,PetscScalar *ar,PetscScalar *ai)
502: {
503: PetscScalar xr;
504: #if !defined(PETSC_USE_COMPLEX)
505: PetscScalar xi, k;
506: #endif
510: xr = *ar;
511: #if !defined(PETSC_USE_COMPLEX)
513: xi = *ai;
515: if (xi != 0.0) {
516: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
517: data->Wa*data->Wa*xi*xi;
518: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
519: data->Wb*data->Wa*(xr*xr + xi*xi))/k;
520: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
521: } else
522: #endif
523: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
524: return(0);
525: }
529: PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi)
530: {
531: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
532: PetscErrorCode ierr;
535: dvd_harm_backtrans(data,&ar,&ai);
536: *br = ar;
537: *bi = ai;
538: return(0);
539: }
543: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
544: {
545: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
546: PetscInt i;
547: PetscErrorCode ierr;
550: for (i=0;i<d->size_H;i++) {
551: dvd_harm_backtrans(data, &d->eigr[i-d->cX_in_H], &d->eigi[i-d->cX_in_H]);
552: }
553: return(0);
554: }