Actual source code: dvdutils.c
slepc-3.7.3 2016-09-29
1: /*
2: SLEPc eigensolver: "davidson"
4: Some utils
6: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
7: SLEPc - Scalable Library for Eigenvalue Problem Computations
8: Copyright (c) 2002-2016, 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;
56: static PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
57: {
59: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
62: /* Free local data */
63: if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
64: PetscFree(d->improvex_precond_data);
65: return(0);
66: }
70: static PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
71: {
73: dvdPCWrapper *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;
76: PCApply(dvdpc->pc,x,Px);
77: return(0);
78: }
82: /*
83: Create a trivial preconditioner
84: */
85: static PetscErrorCode dvd_precond_none(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
86: {
90: VecCopy(x,Px);
91: return(0);
92: }
96: /*
97: Create a static preconditioner from a PC
98: */
99: PetscErrorCode dvd_static_precond_PC(dvdDashboard *d,dvdBlackboard *b,PC pc)
100: {
102: dvdPCWrapper *dvdpc;
103: Mat P;
104: PetscBool t0,t1,t2;
107: /* Setup the step */
108: if (b->state >= DVD_STATE_CONF) {
109: /* If the preconditioner is valid */
110: if (pc) {
111: PetscNewLog(d->eps,&dvdpc);
112: dvdpc->pc = pc;
113: PetscObjectReference((PetscObject)pc);
114: d->improvex_precond_data = dvdpc;
115: d->improvex_precond = dvd_static_precond_PC_0;
117: /* PC saves the matrix associated with the linear system, and it has to
118: be initialize to a valid matrix */
119: PCGetOperatorsSet(pc,NULL,&t0);
120: PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t1);
121: PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t2);
122: if (t0 && !t1) {
123: PCGetOperators(pc,NULL,&P);
124: PetscObjectReference((PetscObject)P);
125: PCSetOperators(pc,P,P);
126: PCSetReusePreconditioner(pc,PETSC_TRUE);
127: MatDestroy(&P);
128: } else if (t2) {
129: PCSetOperators(pc,d->A,d->A);
130: PCSetReusePreconditioner(pc,PETSC_TRUE);
131: } else {
132: d->improvex_precond = dvd_precond_none;
133: }
135: EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_precond_d);
137: /* Else, use no preconditioner */
138: } else d->improvex_precond = dvd_precond_none;
139: }
140: return(0);
141: }
145: static PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
146: {
147: PetscErrorCode ierr;
148: dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
151: /* Compute inv(D - eig)*x */
152: if (dvdjp->diagB == 0) {
153: /* Px <- diagA - l */
154: VecCopy(dvdjp->diagA,Px);
155: VecShift(Px,-d->eigr[i]);
156: } else {
157: /* Px <- diagA - l*diagB */
158: VecWAXPY(Px,-d->eigr[i],dvdjp->diagB,dvdjp->diagA);
159: }
161: /* Px(i) <- x/Px(i) */
162: VecPointwiseDivide(Px,x,Px);
163: return(0);
164: }
168: static PetscErrorCode dvd_jacobi_precond_d(dvdDashboard *d)
169: {
170: PetscErrorCode ierr;
171: dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;
174: if (dvdjp->diagA) { VecDestroy(&dvdjp->diagA); }
175: if (dvdjp->diagB) { VecDestroy(&dvdjp->diagB); }
176: PetscFree(d->improvex_precond_data);
177: return(0);
178: }
182: /*
183: Create the Jacobi preconditioner for Generalized Eigenproblems
184: */
185: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d,dvdBlackboard *b)
186: {
187: PetscErrorCode ierr;
188: dvdJacobiPrecond *dvdjp;
189: PetscBool t;
192: /* Check if the problem matrices support GetDiagonal */
193: MatHasOperation(d->A,MATOP_GET_DIAGONAL,&t);
194: if (t && d->B) {
195: MatHasOperation(d->B,MATOP_GET_DIAGONAL,&t);
196: }
198: /* Setup the step */
199: if (b->state >= DVD_STATE_CONF) {
200: PetscNewLog(d->eps,&dvdjp);
201: if (t) {
202: MatCreateVecs(d->A,&dvdjp->diagA,NULL);
203: MatGetDiagonal(d->A,dvdjp->diagA);
204: if (d->B) {
205: MatCreateVecs(d->B,&dvdjp->diagB,NULL);
206: MatGetDiagonal(d->B,dvdjp->diagB);
207: } else dvdjp->diagB = 0;
208: d->improvex_precond_data = dvdjp;
209: d->improvex_precond = dvd_jacobi_precond_0;
211: EPSDavidsonFLAdd(&d->destroyList,dvd_jacobi_precond_d);
213: /* Else, use no preconditioner */
214: } else {
215: dvdjp->diagA = dvdjp->diagB = 0;
216: d->improvex_precond = dvd_precond_none;
217: }
218: }
219: return(0);
220: }
224: static PetscErrorCode dvd_harm_d(dvdDashboard *d)
225: {
229: /* Free local data */
230: PetscFree(d->calcpairs_W_data);
231: return(0);
232: }
236: static PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh,PetscScalar t)
237: {
239: switch (dvdh->mode) {
240: case DVD_HARM_RR: /* harmonic RR */
241: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 0.0; dvdh->Pb = -1.0;
242: break;
243: case DVD_HARM_RRR: /* relative harmonic RR */
244: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = 0.0;
245: break;
246: case DVD_HARM_REIGS: /* rightmost eigenvalues */
247: dvdh->Wa = 1.0; dvdh->Wb = t; dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
248: break;
249: case DVD_HARM_LEIGS: /* largest eigenvalues */
250: dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0;
251: break;
252: case DVD_HARM_NONE:
253: default:
254: SETERRQ(PETSC_COMM_SELF,1, "Harmonic type not supported");
255: }
257: /* Check the transformation does not change the sign of the imaginary part */
258: #if !defined(PETSC_USE_COMPLEX)
259: if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0) {
260: dvdh->Pa *= -1.0;
261: dvdh->Pb *= -1.0;
262: }
263: #endif
264: return(0);
265: }
269: static PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
270: {
271: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
273: PetscInt l,k;
274: BV BX = d->BX?d->BX:d->eps->V;
277: /* Update the target if it is necessary */
278: if (!data->withTarget) {
279: dvd_harm_transf(data,d->eigr[0]);
280: }
282: /* W(i) <- Wa*AV(i) - Wb*BV(i) */
283: BVGetActiveColumns(d->eps->V,&l,&k);
284: if (k != l+d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
285: BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e);
286: BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
287: BVSetActiveColumns(BX,l+d->V_new_s,l+d->V_new_e);
288: BVCopy(d->AX,d->W);
289: BVScale(d->W,data->Wa);
290: BVMult(d->W,-data->Wb,1.0,BX,NULL);
291: BVSetActiveColumns(d->W,l,k);
292: BVSetActiveColumns(d->AX,l,k);
293: BVSetActiveColumns(BX,l,k);
294: return(0);
295: }
299: static PetscErrorCode dvd_harm_proj(dvdDashboard *d)
300: {
302: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
303: PetscInt i,j,l0,l,k,ld;
304: PetscScalar h,g,*H,*G;
307: BVGetActiveColumns(d->eps->V,&l0,&k);
308: l = l0 + d->V_new_s;
309: k = l0 + d->V_new_e;
310: MatGetSize(d->H,&ld,NULL);
311: MatDenseGetArray(d->H,&H);
312: MatDenseGetArray(d->G,&G);
313: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
314: /* Right part */
315: for (i=l;i<k;i++) {
316: for (j=l0;j<k;j++) {
317: h = H[ld*i+j];
318: g = G[ld*i+j];
319: H[ld*i+j] = data->Pa*h - data->Pb*g;
320: G[ld*i+j] = data->Wa*h - data->Wb*g;
321: }
322: }
323: /* Left part */
324: for (i=l0;i<l;i++) {
325: for (j=l;j<k;j++) {
326: h = H[ld*i+j];
327: g = G[ld*i+j];
328: H[ld*i+j] = data->Pa*h - data->Pb*g;
329: G[ld*i+j] = data->Wa*h - data->Wb*g;
330: }
331: }
332: MatDenseRestoreArray(d->H,&H);
333: MatDenseRestoreArray(d->G,&G);
334: return(0);
335: }
339: PetscErrorCode dvd_harm_updateproj(dvdDashboard *d)
340: {
342: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
343: PetscInt i,j,l,k,ld;
344: PetscScalar h,g,*H,*G;
347: BVGetActiveColumns(d->eps->V,&l,&k);
348: k = l + d->V_tra_s;
349: MatGetSize(d->H,&ld,NULL);
350: MatDenseGetArray(d->H,&H);
351: MatDenseGetArray(d->G,&G);
352: /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
353: /* Right part */
354: for (i=l;i<k;i++) {
355: for (j=0;j<l;j++) {
356: h = H[ld*i+j];
357: g = G[ld*i+j];
358: H[ld*i+j] = data->Pa*h - data->Pb*g;
359: G[ld*i+j] = data->Wa*h - data->Wb*g;
360: }
361: }
362: /* Lower triangular part */
363: for (i=0;i<l;i++) {
364: for (j=l;j<k;j++) {
365: h = H[ld*i+j];
366: g = G[ld*i+j];
367: H[ld*i+j] = data->Pa*h - data->Pb*g;
368: G[ld*i+j] = data->Wa*h - data->Wb*g;
369: }
370: }
371: MatDenseRestoreArray(d->H,&H);
372: MatDenseRestoreArray(d->G,&G);
373: return(0);
374: }
378: static PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data,PetscScalar *ar,PetscScalar *ai)
379: {
380: PetscScalar xr;
381: #if !defined(PETSC_USE_COMPLEX)
382: PetscScalar xi, k;
383: #endif
387: xr = *ar;
388: #if !defined(PETSC_USE_COMPLEX)
390: xi = *ai;
392: if (xi != 0.0) {
393: k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) + data->Wa*data->Wa*xi*xi;
394: *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr + data->Wb*data->Wa*(xr*xr + xi*xi))/k;
395: *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
396: } else
397: #endif
398: *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
399: return(0);
400: }
404: static PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi)
405: {
406: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
410: dvd_harm_backtrans(data,&ar,&ai);
411: *br = ar;
412: *bi = ai;
413: return(0);
414: }
418: static PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
419: {
420: dvdHarmonic *data = (dvdHarmonic*)d->calcpairs_W_data;
421: PetscInt i,l,k;
425: BVGetActiveColumns(d->eps->V,&l,&k);
426: for (i=0;i<k-l;i++) {
427: dvd_harm_backtrans(data,&d->eigr[i],&d->eigi[i]);
428: }
429: return(0);
430: }
434: PetscErrorCode dvd_harm_conf(dvdDashboard *d,dvdBlackboard *b,HarmType_t mode,PetscBool fixedTarget,PetscScalar t)
435: {
437: dvdHarmonic *dvdh;
440: /* Set the problem to GNHEP:
441: d->G maybe is upper triangular due to biorthogonality of V and W */
442: d->sEP = d->sA = d->sB = 0;
444: /* Setup the step */
445: if (b->state >= DVD_STATE_CONF) {
446: PetscNewLog(d->eps,&dvdh);
447: dvdh->withTarget = fixedTarget;
448: dvdh->mode = mode;
449: if (fixedTarget) dvd_harm_transf(dvdh, t);
450: d->calcpairs_W_data = dvdh;
451: d->calcpairs_W = dvd_harm_updateW;
452: d->calcpairs_proj_trans = dvd_harm_proj;
453: d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
454: d->calcpairs_eig_backtrans = dvd_harm_eig_backtrans;
456: EPSDavidsonFLAdd(&d->destroyList,dvd_harm_d);
457: }
458: return(0);
459: }