Actual source code: dvd_calcpairs.c
1: /*
2: SLEPc eigensolver: "davidson"
4: Step: calc the best eigenpairs in the subspace V.
6: For that, performs these steps:
7: 1) Update W <- A * V
8: 2) Update H <- V' * W
9: 3) Obtain eigenpairs of H
10: 4) Select some eigenpairs
11: 5) Compute the Ritz pairs of the selected ones
13: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
14: SLEPc - Scalable Library for Eigenvalue Problem Computations
15: Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
17: This file is part of SLEPc.
18:
19: SLEPc is free software: you can redistribute it and/or modify it under the
20: terms of version 3 of the GNU Lesser General Public License as published by
21: the Free Software Foundation.
23: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
24: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
25: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
26: more details.
28: You should have received a copy of the GNU Lesser General Public License
29: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
30: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31: */
33: #include davidson.h
34: #include <slepcblaslapack.h>
36: PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d);
37: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d);
38: PetscErrorCode dvd_calcpairs_projeig_eig(dvdDashboard *d);
39: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d);
40: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d);
41: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n);
42: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n);
43: PetscErrorCode dvd_calcpairs_X(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
44: Vec *X);
45: PetscErrorCode dvd_calcpairs_Y(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
46: Vec *Y);
47: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
48: Vec *R);
49: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
50: PetscInt r_e, Vec *R);
51: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
52: DvdMult_copy_func **sr);
53: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d);
54: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
55: DvdMult_copy_func **sr);
56: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d);
57: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d);
58: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
59: DvdMult_copy_func **sr);
60: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d);
61: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
62: DvdMult_copy_func **sr);
63: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,Vec *real_BV,PetscInt *size_cX,Vec **BV,PetscInt *size_BV,PetscInt *max_size_BV,PetscBool BV_shift,PetscInt *cX_in_proj,PetscScalar *MTX);
65: /**** Control routines ********************************************************/
68: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d, dvdBlackboard *b,
69: orthoV_type_t orth, IP ipI,
70: PetscInt cX_proj, PetscBool harm)
71: {
72: PetscErrorCode ierr;
73: PetscInt i;
74: PetscBool std_probl,her_probl;
78: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
79: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
81: /* Setting configuration constrains */
82: #if !defined(PETSC_USE_COMPLEX)
83: /* if the last converged eigenvalue is complex its conjugate pair is also
84: converged */
85: b->max_nev = PetscMax(b->max_nev, d->nev+(her_probl && !d->B?0:1));
86: #else
87: b->max_nev = PetscMax(b->max_nev, d->nev);
88: #endif
89: b->max_size_proj = PetscMax(b->max_size_proj, b->max_size_V+cX_proj);
90: d->size_real_V = b->max_size_V+b->max_nev;
91: d->W_shift = d->B?PETSC_TRUE:PETSC_FALSE;
92: d->size_real_W = harm?(b->max_size_V+(d->W_shift?b->max_nev:b->max_size_cP)):0;
93: d->size_real_AV = b->max_size_V+b->max_size_cP;
94: d->size_BDS = 0;
95: if (d->B && her_probl && (orth == DVD_ORTHOV_I || orth == DVD_ORTHOV_BOneMV)) {
96: d->size_real_BV = b->size_V; d->BV_shift = PETSC_TRUE;
97: if (orth == DVD_ORTHOV_BOneMV) d->size_BDS = d->eps->nds;
98: } else if (d->B) {
99: d->size_real_BV = b->max_size_V + b->max_size_P; d->BV_shift = PETSC_FALSE;
100: } else {
101: d->size_real_BV = 0; d->BV_shift = PETSC_FALSE;
102: }
103: b->own_vecs+= d->size_real_V + d->size_real_W + d->size_real_AV +
104: d->size_real_BV + d->size_BDS;
105: b->own_scalars+= b->max_size_proj*b->max_size_proj*2*(std_probl?1:2) +
106: /* H, G?, S, T? */
107: b->max_size_proj*b->max_size_proj*(std_probl?1:2) +
108: /* pX, pY? */
109: b->max_nev*b->max_nev*(her_probl?0:(!d->B?1:2));
110: /* cS?, cT? */
111: b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X);
112: /* updateV0 */
113: b->max_size_auxS = PetscMax(PetscMax(
114: b->max_size_auxS,
115: b->max_size_X*b->max_size_proj*2*(!d->B?1:2) + /* updateAV1,BV1 */
116: b->max_size_X*b->max_nev*(her_probl?0:(!d->B?1:2)) + /* updateV0,W0 */
117: /* SlepcReduction: in */
118: PetscMax(
119: b->max_size_X*b->max_size_proj*2*(!d->B?1:2) + /* updateAV1,BV1 */
120: b->max_size_X*b->max_nev*(her_probl?0:(!d->B?1:2)), /* updateV0,W0 */
121: /* SlepcReduction: out */
122: PetscMax(
123: b->max_size_proj*b->max_size_proj, /* updateAV0,BV0 */
124: b->max_size_proj+b->max_nev))), /* dvd_orth */
125: std_probl?0:(b->max_size_proj*11+16) /* projeig */);
126: #if defined(PETSC_USE_COMPLEX)
127: b->max_size_auxS = PetscMax(b->max_size_auxS, b->max_size_V);
128: /* dvd_calcpairs_projeig_eig */
129: #endif
131: /* Setup the step */
132: if (b->state >= DVD_STATE_CONF) {
133: d->max_cX_in_proj = cX_proj;
134: d->max_size_P = b->max_size_P;
135: d->real_V = b->free_vecs; b->free_vecs+= d->size_real_V;
136: if (harm) {
137: d->real_W = b->free_vecs; b->free_vecs+= d->size_real_W;
138: } else {
139: d->real_W = PETSC_NULL;
140: }
141: d->real_AV = d->AV = b->free_vecs; b->free_vecs+= d->size_real_AV;
142: d->max_size_proj = b->max_size_proj;
143: d->real_H = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
144: d->ldH = b->max_size_proj;
145: d->pX = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
146: d->S = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
147: if (!her_probl) {
148: d->cS = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
149: d->max_size_cS = d->ldcS = b->max_nev;
150: } else {
151: d->cS = PETSC_NULL;
152: d->max_size_cS = d->ldcS = 0;
153: }
154: d->ipV = ipI;
155: d->ipW = ipI;
156: if (orth == DVD_ORTHOV_BOneMV) {
157: d->BDS = b->free_vecs; b->free_vecs+= d->eps->nds;
158: for (i=0; i<d->eps->nds; i++) {
159: MatMult(d->B, d->eps->DS[i], d->BDS[i]);
160: }
161: } else
162: d->BDS = PETSC_NULL;
163: if (d->B) {
164: d->real_BV = b->free_vecs; b->free_vecs+= d->size_real_BV;
165: } else {
166: d->size_real_BV = 0;
167: d->real_BV = PETSC_NULL;
168: d->BV_shift = PETSC_FALSE;
169: }
170: if (!std_probl) {
171: d->real_G = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
172: d->T = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
173: d->pY = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
174: } else {
175: d->real_G = PETSC_NULL;
176: d->T = PETSC_NULL;
177: d->pY = PETSC_NULL;
178: }
179: if (d->B && !her_probl) {
180: d->cT = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
181: d->ldcT = b->max_nev;
182: } else {
183: d->cT = PETSC_NULL;
184: d->ldcT = 0;
185: }
187: d->calcPairs = dvd_calcpairs_proj;
188: d->calcpairs_residual = dvd_calcpairs_res_0;
189: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
190: d->calcpairs_selectPairs = PETSC_NULL;
191: d->ipI = ipI;
192: DVD_FL_ADD(d->startList, dvd_calcpairs_qz_start);
193: }
195: return(0);
196: }
201: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
202: {
203: PetscBool std_probl, her_probl;
204: PetscInt i;
208: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
209: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
211: d->size_V = 0;
212: d->V = d->real_V;
213: d->cX = d->real_V;
214: d->size_cX = 0;
215: d->max_size_V = d->size_real_V;
216: d->W = d->real_W;
217: d->max_size_W = d->size_real_W;
218: d->size_W = 0;
219: d->size_AV = 0;
220: d->AV = d->real_AV;
221: d->max_size_AV = d->size_real_AV;
222: d->size_H = 0;
223: d->H = d->real_H;
224: if (d->cS) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cS[i] = 0.0;
225: d->size_BV = 0;
226: d->BV = d->real_BV;
227: d->max_size_BV = d->size_real_BV;
228: d->size_G = 0;
229: d->G = d->real_G;
230: if (d->cT) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cT[i] = 0.0;
231: d->cY = d->B && !her_probl ? d->W : PETSC_NULL;
232: d->BcX = d->orthoV_type == DVD_ORTHOV_I && d->B && her_probl ? d->BcX : PETSC_NULL;
233: d->size_cY = 0;
234: d->size_BcX = 0;
235: d->cX_in_V = d->cX_in_H = d->cX_in_G = d->cX_in_W = d->cX_in_AV = d->cX_in_BV = 0;
237: return(0);
238: }
243: PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
244: {
245: PetscErrorCode ierr;
246: DvdReduction r;
247: static const PetscInt MAX_OPS = 7;
248: DvdReductionChunk
249: ops[MAX_OPS];
250: DvdMult_copy_func
251: sr[MAX_OPS], *sr0 = sr;
252: PetscInt size_in, i;
253: PetscScalar *in = d->auxS, *out;
254: PetscBool stdp;
258: stdp = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
259: size_in =
260: (d->size_cX+d->V_tra_s)*d->V_tra_s*(d->cT?2:(d->cS?1:0)) + /* updateV0,W0 */
261: d->size_H*d->V_tra_e*(!stdp?2:1) + /* updateAV1,BV1 */
262: (d->size_H*(d->V_new_e-d->V_new_s)*2+
263: (d->V_new_e-d->V_new_s)*(d->V_new_e-d->V_new_s))*(!stdp?2:1); /* updateAV0,BV0 */
264:
265: out = in+size_in;
267: /* Check consistency */
268: if (2*size_in > d->size_auxS) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
270: /* Prepare reductions */
271: SlepcAllReduceSumBegin(ops, MAX_OPS, in, out, size_in, &r,
272: ((PetscObject)d->V[0])->comm);
273: /* Allocate size_in */
274: d->auxS+= size_in;
275: d->size_auxS-= size_in;
277: /* Update AV, BV, W and the projected matrices */
278: /* 1. S <- S*MT */
279: dvd_calcpairs_updateV0(d, &r, &sr0);
280: dvd_calcpairs_updateW0(d, &r, &sr0);
281: dvd_calcpairs_updateAV0(d);
282: dvd_calcpairs_updateBV0(d);
283: /* 2. V <- orth(V, V_new) */
284: dvd_calcpairs_updateV1(d);
285: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
286: /* Check consistency */
287: if (d->size_AV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
288: for (i=d->V_new_s; i<d->V_new_e; i++) {
289: MatMult(d->A, d->V[i], d->AV[i]);
290: }
291: d->size_AV = d->V_new_e;
292: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
293: if (d->B && d->orthoV_type != DVD_ORTHOV_BOneMV) {
294: /* Check consistency */
295: if (d->size_BV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
296: for (i=d->V_new_s; i<d->V_new_e; i++) {
297: MatMult(d->B, d->V[i], d->BV[i]);
298: }
299: d->size_BV = d->V_new_e;
300: }
301: /* 5 <- W <- [W f(AV,BV)] */
302: dvd_calcpairs_updateW1(d);
303: dvd_calcpairs_updateAV1(d, &r, &sr0);
304: dvd_calcpairs_updateBV1(d, &r, &sr0);
306: /* Deallocate size_in */
307: d->auxS-= size_in;
308: d->size_auxS+= size_in;
310: /* Do reductions */
311: SlepcAllReduceSumEnd(&r);
313: /* Perform the transformation on the projected problem */
314: if (d->calcpairs_proj_trans) {
315: d->calcpairs_proj_trans(d);
316: }
318: d->V_tra_s = d->V_tra_e = 0;
319: d->V_new_s = d->V_new_e;
321: /* Solve the projected problem */
322: if (DVD_IS(d->sEP, DVD_EP_STD)) {
323: if (DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
324: dvd_calcpairs_projeig_eig(d);
325: } else {
326: dvd_calcpairs_projeig_qz_std(d);
327: }
328: } else {
329: dvd_calcpairs_projeig_qz_gen(d);
330: }
332: /* Check consistency */
333: if (d->size_V != d->V_new_e || d->size_V+d->cX_in_H != d->size_H || d->cX_in_V != d->cX_in_H ||
334: d->size_V != d->size_AV || d->cX_in_H != d->cX_in_AV ||
335: (DVD_ISNOT(d->sEP, DVD_EP_STD) && (
336: d->size_V+d->cX_in_G != d->size_G || d->cX_in_H != d->cX_in_G ||
337: d->size_H != d->size_G || (d->BV && (
338: d->size_V != d->size_BV || d->cX_in_H != d->cX_in_BV)))) ||
339: (d->W && d->size_W != d->size_V)) {
340: SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!");
341: }
343: return(0);
344: }
346: /**** Basic routines **********************************************************/
350: /* auxV: V_tra_s, DvdMult_copy_func: 1 */
351: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
352: DvdMult_copy_func **sr)
353: {
354: PetscErrorCode ierr;
355: PetscInt rm;
359: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
361: /* cX <- [cX V*MT(0:V_tra_s-1)], V <- V*MT(V_tra_s:V_tra_e) */
362: dvd_calcpairs_updateBV0_gen(d,d->real_V,&d->size_cX,&d->V,&d->size_V,&d->max_size_V,PETSC_TRUE,&d->cX_in_V,d->MTX);
364: /* Udpate cS for standard problems */
365: if (d->cS && !d->cT && !d->cY && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
366: /* Check consistency */
367: if (d->size_cS+d->V_tra_s != d->size_cX) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
369: /* auxV <- AV * MTX(0:V_tra_e-1) */
370: SlepcUpdateVectorsZ(d->auxV, 0.0, 1.0, d->AV-d->cX_in_AV, d->size_AV+d->cX_in_AV, d->MTX, d->ldMTX, d->size_MT, d->V_tra_s-d->max_cX_in_proj);
372: /* cS(:, size_cS:) <- cX' * auxV */
373: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
374: VecsMultS(&d->cS[d->ldcS*d->size_cS], 0, d->ldcS, d->cX, 0, d->size_cX-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
375: d->size_cS+= d->V_tra_s-rm;
376: }
377:
378: return(0);
379: }
384: /* auxS: size_cX+V_new_e+1 */
385: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d)
386: {
387: PetscErrorCode ierr;
388: Vec *cX = d->BcX? d->BcX : ( (d->cY && !d->W)? d->cY : d->cX );
392: if (d->V_new_s == d->V_new_e) { return(0); }
394: /* Check consistency */
395: if (d->size_V != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
397: /* V <- gs([cX V(0:V_new_s-1)], V(V_new_s:V_new_e-1)) */
398: if (d->orthoV_type == DVD_ORTHOV_BOneMV) {
399: dvd_BorthV(d->ipV, d->eps->DS, d->BDS, d->eps->nds, d->cX, d->real_BV,
400: d->size_cX, d->V, d->BV, d->V_new_s, d->V_new_e,
401: d->auxS, d->eps->rand);
402: d->size_BV = d->V_new_e;
403: } else {
404: dvd_orthV(d->ipV, d->eps->DS, d->eps->nds, cX, d->size_cX, d->V,
405: d->V_new_s, d->V_new_e, d->auxS, d->eps->rand);
406:
407: }
408: d->size_V = d->V_new_e;
410: return(0);
411: }
415: /* auxV: V_tra_s, DvdMult_copy_func: 2 */
416: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
417: DvdMult_copy_func **sr)
418: {
419: PetscErrorCode ierr;
420: PetscInt rm;
424: if (!d->W || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }
426: /* cY <- [cY W*MTY(0:V_tra_s-1)], W <- W*MTY(V_tra_s:V_tra_e) */
427: dvd_calcpairs_updateBV0_gen(d,d->real_W,&d->size_cY,&d->W,&d->size_W,&d->max_size_W,d->W_shift,&d->cX_in_W,d->MTY);
429: /* Udpate cS and cT */
430: if (d->cY && d->cT && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
431: /* Check consistency */
432: if (d->size_cS+d->V_tra_s != d->size_cY) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
434: /* auxV <- AV * MTX(0:V_tra_e-1) */
435: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
436: SlepcUpdateVectorsZ(d->auxV, 0.0, 1.0, d->AV-d->cX_in_H, d->size_AV-d->cX_in_H, d->MTX, d->ldMTX, d->size_MT, d->V_tra_s-d->max_cX_in_proj);
438: /* cS(:, size_cS:) <- cY' * auxV */
439: VecsMultS(&d->cS[d->ldcS*d->size_cS], 0, d->ldcS, d->cY, 0, d->size_cY-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
441: /* auxV <- BV * MTX(0:V_tra_e-1) */
442: SlepcUpdateVectorsZ(d->auxV, 0.0, 1.0, d->BV-d->cX_in_H, d->size_BV-d->cX_in_H, d->MTX, d->ldMTX, d->size_MT, d->V_tra_s-d->max_cX_in_proj);
444: /* cT(:, size_cS:) <- cY' * auxV */
445: VecsMultS(&d->cT[d->ldcS*d->size_cS], 0, d->ldcS, d->cY, 0, d->size_cY-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
446: d->size_cS+= d->V_tra_s-rm;
447: d->size_cT+= d->V_tra_s-rm;
448: }
450: return(0);
451: }
456: /* auxS: size_cX+V_new_e+1 */
457: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d)
458: {
459: PetscErrorCode ierr;
463: if (!d->W || (d->V_new_s == d->V_new_e)) { return(0); }
465: /* Check consistency */
466: if (d->size_W != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
468: /* Update W */
469: d->calcpairs_W(d);
471: /* W <- gs([cY W(0:V_new_s-1)], W(V_new_s:V_new_e-1)) */
472: dvd_orthV(d->ipW, PETSC_NULL, 0, d->cY, d->size_cY, d->W, d->V_new_s,
473: d->V_new_e, d->auxS, d->eps->rand);
474:
475: d->size_W = d->V_new_e;
477: return(0);
478: }
482: /* auxS: size_H*(V_tra_e-V_tra_s) */
483: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d)
484: {
485: PetscErrorCode ierr;
486: PetscScalar *MTY = d->W?d->MTY:d->MTX;
487: PetscInt cMT, tra_s, rm, cp;
491: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
493: /* AV(V_tra_s-cp-1:) = cAV*MTX(V_tra_s:) */
494: dvd_calcpairs_updateBV0_gen(d,d->real_AV,PETSC_NULL,&d->AV,&d->size_AV,&d->max_size_AV,PETSC_FALSE,&d->cX_in_AV,d->MTX);
495: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
496: cMT = d->V_tra_e - tra_s;
497: rm = d->V_tra_s - tra_s;
498: cp = PetscMin(d->max_cX_in_proj - rm, d->cX_in_H);
500: /* Update H <- MTY(tra_s)' * (H * MTX(tra_s:)) */
501: SlepcDenseMatProdTriang(d->auxS, 0, d->ldH, d->H, d->sH, d->ldH, d->size_H, d->size_H, PETSC_FALSE, &d->MTX[d->ldMTX*tra_s], 0, d->ldMTX, d->size_MT, cMT, PETSC_FALSE);
502: SlepcDenseMatProdTriang(d->H, d->sH, d->ldH, &MTY[d->ldMTX*tra_s], 0, d->ldMTX, d->size_MT, cMT, PETSC_TRUE, d->auxS, 0, d->ldH, d->size_H, cMT, PETSC_FALSE);
503: d->size_H = cMT;
504: d->cX_in_H = cp+rm;
506: return(0);
507: }
512: /* DvdMult_copy_func: 2 */
513: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
514: DvdMult_copy_func **sr)
515: {
516: PetscErrorCode ierr;
517: Vec *W = d->W?d->W:d->V;
521: if (d->V_new_s == d->V_new_e) { return(0); }
523: /* Check consistency */
524: if (d->size_H != d->V_new_s+d->cX_in_H || d->size_V != d->V_new_e) {
525: SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!");
526: }
528: /* H = [H W(old)'*AV(new);
529: W(new)'*AV(old) W(new)'*AV(new) ],
530: being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
531: VecsMultS(d->H,d->sH,d->ldH,W-d->cX_in_H,d->V_new_s+d->cX_in_H, d->V_new_e+d->cX_in_H, d->AV-d->cX_in_H,d->V_new_s+d->cX_in_H,d->V_new_e+d->cX_in_H, r, (*sr)++);
532: d->size_H = d->V_new_e+d->cX_in_H;
534: return(0);
535: }
539: /* auxS: max(BcX*(size_cX+V_new_e+1), size_G*(V_tra_e-V_tra_s)) */
540: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d)
541: {
542: PetscErrorCode ierr;
543: PetscScalar *MTY = d->W?d->MTY:d->MTX;
544: PetscInt cMT, rm, cp, tra_s, i;
545: PetscBool lindep;
546: PetscReal norm;
550: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
552: /* BV <- BV*MTX */
553: dvd_calcpairs_updateBV0_gen(d,d->real_BV,PETSC_NULL,&d->BV,&d->size_BV,&d->max_size_BV,d->BV_shift,&d->cX_in_BV,d->MTX);
555: /* If BcX, BcX <- orth(BcX) */
556: if (d->BcX) {
557: for (i=0; i<d->V_tra_s; i++) {
558: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_BcX+i, PETSC_NULL,
559: d->BcX, d->BcX[d->size_BcX+i], PETSC_NULL,
560: &norm, &lindep);
561: if(lindep) SETERRQ(((PetscObject)d->ipI)->comm,1, "Error during orth(BcX, B*cX(new))");
562: VecScale(d->BcX[d->size_BcX+i], 1.0/norm);
563: }
564: d->size_BcX+= d->V_tra_s;
565: }
567: /* Update G <- MTY' * (G * MTX) */
568: if (d->G) {
569: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
570: cMT = d->V_tra_e - tra_s;
571: rm = d->V_tra_s - tra_s;
572: cp = PetscMin(d->max_cX_in_proj - rm, d->cX_in_G);
573: SlepcDenseMatProdTriang(d->auxS, 0, d->ldH, d->G, d->sG, d->ldH, d->size_G, d->size_G, PETSC_FALSE, &d->MTX[d->ldMTX*tra_s], 0, d->ldMTX, d->size_MT, cMT, PETSC_FALSE);
574: SlepcDenseMatProdTriang(d->G, d->sG, d->ldH, &MTY[d->ldMTX*tra_s], 0, d->ldMTX, d->size_MT, cMT, PETSC_TRUE, d->auxS, 0, d->ldH, d->size_G, cMT, PETSC_FALSE);
575: d->size_G = cMT;
576: d->cX_in_G = cp+rm;
577: }
579: return(0);
580: }
585: /* DvdMult_copy_func: 2 */
586: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
587: DvdMult_copy_func **sr)
588: {
589: PetscErrorCode ierr;
590: Vec *W = d->W?d->W:d->V, *BV = d->BV?d->BV:d->V;
594: if (!d->G || d->V_new_s == d->V_new_e) { return(0); }
596: /* G = [G W(old)'*BV(new);
597: W(new)'*BV(old) W(new)'*BV(new) ],
598: being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
599: VecsMultS(d->G,d->sG,d->ldH,W-d->cX_in_G,d->V_new_s+d->cX_in_G,d->V_new_e+d->cX_in_G,BV-d->cX_in_G,d->V_new_s+d->cX_in_G,d->V_new_e+d->cX_in_G,r,(*sr)++);
600: d->size_G = d->V_new_e+d->cX_in_G;
602: return(0);
603: }
605: /* in complex, d->size_H real auxiliar values are needed */
608: PetscErrorCode dvd_calcpairs_projeig_eig(dvdDashboard *d)
609: {
610: PetscErrorCode ierr;
611: PetscReal *w;
612: #if defined(PETSC_USE_COMPLEX)
613: PetscInt i;
614: #endif
618: /* S <- H */
619: d->ldS = d->ldpX = d->size_H;
620: SlepcDenseCopyTriang(d->S, DVD_MAT_LTRIANG, d->size_H, d->H, d->sH, d->ldH,
621: d->size_H, d->size_H);
623: /* S = pX' * L * pX */
624: #if !defined(PETSC_USE_COMPLEX)
625: w = d->eigr-d->cX_in_H;
626: #else
627: w = (PetscReal*)d->auxS;
628: #endif
629: EPSDenseHEP(d->size_H, d->S, d->ldS, w, d->pX);
630: #if defined(PETSC_USE_COMPLEX)
631: for (i=0; i<d->size_H; i++) d->eigr[i-d->cX_in_H] = w[i];
632: #endif
634: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_eig;
636: return(0);
637: }
642: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d)
643: {
644: PetscErrorCode ierr;
648: /* S <- H */
649: d->ldS = d->ldpX = d->size_H;
650: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
651: d->size_H, d->size_H);
653: /* S = pX' * H * pX */
654: EPSDenseHessenberg(d->size_H, 0, d->S, d->ldS, d->pX);
655: EPSDenseSchur(d->size_H, 0, d->S, d->ldS, d->pX, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
656:
658: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;
660: return(0);
661: }
665: /*
666: auxS(dgges) = size_H (beta) + 8*size_H+16 (work)
667: auxS(zgges) = size_H (beta) + 1+2*size_H (work) + 8*size_H (rwork)
668: */
669: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d)
670: {
671: #if defined(SLEPC_MISSING_LAPACK_GGES)
673: SETERRQ(((PetscObject)(d->eps))->comm,PETSC_ERR_SUP,"GGES - Lapack routine is unavailable.");
674: #else
675: PetscErrorCode ierr;
676: PetscScalar *beta = d->auxS;
677: #if !defined(PETSC_USE_COMPLEX)
678: PetscScalar *auxS = beta + d->size_H;
679: PetscBLASInt n_auxS = d->size_auxS - d->size_H;
680: #else
681: PetscReal *auxR = (PetscReal*)(beta + d->size_H);
682: PetscScalar *auxS = (PetscScalar*)(auxR+8*d->size_H);
683: PetscBLASInt n_auxS = d->size_auxS - 9*d->size_H;
684: #endif
685: PetscInt i;
686: PetscBLASInt info,n,a;
689: /* S <- H, T <- G */
690: d->ldS = d->ldT = d->ldpX = d->ldpY = d->size_H;
691: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
692: d->size_H, d->size_H);
693: SlepcDenseCopyTriang(d->T, 0, d->size_H, d->G, d->sG, d->ldH,
694: d->size_H, d->size_H);
696: /* S = Z'*H*Q, T = Z'*G*Q */
697: n = d->size_H;
698: #if !defined(PETSC_USE_COMPLEX)
699: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
700: &a, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
701: auxS, &n_auxS, PETSC_NULL, &info);
702: #else
703: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
704: &a, d->eigr-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
705: auxS, &n_auxS, auxR, PETSC_NULL, &info);
706: #endif
707: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack GGES %d", info);
709: /* eigr[i] <- eigr[i] / beta[i] */
710: for (i=0; i<d->size_H; i++)
711: d->eigr[i-d->cX_in_H] /= beta[i],
712: d->eigi[i-d->cX_in_H] /= beta[i];
714: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;
716: return(0);
717: #endif
718: }
722: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n)
723: {
724: PetscErrorCode ierr;
728: EPSSortDenseHEP(d->eps, d->size_H, 0, d->eigr-d->cX_in_H, d->pX, d->ldpX);
729:
731: if (d->calcpairs_eigs_trans) {
732: d->calcpairs_eigs_trans(d);
733: }
735: return(0);
736: }
741: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n)
742: {
743: PetscErrorCode ierr;
744: #if defined(PETSC_USE_COMPLEX)
745: PetscScalar s;
746: PetscInt i, j;
747: #endif
750: if ((d->ldpX != d->size_H) ||
751: ( d->T &&
752: ((d->ldS != d->ldT) || (d->ldpX != d->ldpY) ||
753: (d->ldpX != d->size_H)) ) ) {
754: SETERRQ(PETSC_COMM_SELF,1, "Error before ordering eigenpairs");
755: }
757: if (d->T) {
758: EPSSortDenseSchurGeneralized(d->eps, d->size_H, 0, n, d->S, d->T,
759: d->ldS, d->pY, d->pX, d->eigr-d->cX_in_H,
760: d->eigi-d->cX_in_H);
761: } else {
762: EPSSortDenseSchur(d->eps, d->size_H, 0, d->S, d->ldS, d->pX,
763: d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
764: }
766: if (d->calcpairs_eigs_trans) {
767: d->calcpairs_eigs_trans(d);
768: }
770: /* Some functions need the diagonal elements in T be real */
771: #if defined(PETSC_USE_COMPLEX)
772: if (d->T) for(i=0; i<d->size_H; i++)
773: if (PetscImaginaryPart(d->T[d->ldT*i+i]) != 0.0) {
774: s = PetscConj(d->T[d->ldT*i+i])/PetscAbsScalar(d->T[d->ldT*i+i]);
775: for(j=0; j<=i; j++)
776: d->T[d->ldT*i+j] = PetscRealPart(d->T[d->ldT*i+j]*s),
777: d->S[d->ldS*i+j]*= s;
778: for(j=0; j<d->size_H; j++) d->pX[d->ldpX*i+j]*= s;
779: }
780: #endif
782: return(0);
783: }
788: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
789: the norm, where i = r_s..r_e
790: */
791: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
792: Vec *R)
793: {
794: PetscInt i;
795: PetscErrorCode ierr;
796: Vec *BV = d->BV?d->BV:d->V;
800: for (i=r_s; i<r_e; i++) {
801: /* nX(i) <- ||X(i)|| */
802: if (d->correctXnorm) {
803: /* R(i) <- V*pX(i) */
804: SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
805: &d->V[-d->cX_in_H], d->size_V+d->cX_in_H,
806: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
807: VecNorm(R[i-r_s], NORM_2, &d->nX[i]);
808: } else
809: d->nX[i] = 1.0;
811: /* R(i-r_s) <- AV*pX(i) */
812: SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
813: &d->AV[-d->cX_in_H], d->size_AV+d->cX_in_H,
814: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
816: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BV*pX(i) */
817: SlepcUpdateVectorsZ(&R[i-r_s], 1.0, -d->eigr[i+d->cX_in_H],
818: &BV[-d->cX_in_H], d->size_V+d->cX_in_H,
819: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
820: }
822: d->calcpairs_proj_res(d, r_s, r_e, R);
824: return(0);
825: }
829: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
830: PetscInt r_e, Vec *R)
831: {
832: PetscInt i;
833: PetscErrorCode ierr;
834: PetscBool lindep;
835: Vec *cX;
839: /* If exists the BcX, R <- orth(BcX, R), nR[i] <- ||R[i]|| */
840: if (d->BcX)
841: cX = d->BcX;
843: /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
844: else if (d->cY)
845: cX = d->cY;
847: /* If fany configurations, R <- orth(cX, R), nR[i] <- ||R[i]|| */
848: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN)))
849: cX = d->cX;
851: /* Otherwise, nR[i] <- ||R[i]|| */
852: else
853: cX = PETSC_NULL;
855: if (cX) for (i=0; i<r_e-r_s; i++) {
856: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX, PETSC_NULL,
857: cX, R[i], PETSC_NULL, &d->nR[r_s+i], &lindep);
858:
859: if(lindep || (d->nR[r_s+i] < PETSC_MACHINE_EPSILON)) {
860: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %G!\n",r_s+i,d->nR[r_s+i]);
861: }
863: } else for(i=0; i<r_e-r_s; i++) {
864: VecNorm(R[i], NORM_2, &d->nR[r_s+i]);
865: }
867: return(0);
868: }
870: /**** Pattern routines ********************************************************/
872: /* BV <- BV*MTX */
875: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,Vec *real_BV,PetscInt *size_cBV,Vec **BV,PetscInt *size_BV,PetscInt *max_size_BV,PetscBool BV_shift,PetscInt *cX_in_proj,PetscScalar *MTX)
876: {
877: PetscErrorCode ierr;
878: PetscInt cMT, rm, cp, tra_s, i;
879: Vec *nBV;
883: if (!real_BV || !*BV || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }
885: if (d->V_tra_s > d->max_cX_in_proj && !BV_shift) {
886: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
887: cMT = d->V_tra_e - tra_s;
888: rm = d->V_tra_s - tra_s;
889: cp = PetscMin(d->max_cX_in_proj - rm, *cX_in_proj);
890: nBV = real_BV+d->max_cX_in_proj;
891: /* BV(-cp-rm:-1-rm) <- BV(-cp:-1) */
892: for (i=-cp; i<0; i++) {
893: VecCopy((*BV)[i], nBV[i-rm]);
894: }
895: /* BV(-rm:) <- BV*MTX(tra_s:V_tra_e-1) */
896: SlepcUpdateVectorsZ(&nBV[-rm], 0.0, 1.0, *BV-*cX_in_proj, *size_BV+*cX_in_proj, &MTX[d->ldMTX*tra_s], d->ldMTX, d->size_MT, cMT);
897: *size_BV = d->V_tra_e - d->V_tra_s;
898: *max_size_BV-= nBV - *BV;
899: *BV = nBV;
900: if (cX_in_proj && d->max_cX_in_proj>0) *cX_in_proj = cp+rm;
901: } else if (d->V_tra_s <= d->max_cX_in_proj || BV_shift) {
902: /* [BcX BV] <- [BcX BV*MTX] */
903: SlepcUpdateVectorsZ(*BV-*cX_in_proj, 0.0, 1.0, *BV-*cX_in_proj, *size_BV+*cX_in_proj, MTX, d->ldMTX, d->size_MT, d->V_tra_e);
904: *BV+= d->V_tra_s-*cX_in_proj;
905: *max_size_BV-= d->V_tra_s-*cX_in_proj;
906: *size_BV = d->V_tra_e - d->V_tra_s;
907: if (size_cBV && BV_shift) *size_cBV = *BV - real_BV;
908: if (d->max_cX_in_proj>0) *cX_in_proj = PetscMin(*BV - real_BV, d->max_cX_in_proj);
909: } else { /* !BV_shift */
910: /* BV <- BV*MTX(V_tra_s:) */
911: SlepcUpdateVectorsZ(*BV, 0.0, 1.0, *BV, *size_BV,
912: &MTX[d->V_tra_s*d->ldMTX], d->ldMTX, d->size_MT, d->V_tra_e-d->V_tra_s);
913:
914: *size_BV = d->V_tra_e - d->V_tra_s;
915: }
917: return(0);
918: }