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-2012, 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_qz_d(dvdDashboard *d);
39: PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d);
40: PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d, PetscInt n);
41: PetscErrorCode dvd_calcpairs_X(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
42: Vec *X);
43: PetscErrorCode dvd_calcpairs_Y(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
44: Vec *Y);
45: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
46: Vec *R);
47: PetscErrorCode dvd_calcpairs_eig_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R);
48: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
49: PetscInt r_e, Vec *R);
50: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
51: DvdMult_copy_func **sr);
52: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d);
53: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
54: DvdMult_copy_func **sr);
55: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d);
56: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d);
57: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
58: DvdMult_copy_func **sr);
59: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d);
60: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
61: DvdMult_copy_func **sr);
62: 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,DSMatType MT);
64: /**** Control routines ********************************************************/
67: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d, dvdBlackboard *b,
68: EPSOrthType orth, IP ipI,
69: PetscInt cX_proj, PetscBool harm)
70: {
71: PetscErrorCode ierr;
72: PetscInt i,max_cS;
73: PetscBool std_probl,her_probl,ind_probl,her_ind_probl;
74: const DSType dstype;
75: const char *prefix;
76: PetscErrorCode (*f)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
77: void *ctx;
81: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
82: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
83: ind_probl = DVD_IS(d->sEP, DVD_EP_INDEFINITE)?PETSC_TRUE:PETSC_FALSE;
84: her_ind_probl = (her_probl || ind_probl)? PETSC_TRUE:PETSC_FALSE;
86: /* Setting configuration constrains */
87: #if !defined(PETSC_USE_COMPLEX)
88: /* if the last converged eigenvalue is complex its conjugate pair is also
89: converged */
90: b->max_nev = PetscMax(b->max_nev, d->nev+(her_probl && !d->B?0:1));
91: #else
92: b->max_nev = PetscMax(b->max_nev, d->nev);
93: #endif
94: b->max_size_proj = PetscMax(b->max_size_proj, b->max_size_V+cX_proj);
95: d->size_real_V = b->max_size_V+b->max_nev;
96: d->W_shift = d->B?PETSC_TRUE:PETSC_FALSE;
97: d->size_real_W = harm?(b->max_size_V+(d->W_shift?b->max_nev:b->max_size_cP)):0;
98: d->size_real_AV = b->max_size_V+b->max_size_cP;
99: d->size_BDS = 0;
100: if (d->B && her_ind_probl && (orth == EPS_ORTH_I || orth == EPS_ORTH_BOPT)) {
101: d->size_real_BV = b->size_V; d->BV_shift = PETSC_TRUE;
102: if (orth == EPS_ORTH_BOPT) d->size_BDS = d->eps->nds;
103: } else if (d->B) {
104: d->size_real_BV = b->max_size_V + b->max_size_P; d->BV_shift = PETSC_FALSE;
105: } else {
106: d->size_real_BV = 0; d->BV_shift = PETSC_FALSE;
107: }
108: b->own_vecs+= d->size_real_V + d->size_real_W + d->size_real_AV +
109: d->size_real_BV + d->size_BDS;
110: b->own_scalars+= b->max_size_proj*b->max_size_proj*2*(std_probl?1:2) +
111: /* H, G?, S, T? */
112: b->max_nev*b->max_nev*(her_ind_probl?0:(!d->B?1:2)) +
113: /* cS?, cT? */
114: FromRealToScalar(d->size_real_V)*(ind_probl?1:0) + /* nBV */
115: FromRealToScalar(b->max_size_proj)*(ind_probl?1:0) + /* nBpX */
116: (d->eps->arbit_func? b->size_V*2 : 0); /* rr, ri */
117: b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X);
118: /* updateV0 */
119: max_cS = PetscMax(b->max_size_X,cX_proj);
120: b->max_size_auxS = PetscMax(PetscMax(
121: b->max_size_auxS,
122: b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
123: max_cS*b->max_nev*(her_ind_probl?0:(!d->B?1:2)) + /* updateV0,W0 */
124: /* SlepcReduction: in */
125: PetscMax(
126: b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
127: max_cS*b->max_nev*(her_ind_probl?0:(!d->B?1:2)), /* updateV0,W0 */
128: /* SlepcReduction: out */
129: PetscMax(
130: b->max_size_proj*b->max_size_proj, /* updateAV0,BV0 */
131: b->max_size_proj+b->max_nev))), /* dvd_orth */
132: std_probl?0:(b->max_size_proj*11+16) /* projeig */);
133: #if defined(PETSC_USE_COMPLEX)
134: b->max_size_auxS = PetscMax(b->max_size_auxS, b->max_size_V);
135: /* dvd_calcpairs_projeig_eig */
136: #endif
138: /* Setup the step */
139: if (b->state >= DVD_STATE_CONF) {
140: d->max_cX_in_proj = cX_proj;
141: d->max_size_P = b->max_size_P;
142: d->real_V = b->free_vecs; b->free_vecs+= d->size_real_V;
143: if (harm) {
144: d->real_W = b->free_vecs; b->free_vecs+= d->size_real_W;
145: } else {
146: d->real_W = PETSC_NULL;
147: }
148: d->real_AV = d->AV = b->free_vecs; b->free_vecs+= d->size_real_AV;
149: d->max_size_proj = b->max_size_proj;
150: d->real_H = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
151: d->ldH = b->max_size_proj;
152: d->S = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
153: if (!her_ind_probl) {
154: d->cS = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
155: d->max_size_cS = d->ldcS = b->max_nev;
156: } else {
157: d->cS = PETSC_NULL;
158: d->max_size_cS = d->ldcS = 0;
159: d->orthoV_type = orth;
160: if (ind_probl) {
161: d->real_nBV = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(d->size_real_V);
162: d->nBpX = (PetscReal*)b->free_scalars; b->free_scalars+= FromRealToScalar(d->max_size_proj);
163: } else
164: d->real_nBV = d->nBDS = d->nBpX = PETSC_NULL;
165: }
166: d->ipV = ipI;
167: d->ipW = ipI;
168: if (orth == EPS_ORTH_BOPT) {
169: d->BDS = b->free_vecs; b->free_vecs+= d->eps->nds;
170: for (i=0; i<d->eps->nds; i++) {
171: MatMult(d->B, d->eps->defl[i], d->BDS[i]);
172: }
173: } else
174: d->BDS = PETSC_NULL;
175: if (d->B) {
176: d->real_BV = b->free_vecs; b->free_vecs+= d->size_real_BV;
177: } else {
178: d->size_real_BV = 0;
179: d->real_BV = PETSC_NULL;
180: d->BV_shift = PETSC_FALSE;
181: }
182: if (!std_probl) {
183: d->real_G = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
184: d->T = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
185: } else {
186: d->real_G = PETSC_NULL;
187: d->T = PETSC_NULL;
188: }
189: if (d->B && !her_ind_probl) {
190: d->cT = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
191: d->ldcT = b->max_nev;
192: } else {
193: d->cT = PETSC_NULL;
194: d->ldcT = 0;
195: }
196: if (d->eps->arbit_func) {
197: d->eps->rr = b->free_scalars; b->free_scalars+= b->size_V;
198: d->eps->ri = b->free_scalars; b->free_scalars+= b->size_V;
199: } else {
200: d->eps->rr = PETSC_NULL;
201: d->eps->ri = PETSC_NULL;
202: }
203: /* Create a DS if the method works with Schur decompositions */
204: if (d->cS) {
205: DSCreate(((PetscObject)d->eps->ds)->comm,&d->conv_ps);
206: DSSetType(d->conv_ps,d->cT ? DSGNHEP : DSNHEP);
207: /* Transfer as much as possible options from eps->ds to conv_ps */
208: DSGetOptionsPrefix(d->eps->ds,&prefix);
209: DSSetOptionsPrefix(d->conv_ps,prefix);
210: DSSetFromOptions(d->conv_ps);
211: DSGetEigenvalueComparison(d->eps->ds,&f,&ctx);
212: DSSetEigenvalueComparison(d->conv_ps,f,ctx);
213: DSAllocate(d->conv_ps,b->max_nev);
214: } else {
215: d->conv_ps = PETSC_NULL;
216: }
217: d->calcPairs = dvd_calcpairs_proj;
218: d->calcpairs_residual = dvd_calcpairs_res_0;
219: d->calcpairs_residual_eig = dvd_calcpairs_eig_res_0;
220: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
221: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs;
222: d->ipI = ipI;
223: /* Create and configure a DS for solving the projected problems */
224: if (d->real_W) { /* If we use harmonics */
225: dstype = DSGNHEP;
226: } else {
227: if (ind_probl) {
228: dstype = DSGHIEP;
229: } else if (std_probl) {
230: dstype = her_probl ? DSHEP : DSNHEP;
231: } else {
232: dstype = her_probl ? DSGHEP : DSGNHEP;
233: }
234: }
235: d->ps = d->eps->ds;
236: DSSetType(d->ps,dstype);
237: DSAllocate(d->ps,d->max_size_proj);
239: DVD_FL_ADD(d->startList, dvd_calcpairs_qz_start);
240: DVD_FL_ADD(d->destroyList, dvd_calcpairs_qz_d);
241: }
243: return(0);
244: }
249: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
250: {
251: PetscBool her_probl,ind_probl,her_ind_probl;
252: PetscInt i;
255: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
256: ind_probl = DVD_IS(d->sEP, DVD_EP_INDEFINITE)?PETSC_TRUE:PETSC_FALSE;
257: her_ind_probl = (her_probl || ind_probl)? PETSC_TRUE:PETSC_FALSE;
259: d->size_V = 0;
260: d->V = d->real_V;
261: d->cX = d->real_V;
262: d->size_cX = 0;
263: d->max_size_V = d->size_real_V;
264: d->W = d->real_W;
265: d->max_size_W = d->size_real_W;
266: d->size_W = 0;
267: d->size_AV = 0;
268: d->AV = d->real_AV;
269: d->max_size_AV = d->size_real_AV;
270: d->size_H = 0;
271: d->H = d->real_H;
272: if (d->cS) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cS[i] = 0.0;
273: d->size_BV = 0;
274: d->BV = d->real_BV;
275: d->max_size_BV = d->size_real_BV;
276: d->size_G = 0;
277: d->G = d->real_G;
278: if (d->cT) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cT[i] = 0.0;
279: d->cY = d->B && !her_ind_probl ? d->W : PETSC_NULL;
280: d->BcX = d->orthoV_type == EPS_ORTH_I && d->B && her_probl ? d->BcX : PETSC_NULL;
281: d->size_cY = 0;
282: d->size_BcX = 0;
283: 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;
284: d->nBV = d->nBcX = d->real_nBV;
285: return(0);
286: }
290: PetscErrorCode dvd_calcpairs_qz_d(dvdDashboard *d)
291: {
292: PetscErrorCode ierr;
295: DSDestroy(&d->conv_ps);
296: return(0);
297: }
301: PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
302: {
303: PetscErrorCode ierr;
304: DvdReduction r;
305: #define MAX_OPS 7
306: DvdReductionChunk
307: ops[MAX_OPS];
308: DvdMult_copy_func
309: sr[MAX_OPS], *sr0 = sr;
310: PetscInt size_in, i;
311: PetscScalar *in = d->auxS, *out;
312: PetscBool stdp;
316: stdp = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
317: size_in =
318: (d->size_cX+d->V_tra_s-d->cX_in_H)*d->V_tra_s*(d->cT?2:(d->cS?1:0)) + /* updateV0,W0 */
319: (d->size_H*(d->V_new_e-d->V_new_s)*2+
320: (d->V_new_e-d->V_new_s)*(d->V_new_e-d->V_new_s))*(!stdp?2:1); /* updateAV1,BV1 */
321:
322: out = in+size_in;
324: /* Check consistency */
325: if (2*size_in > d->size_auxS) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
327: /* Prepare reductions */
328: SlepcAllReduceSumBegin(ops, MAX_OPS, in, out, size_in, &r,
329: ((PetscObject)d->V[0])->comm);
330: /* Allocate size_in */
331: d->auxS+= size_in;
332: d->size_auxS-= size_in;
334: /* Update AV, BV, W and the projected matrices */
335: /* 1. S <- S*MT */
336: dvd_calcpairs_updateV0(d, &r, &sr0);
337: dvd_calcpairs_updateW0(d, &r, &sr0);
338: dvd_calcpairs_updateAV0(d);
339: dvd_calcpairs_updateBV0(d);
340: /* 2. V <- orth(V, V_new) */
341: dvd_calcpairs_updateV1(d);
342: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
343: /* Check consistency */
344: if (d->size_AV != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
345: for (i=d->V_new_s; i<d->V_new_e; i++) {
346: MatMult(d->A, d->V[i], d->AV[i]);
347: }
348: d->size_AV = d->V_new_e;
349: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
350: if (d->B && d->orthoV_type != EPS_ORTH_BOPT) {
351: /* Check consistency */
352: if (d->size_BV != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
353: for (i=d->V_new_s; i<d->V_new_e; i++) {
354: MatMult(d->B, d->V[i], d->BV[i]);
355: }
356: d->size_BV = d->V_new_e;
357: }
358: /* 5 <- W <- [W f(AV,BV)] */
359: dvd_calcpairs_updateW1(d);
360: dvd_calcpairs_updateAV1(d, &r, &sr0);
361: dvd_calcpairs_updateBV1(d, &r, &sr0);
363: /* Deallocate size_in */
364: d->auxS-= size_in;
365: d->size_auxS+= size_in;
367: /* Do reductions */
368: SlepcAllReduceSumEnd(&r);
370: /* Perform the transformation on the projected problem */
371: if (d->calcpairs_proj_trans) {
372: d->calcpairs_proj_trans(d);
373: }
375: d->V_tra_s = d->V_tra_e = 0;
376: d->V_new_s = d->V_new_e;
378: /* Solve the projected problem */
379: if (d->size_H>0) {
380: dvd_calcpairs_projeig_solve(d);
381: }
383: /* Check consistency */
384: 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 ||
385: d->size_V != d->size_AV || d->cX_in_H != d->cX_in_AV ||
386: (DVD_ISNOT(d->sEP, DVD_EP_STD) && (
387: d->size_V+d->cX_in_G != d->size_G || d->cX_in_H != d->cX_in_G ||
388: d->size_H != d->size_G || (d->BV && (
389: d->size_V != d->size_BV || d->cX_in_H != d->cX_in_BV)))) ||
390: (d->W && d->size_W != d->size_V)) {
391: SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
392: }
394: return(0);
395: #undef MAX_OPS
396: }
398: /**** Basic routines **********************************************************/
402: /* auxV: V_tra_s, DvdMult_copy_func: 1 */
403: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
404: DvdMult_copy_func **sr)
405: {
406: PetscErrorCode ierr;
407: PetscInt rm,i,ld;
408: PetscScalar *pQ;
412: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
414: /* Update nBcX and nBV */
415: if (d->nBcX && d->nBpX && d->nBV) {
416: d->nBV+= d->V_tra_s;
417: for (i=0; i<d->V_tra_s; i++) d->nBcX[d->size_cX+i] = d->nBpX[i];
418: for (i=d->V_tra_s; i<d->V_tra_e; i++) d->nBV[i-d->V_tra_s] = d->nBpX[i];
419: }
421: /* cX <- [cX V*MT(0:V_tra_s-1)], V <- V*MT(V_tra_s:V_tra_e) */
422: 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,DS_MAT_Q);
424: /* Udpate cS for standard problems */
425: if (d->cS && !d->cT && !d->cY && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
426: /* Check consistency */
427: if (d->size_cS+d->V_tra_s != d->size_cX) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
429: /* auxV <- AV * ps.Q(0:V_tra_e-1) */
430: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
431: DSGetLeadingDimension(d->ps,&ld);
432: DSGetArray(d->ps,DS_MAT_Q,&pQ);
433: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->AV-d->cX_in_AV,d->size_AV+d->cX_in_AV,pQ,ld,d->size_MT,d->V_tra_s-rm);
434: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
436: /* cS(:, size_cS:) <- cX' * auxV */
437: 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)++);
438: d->size_cS+= d->V_tra_s-rm;
439: }
441: return(0);
442: }
447: /* auxS: size_cX+V_new_e+1 */
448: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d)
449: {
450: PetscErrorCode ierr;
451: Vec *cX = d->BcX? d->BcX : ( (d->cY && !d->W)? d->cY : d->cX );
455: if (d->V_new_s == d->V_new_e) { return(0); }
457: /* Check consistency */
458: if (d->size_V != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
460: /* V <- gs([cX V(0:V_new_s-1)], V(V_new_s:V_new_e-1)) */
461: if (d->orthoV_type == EPS_ORTH_BOPT) {
462: dvd_BorthV_faster(d->ipV,d->eps->defl,d->BDS,d->nBDS,d->eps->nds,d->cX,d->real_BV,d->nBcX,d->size_cX,d->V,d->BV,d->nBV,d->V_new_s,d->V_new_e,d->auxS,d->eps->rand);
463: d->size_BV = d->V_new_e;
464: } else if (DVD_IS(d->sEP, DVD_EP_INDEFINITE)) {
465: dvd_BorthV_stable(d->ipV,d->eps->defl,d->nBDS,d->eps->nds,d->cX,d->nBcX,d->size_cX,d->V,d->nBV,d->V_new_s,d->V_new_e,d->auxS,d->eps->rand);
466: } else {
467: dvd_orthV(d->ipV,d->eps->defl,d->eps->nds,cX,d->size_cX,d->V,d->V_new_s,d->V_new_e,d->auxS,d->eps->rand);
468: }
469: d->size_V = d->V_new_e;
471: return(0);
472: }
476: /* auxV: V_tra_s, DvdMult_copy_func: 2 */
477: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
478: DvdMult_copy_func **sr)
479: {
480: PetscErrorCode ierr;
481: PetscInt rm,ld;
482: PetscScalar *pQ;
486: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
488: /* cY <- [cY W*ps.Z(0:V_tra_s-1)], W <- W*ps.Z(V_tra_s:V_tra_e) */
489: 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,DS_MAT_Z);
491: /* Udpate cS and cT */
492: if (d->cT && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
493: /* Check consistency */
494: if (d->size_cS+d->V_tra_s != d->size_cX || (d->W && d->size_cY != d->size_cX)) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
496: DSGetLeadingDimension(d->ps,&ld);
497: DSGetArray(d->ps,DS_MAT_Q,&pQ);
498: /* auxV <- AV * ps.Q(0:V_tra_e-1) */
499: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
500: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->AV-d->cX_in_H,d->size_AV-d->cX_in_H,pQ,ld,d->size_MT,d->V_tra_s-rm);
502: /* cS(:, size_cS:) <- cY' * auxV */
503: VecsMultS(&d->cS[d->ldcS*d->size_cS], 0, d->ldcS, d->cY?d->cY:d->cX, 0, d->size_cX-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
505: /* auxV <- BV * ps.Q(0:V_tra_e-1) */
506: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->BV-d->cX_in_H,d->size_BV-d->cX_in_H,pQ,ld,d->size_MT,d->V_tra_s-rm);
507: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
509: /* cT(:, size_cS:) <- cY' * auxV */
510: VecsMultS(&d->cT[d->ldcS*d->size_cS], 0, d->ldcS, d->cY?d->cY:d->cX, 0, d->size_cX-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
511:
512: d->size_cS+= d->V_tra_s-rm;
513: d->size_cT+= d->V_tra_s-rm;
514: }
516: return(0);
517: }
522: /* auxS: size_cX+V_new_e+1 */
523: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d)
524: {
525: PetscErrorCode ierr;
526: Vec *cY = d->cY?d->cY:d->cX;
530: if (!d->W || d->V_new_s == d->V_new_e) { return(0); }
532: /* Check consistency */
533: if (d->size_W != d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
535: /* Update W */
536: d->calcpairs_W(d);
538: /* W <- gs([cY W(0:V_new_s-1)], W(V_new_s:V_new_e-1)) */
539: dvd_orthV(d->ipW, PETSC_NULL, 0, cY, d->size_cX, d->W-d->cX_in_W, d->V_new_s+d->cX_in_W, d->V_new_e+d->cX_in_W, d->auxS, d->eps->rand);
540: d->size_W = d->V_new_e;
542: return(0);
543: }
547: /* auxS: size_H*(V_tra_e-V_tra_s) */
548: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d)
549: {
550: PetscErrorCode ierr;
551: PetscInt cMT,tra_s,ld;
552: PetscScalar *pQ,*pZ;
556: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
558: /* AV(V_tra_s-cp-1:) = cAV*ps.Q(V_tra_s:) */
559: 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,DS_MAT_Q);
560: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj,0);
561: cMT = d->V_tra_e - tra_s;
563: /* Update H <- ps.Z(tra_s)' * (H * ps.Q(tra_s:)) */
564: DSGetLeadingDimension(d->ps,&ld);
565: DSGetArray(d->ps,DS_MAT_Q,&pQ);
566: if (d->W) {
567: DSGetArray(d->ps,DS_MAT_Z,&pZ);
568: } else {
569: pZ = pQ;
570: }
571: SlepcDenseMatProdTriang(d->auxS,0,d->ldH,d->H,d->sH,d->ldH,d->size_H,d->size_H,PETSC_FALSE,&pQ[ld*tra_s],0,ld,d->size_MT,cMT,PETSC_FALSE);
572: SlepcDenseMatProdTriang(d->H,d->sH,d->ldH,&pZ[ld*tra_s],0,ld,d->size_MT,cMT,PETSC_TRUE,d->auxS,0,d->ldH,d->size_H,cMT,PETSC_FALSE);
573: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
574: if (d->W) {
575: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
576: }
577: d->size_H = cMT;
578: d->cX_in_H = d->cX_in_AV;
580: return(0);
581: }
586: /* DvdMult_copy_func: 2 */
587: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
588: DvdMult_copy_func **sr)
589: {
590: PetscErrorCode ierr;
591: Vec *W = d->W?d->W:d->V;
595: if (d->V_new_s == d->V_new_e) { return(0); }
597: /* Check consistency */
598: if (d->size_H != d->V_new_s+d->cX_in_H || d->size_V != d->V_new_e) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
600: /* H = [H W(old)'*AV(new);
601: W(new)'*AV(old) W(new)'*AV(new) ],
602: being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
603: 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)++);
604: d->size_H = d->V_new_e+d->cX_in_H;
606: return(0);
607: }
611: /* auxS: max(BcX*(size_cX+V_new_e+1), size_G*(V_tra_e-V_tra_s)) */
612: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d)
613: {
614: PetscErrorCode ierr;
615: PetscInt cMT,tra_s,i,ld;
616: PetscBool lindep;
617: PetscReal norm;
618: PetscScalar *pQ,*pZ;
622: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
624: /* BV <- BV*MT */
625: 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,DS_MAT_Q);
627: /* If BcX, BcX <- orth(BcX) */
628: if (d->BcX) {
629: for (i=0; i<d->V_tra_s; i++) {
630: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_BcX+i, PETSC_NULL,
631: d->BcX, d->BcX[d->size_BcX+i], PETSC_NULL,
632: &norm, &lindep);
633: if(lindep) SETERRQ(PETSC_COMM_SELF,1, "Error during orth(BcX, B*cX(new))");
634: VecScale(d->BcX[d->size_BcX+i], 1.0/norm);
635: }
636: d->size_BcX+= d->V_tra_s;
637: }
639: /* Update G <- ps.Z' * (G * ps.Q) */
640: if (d->G) {
641: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj,0);
642: cMT = d->V_tra_e - tra_s;
643: DSGetLeadingDimension(d->ps,&ld);
644: DSGetArray(d->ps,DS_MAT_Q,&pQ);
645: if (d->W) {
646: DSGetArray(d->ps,DS_MAT_Z,&pZ);
647: } else {
648: pZ = pQ;
649: }
650: SlepcDenseMatProdTriang(d->auxS,0,d->ldH,d->G,d->sG,d->ldH,d->size_G,d->size_G,PETSC_FALSE,&pQ[ld*tra_s],0,ld,d->size_MT,cMT,PETSC_FALSE);
651: SlepcDenseMatProdTriang(d->G,d->sG,d->ldH,&pZ[ld*tra_s],0,ld,d->size_MT,cMT,PETSC_TRUE,d->auxS,0,d->ldH,d->size_G,cMT,PETSC_FALSE);
652: DSRestoreArray(d->ps,DS_MAT_Q,&pQ);
653: if (d->W) {
654: DSRestoreArray(d->ps,DS_MAT_Z,&pZ);
655: }
656: d->size_G = cMT;
657: d->cX_in_G = d->cX_in_V;
658: }
660: return(0);
661: }
666: /* DvdMult_copy_func: 2 */
667: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
668: DvdMult_copy_func **sr)
669: {
670: PetscErrorCode ierr;
671: Vec *W = d->W?d->W:d->V, *BV = d->BV?d->BV:d->V;
675: if (!d->G || d->V_new_s == d->V_new_e) { return(0); }
677: /* G = [G W(old)'*BV(new);
678: W(new)'*BV(old) W(new)'*BV(new) ],
679: being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
680: 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)++);
681: d->size_G = d->V_new_e+d->cX_in_G;
683: return(0);
684: }
686: /* in complex, d->size_H real auxiliar values are needed */
689: PetscErrorCode dvd_calcpairs_projeig_solve(dvdDashboard *d)
690: {
691: PetscErrorCode ierr;
692: PetscScalar *A;
693: PetscInt ld,i;
696: DSSetDimensions(d->ps,d->size_H,PETSC_IGNORE,0,0);
697: DSGetLeadingDimension(d->ps,&ld);
698: DSGetArray(d->ps,DS_MAT_A,&A);
699: SlepcDenseCopyTriang(A,0,ld,d->H,d->sH,d->ldH,d->size_H,d->size_H);
700: DSRestoreArray(d->ps,DS_MAT_A,&A);
701: if (d->G) {
702: DSGetArray(d->ps,DS_MAT_B,&A);
703: SlepcDenseCopyTriang(A,0,ld,d->G,d->sG,d->ldH,d->size_H,d->size_H);
704: DSRestoreArray(d->ps,DS_MAT_B,&A);
705: }
706: /* Set the signature on projected matrix B */
707: if (DVD_IS(d->sEP, DVD_EP_INDEFINITE)) {
708: DSGetArray(d->ps,DS_MAT_B,&A);
709: PetscMemzero(A,sizeof(PetscScalar)*d->size_H*ld);
710: for (i=0; i<d->size_H; i++) {
711: A[i+ld*i] = d->nBV[i];
712: }
713: DSRestoreArray(d->ps,DS_MAT_B,&A);
714: }
715: DSSetState(d->ps,DS_STATE_RAW);
716: DSSolve(d->ps,d->eigr-d->cX_in_H,d->eigi-d->cX_in_H);
717: return(0);
718: }
722: PetscErrorCode dvd_calcpairs_apply_arbitrary_func(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscScalar **rr_,PetscScalar **ri_)
723: {
724: PetscInt i,k,ld;
725: PetscScalar *pX,*rr,*ri,ar,ai;
726: Vec *X = d->auxV,xr,xi;
727: PetscErrorCode ierr;
728: #if !defined(PETSC_USE_COMPLEX)
729: PetscInt j;
730: #endif
731:
733: /* Quick exit without neither arbitrary selection nor harmonic extraction */
734: if (!d->eps->arbit_func && !d->calcpairs_eig_backtrans) {
735: *rr_ = d->eigr-d->cX_in_H;
736: *ri_ = d->eigi-d->cX_in_H;
737: return(0);
738: }
740: /* Quick exit without arbitrary selection, but with harmonic extraction */
741: if (!d->eps->arbit_func && d->calcpairs_eig_backtrans) {
742: *rr_ = rr = d->auxS;
743: *ri_ = ri = d->auxS+r_e-r_s;
744: for (i=r_s; i<r_e; i++) {
745: d->calcpairs_eig_backtrans(d,d->eigr[i],d->eigi[i],&rr[i-r_s],&ri[i-r_s]);
746: }
747: return(0);
748: }
750: DSGetLeadingDimension(d->ps,&ld);
751: *rr_ = rr = d->eps->rr + d->eps->nconv;
752: *ri_ = ri = d->eps->ri + d->eps->nconv;
753: for (i=r_s; i<r_e; i++) {
754: k = i;
755: DSVectors(d->ps,DS_MAT_X,&k,PETSC_NULL);
756: DSNormalize(d->ps,DS_MAT_X,i);
757: DSGetArray(d->ps,DS_MAT_X,&pX);
758: dvd_improvex_compute_X(d,i,k+1,X,pX,ld);
759: DSRestoreArray(d->ps,DS_MAT_X,&pX);
760: #if !defined(PETSC_USE_COMPLEX)
761: if (d->nX[i] != 1.0) {
762: for (j=i; j<k+1; j++) {
763: VecScale(X[j-i],1/d->nX[i]);
764: }
765: }
766: xr = X[0];
767: xi = X[1];
768: if (i == k) {
769: VecZeroEntries(xi);
770: }
771: #else
772: xr = X[0];
773: xi = PETSC_NULL;
774: if (d->nX[i] != 1.0) {
775: VecScale(xr,1.0/d->nX[i]);
776: }
777: #endif
778: if (d->calcpairs_eig_backtrans) {
779: d->calcpairs_eig_backtrans(d,d->eigr[i],d->eigi[i],&ar,&ai);
780: } else {
781: ar = d->eigr[i];
782: ai = d->eigi[i];
783: }
784: (d->eps->arbit_func)(ar,ai,xr,xi,&rr[i-r_s],&ri[i-r_s],d->eps->arbit_ctx);
785: #if !defined(PETSC_USE_COMPLEX)
786: if (i != k) {
787: rr[i+1-r_s] = rr[i-r_s];
788: ri[i+1-r_s] = ri[i-r_s];
789: i++;
790: }
791: #endif
792: }
793: return(0);
794: }
795:
799: PetscErrorCode dvd_calcpairs_selectPairs(dvdDashboard *d, PetscInt n)
800: {
801: PetscInt k;
802: PetscScalar *rr,*ri;
803: PetscErrorCode ierr;
806: n = PetscMin(n,d->size_H-d->cX_in_H);
807: /* Put the best n pairs at the beginning. Useful for restarting */
808: DSSetDimensions(d->ps,PETSC_IGNORE,PETSC_IGNORE,d->cX_in_H,PETSC_IGNORE);
809: dvd_calcpairs_apply_arbitrary_func(d,d->cX_in_H,d->size_H,&rr,&ri);
810: k = n;
811: DSSort(d->ps,d->eigr-d->cX_in_H,d->eigi-d->cX_in_H,rr,ri,&k);
812: /* Put the best pair at the beginning. Useful to check its residual */
813: #if !defined(PETSC_USE_COMPLEX)
814: if (n != 1 && (n != 2 || d->eigi[0] == 0.0))
815: #else
816: if (n != 1)
817: #endif
818: {
819: dvd_calcpairs_apply_arbitrary_func(d,d->cX_in_H,d->size_H,&rr,&ri);
820: k = 1;
821: DSSort(d->ps,d->eigr-d->cX_in_H,d->eigi-d->cX_in_H,rr,ri,&k);
822: }
823: if (d->calcpairs_eigs_trans) {
824: d->calcpairs_eigs_trans(d);
825: }
826: return(0);
827: }
832: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
833: the norm associated to the Schur pair, where i = r_s..r_e
834: */
835: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
836: Vec *R)
837: {
838: PetscInt i,ldpX;
839: PetscScalar *pX;
840: PetscErrorCode ierr;
841: Vec *BV = d->BV?d->BV:d->V;
844: DSGetLeadingDimension(d->ps,&ldpX);
845: DSGetArray(d->ps,DS_MAT_Q,&pX);
846: for (i=r_s; i<r_e; i++) {
847: /* nX(i) <- ||X(i)|| */
848: if (d->correctXnorm) {
849: /* R(i) <- V*pX(i) */
850: SlepcUpdateVectorsZ(&R[i-r_s],0.0,1.0,&d->V[-d->cX_in_H],d->size_V+d->cX_in_H,&pX[ldpX*(i+d->cX_in_H)],ldpX,d->size_H,1);
851: VecNorm(R[i-r_s],NORM_2,&d->nX[i]);
852: } else {
853: d->nX[i] = 1.0;
854: }
855: /* R(i-r_s) <- AV*pX(i) */
856: SlepcUpdateVectorsZ(&R[i-r_s],0.0,1.0,&d->AV[-d->cX_in_H],d->size_AV+d->cX_in_H,&pX[ldpX*(i+d->cX_in_H)],ldpX,d->size_H,1);
857: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BV*pX(i) */
858: SlepcUpdateVectorsZ(&R[i-r_s],1.0,-d->eigr[i+d->cX_in_H],&BV[-d->cX_in_H],d->size_V+d->cX_in_H,&pX[ldpX*(i+d->cX_in_H)],ldpX,d->size_H,1);
859: }
860: DSRestoreArray(d->ps,DS_MAT_Q,&pX);
861: d->calcpairs_proj_res(d, r_s, r_e, R);
862: return(0);
863: }
867: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
868: PetscInt r_e, Vec *R)
869: {
870: PetscInt i;
871: PetscErrorCode ierr;
872: PetscBool lindep;
873: Vec *cX;
877: /* If exists the BcX, R <- orth(BcX, R), nR[i] <- ||R[i]|| */
878: if (d->BcX)
879: cX = d->BcX;
881: /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
882: else if (d->cY)
883: cX = d->cY;
885: /* If fany configurations, R <- orth(cX, R), nR[i] <- ||R[i]|| */
886: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN)))
887: cX = d->cX;
889: /* Otherwise, nR[i] <- ||R[i]|| */
890: else
891: cX = PETSC_NULL;
893: if (cX) {
894: if (cX && d->orthoV_type == EPS_ORTH_BOPT) {
895: Vec auxV;
896: VecDuplicate(d->auxV[0],&auxV);
897: for (i=0; i<r_e-r_s; i++) {
898: IPBOrthogonalize(d->ipV,d->eps->nds,d->eps->defl,d->BDS,d->nBDS,d->size_cX,PETSC_NULL,d->cX,d->real_BV,d->nBcX,R[i],auxV,PETSC_NULL,&d->nR[r_s+i],&lindep);
899: }
900: VecDestroy(&auxV);
901: } else if (DVD_IS(d->sEP, DVD_EP_INDEFINITE)) {
902: for (i=0; i<r_e-r_s; i++) {
903: IPPseudoOrthogonalize(d->ipV,d->size_cX,cX,d->nBcX,R[i],PETSC_NULL,&d->nR[r_s+i],&lindep);
904: }
905: } else {
906: for (i=0; i<r_e-r_s; i++) {
907: IPOrthogonalize(d->ipI,0,PETSC_NULL,d->size_cX,PETSC_NULL,cX,R[i],PETSC_NULL,&d->nR[r_s+i],&lindep);
908: }
909: }
910: if(lindep || (PetscAbs(d->nR[r_s+i]) < PETSC_MACHINE_EPSILON)) {
911: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %G!\n",r_s+i,d->nR[r_s+i]);
912: }
913: }
914: if (!cX || (cX && d->orthoV_type == EPS_ORTH_BOPT)) {
915: for(i=0; i<r_e-r_s; i++) {
916: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
917: }
918: for(i=0; i<r_e-r_s; i++) {
919: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
920: }
921: }
923: return(0);
924: }
928: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
929: the norm associated to the eigenpair, where i = r_s..r_e
930: R, vectors of Vec of size r_e-r_s,
931: auxV, PetscMax(r_e+cX_in_H, 2*(r_e-r_s)) vectors,
932: auxS, auxiliar vector of size (d->size_cX+r_e)^2+6(d->size_cX+r_e)+(r_e-r_s)*d->size_H
933: */
934: PetscErrorCode dvd_calcpairs_eig_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
935: {
936: PetscInt i,size_in,n,ld,ldc,k;
937: PetscErrorCode ierr;
938: Vec *Bx;
939: PetscScalar *cS,*cT,*pcX,*pX,*pX0;
940: DvdReduction r;
941: DvdReductionChunk
942: ops[2];
943: DvdMult_copy_func
944: sr[2];
945: #if !defined(PETSC_USE_COMPLEX)
946: PetscScalar b[8];
947: Vec X[4];
948: #endif
952: /* Quick return */
953: if (!d->cS) { return(0); }
955: size_in = (d->size_cX+r_e)*(d->cX_in_AV+r_e)*(d->cT?2:1);
956: /* Check consistency */
957: if (d->size_auxV < PetscMax(2*(r_e-r_s),d->cX_in_AV+r_e) || d->size_auxS < PetscMax(d->size_H*(r_e-r_s) /* pX0 */, 2*size_in /* SlepcAllReduceSum */ )) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
959: n = d->size_cX+r_e;
960: DSSetDimensions(d->conv_ps,n,PETSC_IGNORE,0,0);
961: DSGetLeadingDimension(d->conv_ps,&ldc);
962: DSGetArray(d->conv_ps,DS_MAT_A,&cS);
963: SlepcDenseCopyTriang(cS,0,ldc,d->cS,0,d->ldcS,d->size_cS,d->size_cS);
964: if (d->cT) {
965: DSGetArray(d->conv_ps,DS_MAT_B,&cT);
966: SlepcDenseCopyTriang(cT,0,ldc,d->cT,0,d->ldcT,d->size_cS,d->size_cS);
967: }
968: DSGetLeadingDimension(d->ps,&ld);
969: DSGetArray(d->ps,DS_MAT_Q,&pX);
970: /* Prepare reductions */
971: SlepcAllReduceSumBegin(ops,2,d->auxS,d->auxS+size_in,size_in,&r,((PetscObject)d->V[0])->comm);
972: /* auxV <- AV * pX(0:r_e+cX_in_H) */
973: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->AV-d->cX_in_AV,d->size_AV+d->cX_in_AV,pX,ld,d->size_H,d->cX_in_AV+r_e);
974: /* cS(:, size_cS:) <- cX' * auxV */
975: VecsMultS(&cS[ldc*d->size_cS],0,ldc,d->cY?d->cY:d->cX,0,d->size_cX+r_e,d->auxV,0,d->cX_in_AV+r_e,&r,&sr[0]);
977: if (d->cT) {
978: /* R <- BV * pX(0:r_e+cX_in_H) */
979: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->BV-d->cX_in_BV,d->size_BV+d->cX_in_BV,pX,ld,d->size_G,d->cX_in_BV+r_e);
980: /* cT(:, size_cS:) <- cX' * auxV */
981: VecsMultS(&cT[ldc*d->size_cT],0,ldc,d->cY?d->cY:d->cX,0,d->size_cY+r_e,d->auxV,0,d->cX_in_BV+r_e,&r,&sr[1]);
982: }
983: /* Do reductions */
984: SlepcAllReduceSumEnd(&r);
986: DSRestoreArray(d->conv_ps,DS_MAT_A,&cS);
987: if (d->cT) {
988: DSRestoreArray(d->conv_ps,DS_MAT_B,&cT);
989: }
990: DSSetState(d->conv_ps,DS_STATE_INTERMEDIATE);
991: /* eig(S,T) */
992: k = d->size_cX+r_s;
993: DSVectors(d->conv_ps,DS_MAT_X,&k,PETSC_NULL);
994: DSNormalize(d->conv_ps,DS_MAT_X,d->size_cX+r_s);
995: /* pX0 <- ps.Q(0:d->cX_in_AV+r_e-1) * conv_ps.X(size_cX-cX_in_H:) */
996: pX0 = d->auxS;
997: DSGetArray(d->conv_ps,DS_MAT_X,&pcX);
998: SlepcDenseMatProd(pX0,d->size_H,0.0,1.0,pX,ld,d->size_H,d->cX_in_AV+r_e,PETSC_FALSE,&pcX[(d->size_cX-d->cX_in_H)*ldc],ldc,n,r_e-r_s,PETSC_FALSE);
999: DSRestoreArray(d->ps,DS_MAT_Q,&pX);
1000: /* auxV <- cX(0:size_cX-cX_in_AV)*conv_ps.X + V*pX0 */
1001: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->cX,d->size_cX-d->cX_in_AV,pcX,ldc,n,r_e-r_s);
1002: DSRestoreArray(d->conv_ps,DS_MAT_X,&pcX);
1003: SlepcUpdateVectorsZ(d->auxV,d->size_cX-d->cX_in_AV==0?0.0:1.0,1.0,d->V-d->cX_in_AV,d->size_V+d->cX_in_AV,pX0,d->size_H,d->size_H,r_e-r_s);
1004: /* R <- A*auxV */
1005: for (i=0; i<r_e-r_s; i++) {
1006: MatMult(d->A,d->auxV[i],R[i]);
1007: }
1008: /* Bx <- B*auxV */
1009: if (d->B) {
1010: Bx = &d->auxV[r_e-r_s];
1011: for (i=0; i<r_e-r_s; i++) {
1012: MatMult(d->B,d->auxV[i],Bx[i]);
1013: }
1014: } else {
1015: Bx = d->auxV;
1016: }
1017: /* R <- (A - eig*B)*V*pX */
1018: for(i=0; i<r_e-r_s; i++) {
1019: #if !defined(PETSC_USE_COMPLEX)
1020: if (d->eigi[r_s+i] != 0.0) {
1021: /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [ 1 0
1022: 0 1
1023: -eigr_i -eigi_i
1024: eigi_i -eigr_i] */
1025: b[0] = b[5] = 1.0;
1026: b[2] = b[7] = -d->eigr[r_s+i];
1027: b[6] = -(b[3] = d->eigi[r_s+i]);
1028: b[1] = b[4] = 0.0;
1029: X[0] = R[i]; X[1] = R[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
1030: SlepcUpdateVectorsD(X,4,1.0,b,4,4,2,d->auxS,d->size_auxS);
1031: i++;
1032: } else
1033: #endif
1034: {
1035: /* R <- Ax -eig*Bx */
1036: VecAXPBY(R[i], -d->eigr[r_s+i], 1.0, Bx[i]);
1037: }
1038: }
1039: /* nR <- ||R|| */
1040: for(i=0; i<r_e-r_s; i++) {
1041: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
1042: }
1043: for(i=0; i<r_e-r_s; i++) {
1044: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
1045: }
1047: return(0);
1048: }
1051: /**** Pattern routines ********************************************************/
1053: /* BV <- BV*MT */
1056: 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,DSMatType mat)
1057: {
1058: PetscErrorCode ierr;
1059: PetscInt cMT,rm,cp,tra_s,i,ld;
1060: Vec *nBV;
1061: PetscScalar *MT;
1065: if (!real_BV || !*BV || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }
1067: DSGetLeadingDimension(d->ps,&ld);
1068: DSGetArray(d->ps,mat,&MT);
1069: if (d->V_tra_s > d->max_cX_in_proj && !BV_shift) {
1070: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
1071: cMT = d->V_tra_e - tra_s;
1072: rm = d->V_tra_s - tra_s;
1073: cp = PetscMin(d->max_cX_in_proj - rm, *cX_in_proj);
1074: nBV = real_BV+d->max_cX_in_proj;
1075: /* BV(-cp-rm:-1-rm) <- BV(-cp:-1) */
1076: for (i=-cp; i<0; i++) {
1077: VecCopy((*BV)[i], nBV[i-rm]);
1078: }
1079: /* BV(-rm:) <- BV*MT(tra_s:V_tra_e-1) */
1080: SlepcUpdateVectorsZ(&nBV[-rm],0.0,1.0,*BV-*cX_in_proj,*size_BV+*cX_in_proj,&MT[ld*tra_s],ld,d->size_MT,cMT);
1081: *size_BV = d->V_tra_e - d->V_tra_s;
1082: *max_size_BV-= nBV - *BV;
1083: *BV = nBV;
1084: if (cX_in_proj && d->max_cX_in_proj>0) *cX_in_proj = cp+rm;
1085: } else if (d->V_tra_s <= d->max_cX_in_proj || BV_shift) {
1086: /* [BcX BV] <- [BcX BV*MT] */
1087: SlepcUpdateVectorsZ(*BV-*cX_in_proj,0.0,1.0,*BV-*cX_in_proj,*size_BV+*cX_in_proj,MT,ld,d->size_MT,d->V_tra_e);
1088: *BV+= d->V_tra_s-*cX_in_proj;
1089: *max_size_BV-= d->V_tra_s-*cX_in_proj;
1090: *size_BV = d->V_tra_e - d->V_tra_s;
1091: if (size_cBV && BV_shift) *size_cBV = *BV - real_BV;
1092: if (d->max_cX_in_proj>0) *cX_in_proj = PetscMin(*BV - real_BV, d->max_cX_in_proj);
1093: } else { /* !BV_shift */
1094: /* BV <- BV*MT(V_tra_s:) */
1095: SlepcUpdateVectorsZ(*BV,0.0,1.0,*BV,*size_BV,&MT[d->V_tra_s*ld],ld,d->size_MT,d->V_tra_e-d->V_tra_s);
1096: *size_BV = d->V_tra_e - d->V_tra_s;
1097: }
1098: DSRestoreArray(d->ps,mat,&MT);
1099: return(0);
1100: }