Actual source code: dvdutils.c
slepc-3.6.1 2015-09-03
1: /*
2: SLEPc eigensolver: "davidson"
4: Some utils
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2015, 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: typedef struct {
29: PC pc;
30: } dvdPCWrapper;
32: /*
33: Configure the harmonics.
34: switch (mode) {
35: DVD_HARM_RR: harmonic RR
36: DVD_HARM_RRR: relative harmonic RR
37: DVD_HARM_REIGS: rightmost eigenvalues
38: DVD_HARM_LEIGS: largest eigenvalues
39: }
40: fixedTarged, if true use the target instead of the best eigenvalue
41: target, the fixed target to be used
42: */
43: typedef struct {
44: PetscScalar Wa,Wb; /* span{W} = span{Wa*AV - Wb*BV} */
45: PetscScalar Pa,Pb; /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
46: PetscBool withTarget;
47: HarmType_t mode;
48: } dvdHarmonic;
50: typedef struct {
51: Vec diagA, diagB;
52: } dvdJacobiPrecond;
54: /*
55: Use of PETSc profiler functions
56: */
58: /* Define stages */
59: #define DVD_STAGE_INITV 0
60: #define DVD_STAGE_NEWITER 1
61: #define DVD_STAGE_CALCPAIRS 2
62: #define DVD_STAGE_IMPROVEX 3
63: #define DVD_STAGE_UPDATEV 4
64: #define DVD_STAGE_ORTHV 5
66: typedef struct {
67: PetscErrorCode (*old_initV)(struct _dvdDashboard*);
68: PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
69: PetscErrorCode (*old_improveX)(struct _dvdDashboard*,PetscInt r_s,PetscInt r_e,PetscInt *size_D);
70: PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
71: PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
72: } DvdProfiler;
74: static PetscLogStage stages[6] = {0,0,0,0,0,0};
78: static PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
79: {
81: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
84: /* Free local data */
85: if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
86: PetscFree(d->improvex_precond_data);
87: return(0);
88: }
92: static PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
93: {
95: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
98: PCApply(dvdpc->pc,x,Px);
99: return(0);
100: }
104: /*
105: Create a trivial preconditioner
106: */
107: static PetscErrorCode dvd_precond_none(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
108: {
112: VecCopy(x,Px);
113: return(0);
114: }
118: /*
119: Create a static preconditioner from a PC
120: */
121: PetscErrorCode dvd_static_precond_PC(dvdDashboard *d,dvdBlackboard *b,PC pc)
122: {
124: dvdPCWrapper *dvdpc;
125: Mat P;
126: PetscBool t0,t1,t2;
129: /* Setup the step */
130: if (b->state >= DVD_STATE_CONF) {
131: /* If the preconditioner is valid */
132: if (pc) {
133: PetscNewLog(d->eps,&dvdpc);
134: dvdpc->pc = pc;
135: PetscObjectReference((PetscObject)pc);
136: d->improvex_precond_data = dvdpc;
137: d->improvex_precond = dvd_static_precond_PC_0;
139: /* PC saves the matrix associated with the linear system, and it has to
140: be initialize to a valid matrix */
141: PCGetOperatorsSet(pc,NULL,&t0);
142: PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t1);
143: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t2);
144: if (t0 && !t1) {
145: PCGetOperators(pc,NULL,&P);
146: PetscObjectReference((PetscObject)P);
147: PCSetOperators(pc,P,P);
148: PCSetReusePreconditioner(pc,PETSC_TRUE);
149: MatDestroy(&P);
150: } else if (t2) {
151: PCSetOperators(pc,d->A,d->A);
152: PCSetReusePreconditioner(pc,PETSC_TRUE);
153: } else {
154: d->improvex_precond = dvd_precond_none;
155: }
157: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_precond_d);
159: /* Else, use no preconditioner */
160: } else d->improvex_precond = dvd_precond_none;
161: }
162: return(0);
163: }
167: static PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
168: {
169: PetscErrorCode ierr;
170: dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
173: /* Compute inv(D - eig)*x */
174: if (dvdjp->diagB == 0) {
175: /* Px <- diagA - l */
176: VecCopy(dvdjp->diagA,Px);
177: VecShift(Px,-d->eigr[i]);
178: } else {
179: /* Px <- diagA - l*diagB */
180: VecWAXPY(Px,-d->eigr[i],dvdjp->diagB,dvdjp->diagA);
181: }
183: /* Px(i) <- x/Px(i) */
184: VecPointwiseDivide(Px,x,Px);
185: return(0);
186: }
190: static PetscErrorCode dvd_jacobi_precond_d(dvdDashboard *d)
191: {
192: PetscErrorCode ierr;
193: dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
196: if (dvdjp->diagA) { VecDestroy(&dvdjp->diagA); }
197: if (dvdjp->diagB) { VecDestroy(&dvdjp->diagB); }
198: PetscFree(d->improvex_precond_data);
199: return(0);
200: }
204: /*
205: Create the Jacobi preconditioner for Generalized Eigenproblems
206: */
207: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d,dvdBlackboard *b)
208: {
209: PetscErrorCode ierr;
210: dvdJacobiPrecond *dvdjp;
211: PetscBool t;
214: /* Check if the problem matrices support GetDiagonal */
215: MatHasOperation(d->A,MATOP_GET_DIAGONAL,&t);
216: if (t && d->B) {
217: MatHasOperation(d->B,MATOP_GET_DIAGONAL,&t);
218: }
220: /* Setup the step */
221: if (b->state >= DVD_STATE_CONF) {
222: PetscNewLog(d->eps,&dvdjp);
223: if (t) {
224: MatCreateVecs(d->A,&dvdjp->diagA,NULL);
225: MatGetDiagonal(d->A,dvdjp->diagA);
226: if (d->B) {
227: MatCreateVecs(d->B,&dvdjp->diagB,NULL);
228: MatGetDiagonal(d->B,dvdjp->diagB);
229: } else dvdjp->diagB = 0;
230: d->improvex_precond_data = dvdjp;
231: d->improvex_precond = dvd_jacobi_precond_0;
233: EPSDavidsonFLAdd(&d->destroyList,dvd_jacobi_precond_d);
235: /* Else, use no preconditioner */
236: } else {
237: dvdjp->diagA = dvdjp->diagB = 0;
238: d->improvex_precond = dvd_precond_none;
239: }
240: }
241: return(0);
242: }
246: PetscErrorCode dvd_prof_init()
247: {
251: if (!stages[0]) {
252: PetscLogStageRegister("Dvd_step_initV",&stages[DVD_STAGE_INITV]);
253: PetscLogStageRegister("Dvd_step_calcPairs",&stages[DVD_STAGE_CALCPAIRS]);
254: PetscLogStageRegister("Dvd_step_improveX",&stages[DVD_STAGE_IMPROVEX]);
255: PetscLogStageRegister("Dvd_step_updateV",&stages[DVD_STAGE_UPDATEV]);
256: PetscLogStageRegister("Dvd_step_orthV",&stages[DVD_STAGE_ORTHV]);
257: }
258: return(0);
259: }
263: PetscErrorCode dvd_initV_prof(dvdDashboard* d)
264: {
265: DvdProfiler *p = (DvdProfiler*)d->prof_data;
269: PetscLogStagePush(stages[DVD_STAGE_INITV]);
270: p->old_initV(d);
271: PetscLogStagePop();
272: return(0);
273: }
277: static PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d)
278: {
279: DvdProfiler *p = (DvdProfiler*)d->prof_data;
283: PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
284: p->old_calcPairs(d);
285: PetscLogStagePop();
286: return(0);
287: }
291: static PetscErrorCode dvd_improveX_prof(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
292: {
293: DvdProfiler *p = (DvdProfiler*)d->prof_data;
297: PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
298: p->old_improveX(d,r_s,r_e,size_D);
299: PetscLogStagePop();
300: return(0);
301: }
305: static PetscErrorCode dvd_updateV_prof(dvdDashboard *d)
306: {
307: DvdProfiler *p = (DvdProfiler*)d->prof_data;
311: PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
312: p->old_updateV(d);
313: PetscLogStagePop();
314: return(0);
315: }
319: static PetscErrorCode dvd_profiler_d(dvdDashboard *d)
320: {
322: DvdProfiler *p = (DvdProfiler*)d->prof_data;
325: /* Free local data */
326: PetscFree(p);
327: return(0);
328: }
332: PetscErrorCode dvd_profiler(dvdDashboard *d,dvdBlackboard *b)
333: {
335: DvdProfiler *p;
338: /* Setup the step */
339: if (b->state >= DVD_STATE_CONF) {
340: PetscFree(d->prof_data);
341: PetscNewLog(d->eps,&p);
342: d->prof_data = p;
343: p->old_initV = d->initV; d->initV = dvd_initV_prof;
344: p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
345: p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
346: p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;
348: EPSDavidsonFLAdd(&d->destroyList,dvd_profiler_d);
349: }
350: return(0);
351: }
355: static PetscErrorCode dvd_harm_d(dvdDashboard *d)
356: {
360: /* Free local data */
361: PetscFree(d->calcpairs_W_data);
362: return(0);
363: }
367: static PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh,PetscScalar t)
368: {
370: switch (dvdh->mode) {
371: case DVD_HARM_RR: /* harmonic RR */
372: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0;
373: break;
374: case DVD_HARM_RRR: /* relative harmonic RR */
375: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0;
376: break;
377: case DVD_HARM_REIGS: /* rightmost eigenvalues */
378: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
379: break;
380: case DVD_HARM_LEIGS: /* largest eigenvalues */
381: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0;
382: break;
383: case DVD_HARM_NONE:
384: default:
385: SETERRQ(PETSC_COMM_SELF,1, "Harmonic type not supported");
386: }
388: /* Check the transformation does not change the sign of the imaginary part */
389: #if !defined(PETSC_USE_COMPLEX)
390: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0) {
391: dvdh->Pa *= -1.0;
392: dvdh->Pb *= -1.0;
393: }
394: #endif
395: return(0);
396: }
400: static PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
401: {
402: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
404: PetscInt l,k,i;
405: BV BX = d->BX?d->BX:d->eps->V;
408: /* Update the target if it is necessary */
409: if (!data->withTarget) {
410: dvd_harm_transf(data,d->eigr[0]);
411: }
413: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
414: BVGetActiveColumns(d->eps->V,&l,&k);
415: if (k != l+d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
416: BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e);
417: BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
418: BVSetActiveColumns(BX,l+d->V_new_s,l+d->V_new_e);
419: BVCopy(d->AX,d->W);
420: /* Work around bug in BVScale
421: BVScale(d->W,data->Wa); */
422: for (i=l+d->V_new_s;i<l+d->V_new_e; ++i) {
423: BVScaleColumn(d->W,i,data->Wa);
424: }
425: BVAXPY(d->W,-data->Wb,BX);
426: BVSetActiveColumns(d->W,l,k);
427: BVSetActiveColumns(d->AX,l,k);
428: BVSetActiveColumns(BX,l,k);
429: return(0);
430: }
434: static PetscErrorCode dvd_harm_proj(dvdDashboard *d)
435: {
437: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
438: PetscInt i,j,l0,l,k,ld;
439: PetscScalar h,g,*H,*G;
442: BVGetActiveColumns(d->eps->V,&l0,&k);
443: l = l0 + d->V_new_s;
444: k = l0 + d->V_new_e;
445: MatGetSize(d->H,&ld,NULL);
446: MatDenseGetArray(d->H,&H);
447: MatDenseGetArray(d->G,&G);
448: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
449: /* Right part */
450: for (i=l;i<k;i++) {
451: for (j=l0;j<k;j++) {
452: h = H[ld*i+j];
453: g = G[ld*i+j];
454: H[ld*i+j] = data->Pa*h - data->Pb*g;
455: G[ld*i+j] = data->Wa*h - data->Wb*g;
456: }
457: }
458: /* Left part */
459: for (i=l0;i<l;i++) {
460: for (j=l;j<k;j++) {
461: h = H[ld*i+j];
462: g = G[ld*i+j];
463: H[ld*i+j] = data->Pa*h - data->Pb*g;
464: G[ld*i+j] = data->Wa*h - data->Wb*g;
465: }
466: }
467: MatDenseRestoreArray(d->H,&H);
468: MatDenseRestoreArray(d->G,&G);
469: return(0);
470: }
474: PetscErrorCode dvd_harm_updateproj(dvdDashboard *d)
475: {
477: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
478: PetscInt i,j,l,k,ld;
479: PetscScalar h,g,*H,*G;
482: BVGetActiveColumns(d->eps->V,&l,&k);
483: k = l + d->V_tra_s;
484: MatGetSize(d->H,&ld,NULL);
485: MatDenseGetArray(d->H,&H);
486: MatDenseGetArray(d->G,&G);
487: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
488: /* Right part */
489: for (i=l;i<k;i++) {
490: for (j=0;j<l;j++) {
491: h = H[ld*i+j];
492: g = G[ld*i+j];
493: H[ld*i+j] = data->Pa*h - data->Pb*g;
494: G[ld*i+j] = data->Wa*h - data->Wb*g;
495: }
496: }
497: /* Lower triangular part */
498: for (i=0;i<l;i++) {
499: for (j=l;j<k;j++) {
500: h = H[ld*i+j];
501: g = G[ld*i+j];
502: H[ld*i+j] = data->Pa*h - data->Pb*g;
503: G[ld*i+j] = data->Wa*h - data->Wb*g;
504: }
505: }
506: MatDenseRestoreArray(d->H,&H);
507: MatDenseRestoreArray(d->G,&G);
508: return(0);
509: }
513: static PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data,PetscScalar *ar,PetscScalar *ai)
514: {
515: PetscScalar xr;
516: #if !defined(PETSC_USE_COMPLEX)
517: PetscScalar xi, k;
518: #endif
522: xr = *ar;
523: #if !defined(PETSC_USE_COMPLEX)
525: xi = *ai;
527: if (xi != 0.0) {
528: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) + data->Wa*data->Wa*xi*xi;
529: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr + data->Wb*data->Wa*(xr*xr + xi*xi))/k;
530: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
531: } else
532: #endif
533: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
534: return(0);
535: }
539: static PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi)
540: {
541: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
545: dvd_harm_backtrans(data,&ar,&ai);
546: *br = ar;
547: *bi = ai;
548: return(0);
549: }
553: static PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
554: {
555: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
556: PetscInt i,l,k;
560: BVGetActiveColumns(d->eps->V,&l,&k);
561: for (i=0;i<k-l;i++) {
562: dvd_harm_backtrans(data,&d->eigr[i],&d->eigi[i]);
563: }
564: return(0);
565: }
569: PetscErrorCode dvd_harm_conf(dvdDashboard *d,dvdBlackboard *b,HarmType_t mode,PetscBool fixedTarget,PetscScalar t)
570: {
572: dvdHarmonic *dvdh;
575: /* Set the problem to GNHEP:
576: d->G maybe is upper triangular due to biorthogonality of V and W */
577: d->sEP = d->sA = d->sB = 0;
579: /* Setup the step */
580: if (b->state >= DVD_STATE_CONF) {
581: PetscNewLog(d->eps,&dvdh);
582: dvdh->withTarget = fixedTarget;
583: dvdh->mode = mode;
584: if (fixedTarget) dvd_harm_transf(dvdh, t);
585: d->calcpairs_W_data = dvdh;
586: d->calcpairs_W = dvd_harm_updateW;
587: d->calcpairs_proj_trans = dvd_harm_proj;
588: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
589: d->calcpairs_eig_backtrans = dvd_harm_eig_backtrans;
591: EPSDavidsonFLAdd(&d->destroyList,dvd_harm_d);
592: }
593: return(0);
594: }
598: /*
599: H = [H Y(old)'*X(new);
600: Y(new)'*X(old) Y(new)'*X(new) ],
601: being old=0:l-1, new=l:k-1
602: */
603: PetscErrorCode BVMultS(BV X,BV Y,PetscScalar *H,PetscInt ldh)
604: {
606: PetscInt j,lx,ly,kx,ky;
607: PetscScalar *array;
608: Mat M;
611: BVGetActiveColumns(X,&lx,&kx);
612: BVGetActiveColumns(Y,&ly,&ky);
613: MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
614: BVMatProject(X,NULL,Y,M);
615: MatDenseGetArray(M,&array);
616: /* upper part */
617: for (j=lx;j<kx;j++) {
618: PetscMemcpy(&H[ldh*j],&array[j*ky],ly*sizeof(PetscScalar));
619: }
620: /* lower part */
621: for (j=0;j<kx;j++) {
622: PetscMemcpy(&H[ldh*j+ly],&array[j*ky+ly],(ky-ly)*sizeof(PetscScalar));
623: }
624: MatDenseRestoreArray(M,&array);
625: MatDestroy(&M);
626: return(0);
627: }