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_eig_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R);
50: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
51: PetscInt r_e, Vec *R);
52: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
53: DvdMult_copy_func **sr);
54: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d);
55: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
56: DvdMult_copy_func **sr);
57: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d);
58: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d);
59: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
60: DvdMult_copy_func **sr);
61: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d);
62: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
63: DvdMult_copy_func **sr);
64: 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);
66: /**** Control routines ********************************************************/
69: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d, dvdBlackboard *b,
70: orthoV_type_t orth, IP ipI,
71: PetscInt cX_proj, PetscBool harm)
72: {
73: PetscErrorCode ierr;
74: PetscInt i,max_cS;
75: PetscBool std_probl,her_probl;
79: std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
80: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
82: /* Setting configuration constrains */
83: #if !defined(PETSC_USE_COMPLEX)
84: /* if the last converged eigenvalue is complex its conjugate pair is also
85: converged */
86: b->max_nev = PetscMax(b->max_nev, d->nev+(her_probl && !d->B?0:1));
87: #else
88: b->max_nev = PetscMax(b->max_nev, d->nev);
89: #endif
90: b->max_size_proj = PetscMax(b->max_size_proj, b->max_size_V+cX_proj);
91: d->size_real_V = b->max_size_V+b->max_nev;
92: d->W_shift = d->B?PETSC_TRUE:PETSC_FALSE;
93: d->size_real_W = harm?(b->max_size_V+(d->W_shift?b->max_nev:b->max_size_cP)):0;
94: d->size_real_AV = b->max_size_V+b->max_size_cP;
95: d->size_BDS = 0;
96: if (d->B && her_probl && (orth == DVD_ORTHOV_I || orth == DVD_ORTHOV_BOneMV)) {
97: d->size_real_BV = b->size_V; d->BV_shift = PETSC_TRUE;
98: if (orth == DVD_ORTHOV_BOneMV) d->size_BDS = d->eps->nds;
99: } else if (d->B) {
100: d->size_real_BV = b->max_size_V + b->max_size_P; d->BV_shift = PETSC_FALSE;
101: } else {
102: d->size_real_BV = 0; d->BV_shift = PETSC_FALSE;
103: }
104: b->own_vecs+= d->size_real_V + d->size_real_W + d->size_real_AV +
105: d->size_real_BV + d->size_BDS;
106: b->own_scalars+= b->max_size_proj*b->max_size_proj*2*(std_probl?1:2) +
107: /* H, G?, S, T? */
108: b->max_size_proj*b->max_size_proj*(std_probl?1:2) +
109: /* pX, pY? */
110: b->max_nev*b->max_nev*(her_probl?0:(!d->B?1:2));
111: /* cS?, cT? */
112: b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X);
113: /* updateV0 */
114: max_cS = PetscMax(b->max_size_X,cX_proj);
115: b->max_size_auxS = PetscMax(PetscMax(
116: b->max_size_auxS,
117: b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
118: max_cS*b->max_nev*(her_probl?0:(!d->B?1:2)) + /* updateV0,W0 */
119: /* SlepcReduction: in */
120: PetscMax(
121: b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
122: max_cS*b->max_nev*(her_probl?0:(!d->B?1:2)), /* updateV0,W0 */
123: /* SlepcReduction: out */
124: PetscMax(
125: b->max_size_proj*b->max_size_proj, /* updateAV0,BV0 */
126: b->max_size_proj+b->max_nev))), /* dvd_orth */
127: std_probl?0:(b->max_size_proj*11+16) /* projeig */);
128: #if defined(PETSC_USE_COMPLEX)
129: b->max_size_auxS = PetscMax(b->max_size_auxS, b->max_size_V);
130: /* dvd_calcpairs_projeig_eig */
131: #endif
133: /* Setup the step */
134: if (b->state >= DVD_STATE_CONF) {
135: d->max_cX_in_proj = cX_proj;
136: d->max_size_P = b->max_size_P;
137: d->real_V = b->free_vecs; b->free_vecs+= d->size_real_V;
138: if (harm) {
139: d->real_W = b->free_vecs; b->free_vecs+= d->size_real_W;
140: } else {
141: d->real_W = PETSC_NULL;
142: }
143: d->real_AV = d->AV = b->free_vecs; b->free_vecs+= d->size_real_AV;
144: d->max_size_proj = b->max_size_proj;
145: d->real_H = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
146: d->ldH = b->max_size_proj;
147: d->pX = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
148: d->S = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
149: if (!her_probl) {
150: d->cS = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
151: d->max_size_cS = d->ldcS = b->max_nev;
152: } else {
153: d->cS = PETSC_NULL;
154: d->max_size_cS = d->ldcS = 0;
155: }
156: d->ipV = ipI;
157: d->ipW = ipI;
158: if (orth == DVD_ORTHOV_BOneMV) {
159: d->BDS = b->free_vecs; b->free_vecs+= d->eps->nds;
160: for (i=0; i<d->eps->nds; i++) {
161: MatMult(d->B, d->eps->DS[i], d->BDS[i]);
162: }
163: } else
164: d->BDS = PETSC_NULL;
165: if (d->B) {
166: d->real_BV = b->free_vecs; b->free_vecs+= d->size_real_BV;
167: } else {
168: d->size_real_BV = 0;
169: d->real_BV = PETSC_NULL;
170: d->BV_shift = PETSC_FALSE;
171: }
172: if (!std_probl) {
173: d->real_G = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
174: d->T = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
175: d->pY = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
176: } else {
177: d->real_G = PETSC_NULL;
178: d->T = PETSC_NULL;
179: d->pY = PETSC_NULL;
180: }
181: if (d->B && !her_probl) {
182: d->cT = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
183: d->ldcT = b->max_nev;
184: } else {
185: d->cT = PETSC_NULL;
186: d->ldcT = 0;
187: }
189: d->calcPairs = dvd_calcpairs_proj;
190: d->calcpairs_residual = dvd_calcpairs_res_0;
191: d->calcpairs_residual_eig = dvd_calcpairs_eig_res_0;
192: d->calcpairs_proj_res = dvd_calcpairs_proj_res;
193: d->calcpairs_selectPairs = PETSC_NULL;
194: d->ipI = ipI;
195: DVD_FL_ADD(d->startList, dvd_calcpairs_qz_start);
196: }
198: return(0);
199: }
204: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
205: {
206: PetscBool her_probl;
207: PetscInt i;
210: her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;
212: d->size_V = 0;
213: d->V = d->real_V;
214: d->cX = d->real_V;
215: d->size_cX = 0;
216: d->max_size_V = d->size_real_V;
217: d->W = d->real_W;
218: d->max_size_W = d->size_real_W;
219: d->size_W = 0;
220: d->size_AV = 0;
221: d->AV = d->real_AV;
222: d->max_size_AV = d->size_real_AV;
223: d->size_H = 0;
224: d->H = d->real_H;
225: if (d->cS) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cS[i] = 0.0;
226: d->size_BV = 0;
227: d->BV = d->real_BV;
228: d->max_size_BV = d->size_real_BV;
229: d->size_G = 0;
230: d->G = d->real_G;
231: if (d->cT) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cT[i] = 0.0;
232: d->cY = d->B && !her_probl ? d->W : PETSC_NULL;
233: d->BcX = d->orthoV_type == DVD_ORTHOV_I && d->B && her_probl ? d->BcX : PETSC_NULL;
234: d->size_cY = 0;
235: d->size_BcX = 0;
236: 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: #define 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->cX_in_H)*d->V_tra_s*(d->cT?2:(d->cS?1:0)) + /* updateV0,W0 */
261: (d->size_H*(d->V_new_e-d->V_new_s)*2+
262: (d->V_new_e-d->V_new_s)*(d->V_new_e-d->V_new_s))*(!stdp?2:1); /* updateAV1,BV1 */
263:
264: out = in+size_in;
266: /* Check consistency */
267: if (2*size_in > d->size_auxS) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
269: /* Prepare reductions */
270: SlepcAllReduceSumBegin(ops, MAX_OPS, in, out, size_in, &r,
271: ((PetscObject)d->V[0])->comm);
272: /* Allocate size_in */
273: d->auxS+= size_in;
274: d->size_auxS-= size_in;
276: /* Update AV, BV, W and the projected matrices */
277: /* 1. S <- S*MT */
278: dvd_calcpairs_updateV0(d, &r, &sr0);
279: dvd_calcpairs_updateW0(d, &r, &sr0);
280: dvd_calcpairs_updateAV0(d);
281: dvd_calcpairs_updateBV0(d);
282: /* 2. V <- orth(V, V_new) */
283: dvd_calcpairs_updateV1(d);
284: /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
285: /* Check consistency */
286: if (d->size_AV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
287: for (i=d->V_new_s; i<d->V_new_e; i++) {
288: MatMult(d->A, d->V[i], d->AV[i]);
289: }
290: d->size_AV = d->V_new_e;
291: /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
292: if (d->B && d->orthoV_type != DVD_ORTHOV_BOneMV) {
293: /* Check consistency */
294: if (d->size_BV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
295: for (i=d->V_new_s; i<d->V_new_e; i++) {
296: MatMult(d->B, d->V[i], d->BV[i]);
297: }
298: d->size_BV = d->V_new_e;
299: }
300: /* 5 <- W <- [W f(AV,BV)] */
301: dvd_calcpairs_updateW1(d);
302: dvd_calcpairs_updateAV1(d, &r, &sr0);
303: dvd_calcpairs_updateBV1(d, &r, &sr0);
305: /* Deallocate size_in */
306: d->auxS-= size_in;
307: d->size_auxS+= size_in;
309: /* Do reductions */
310: SlepcAllReduceSumEnd(&r);
312: /* Perform the transformation on the projected problem */
313: if (d->calcpairs_proj_trans) {
314: d->calcpairs_proj_trans(d);
315: }
317: d->V_tra_s = d->V_tra_e = 0;
318: d->V_new_s = d->V_new_e;
320: /* Solve the projected problem */
321: if (d->size_H>0) {
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: }
331: }
333: /* Check consistency */
334: 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 ||
335: d->size_V != d->size_AV || d->cX_in_H != d->cX_in_AV ||
336: (DVD_ISNOT(d->sEP, DVD_EP_STD) && (
337: d->size_V+d->cX_in_G != d->size_G || d->cX_in_H != d->cX_in_G ||
338: d->size_H != d->size_G || (d->BV && (
339: d->size_V != d->size_BV || d->cX_in_H != d->cX_in_BV)))) ||
340: (d->W && d->size_W != d->size_V)) {
341: SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!");
342: }
344: return(0);
345: #undef MAX_OPS
346: }
348: /**** Basic routines **********************************************************/
352: /* auxV: V_tra_s, DvdMult_copy_func: 1 */
353: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
354: DvdMult_copy_func **sr)
355: {
356: PetscErrorCode ierr;
357: PetscInt rm;
361: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
363: /* cX <- [cX V*MT(0:V_tra_s-1)], V <- V*MT(V_tra_s:V_tra_e) */
364: 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);
366: /* Udpate cS for standard problems */
367: if (d->cS && !d->cT && !d->cY && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
368: /* Check consistency */
369: if (d->size_cS+d->V_tra_s != d->size_cX) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
371: /* auxV <- AV * MTX(0:V_tra_e-1) */
372: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
373: 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-rm);
375: /* cS(:, size_cS:) <- cX' * auxV */
376: 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)++);
377: d->size_cS+= d->V_tra_s-rm;
378: }
379:
380: return(0);
381: }
386: /* auxS: size_cX+V_new_e+1 */
387: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d)
388: {
389: PetscErrorCode ierr;
390: Vec *cX = d->BcX? d->BcX : ( (d->cY && !d->W)? d->cY : d->cX );
394: if (d->V_new_s == d->V_new_e) { return(0); }
396: /* Check consistency */
397: if (d->size_V != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
399: /* V <- gs([cX V(0:V_new_s-1)], V(V_new_s:V_new_e-1)) */
400: if (d->orthoV_type == DVD_ORTHOV_BOneMV) {
401: dvd_BorthV(d->ipV, d->eps->DS, d->BDS, d->eps->nds, d->cX, d->real_BV,
402: d->size_cX, d->V, d->BV, d->V_new_s, d->V_new_e,
403: d->auxS, d->eps->rand);
404: d->size_BV = d->V_new_e;
405: } else {
406: dvd_orthV(d->ipV, d->eps->DS, d->eps->nds, cX, d->size_cX, d->V,
407: d->V_new_s, d->V_new_e, d->auxS, d->eps->rand);
408:
409: }
410: d->size_V = d->V_new_e;
412: return(0);
413: }
417: /* auxV: V_tra_s, DvdMult_copy_func: 2 */
418: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
419: DvdMult_copy_func **sr)
420: {
421: PetscErrorCode ierr;
422: PetscInt rm;
426: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
428: /* cY <- [cY W*MTY(0:V_tra_s-1)], W <- W*MTY(V_tra_s:V_tra_e) */
429: 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);
431: /* Udpate cS and cT */
432: if (d->cT && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
433: /* Check consistency */
434: 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!"); }
436: /* auxV <- AV * MTX(0:V_tra_e-1) */
437: rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
438: 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-rm);
440: /* cS(:, size_cS:) <- cY' * auxV */
441: 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)++);
443: /* auxV <- BV * MTX(0:V_tra_e-1) */
444: 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-rm);
446: /* cT(:, size_cS:) <- cY' * auxV */
447: 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)++);
448:
449: d->size_cS+= d->V_tra_s-rm;
450: d->size_cT+= d->V_tra_s-rm;
451: }
453: return(0);
454: }
459: /* auxS: size_cX+V_new_e+1 */
460: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d)
461: {
462: PetscErrorCode ierr;
463: Vec *cY = d->cY?d->cY:d->cX;
467: if (!d->W || d->V_new_s == d->V_new_e) { return(0); }
469: /* Check consistency */
470: if (d->size_W != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
472: /* Update W */
473: d->calcpairs_W(d);
475: /* W <- gs([cY W(0:V_new_s-1)], W(V_new_s:V_new_e-1)) */
476: 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);
477: d->size_W = d->V_new_e;
479: return(0);
480: }
484: /* auxS: size_H*(V_tra_e-V_tra_s) */
485: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d)
486: {
487: PetscErrorCode ierr;
488: PetscScalar *MTY = d->W?d->MTY:d->MTX;
489: PetscInt cMT, tra_s;
493: if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }
495: /* AV(V_tra_s-cp-1:) = cAV*MTX(V_tra_s:) */
496: 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);
497: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
498: cMT = d->V_tra_e - tra_s;
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 = d->cX_in_AV;
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, 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: 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);
572: 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);
573: d->size_G = cMT;
574: d->cX_in_G = d->cX_in_V;
575: }
577: return(0);
578: }
583: /* DvdMult_copy_func: 2 */
584: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
585: DvdMult_copy_func **sr)
586: {
587: PetscErrorCode ierr;
588: Vec *W = d->W?d->W:d->V, *BV = d->BV?d->BV:d->V;
592: if (!d->G || d->V_new_s == d->V_new_e) { return(0); }
594: /* G = [G W(old)'*BV(new);
595: W(new)'*BV(old) W(new)'*BV(new) ],
596: being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
597: 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)++);
598: d->size_G = d->V_new_e+d->cX_in_G;
600: return(0);
601: }
603: /* in complex, d->size_H real auxiliar values are needed */
606: PetscErrorCode dvd_calcpairs_projeig_eig(dvdDashboard *d)
607: {
608: PetscErrorCode ierr;
609: PetscReal *w;
610: #if defined(PETSC_USE_COMPLEX)
611: PetscInt i;
612: #endif
616: /* S <- H */
617: d->ldS = d->ldpX = d->size_H;
618: SlepcDenseCopyTriang(d->S, DVD_MAT_LTRIANG, d->size_H, d->H, d->sH, d->ldH,
619: d->size_H, d->size_H);
621: /* S = pX' * L * pX */
622: #if !defined(PETSC_USE_COMPLEX)
623: w = d->eigr-d->cX_in_H;
624: #else
625: w = (PetscReal*)d->auxS;
626: #endif
627: EPSDenseHEP(d->size_H, d->S, d->ldS, w, d->pX);
628: #if defined(PETSC_USE_COMPLEX)
629: for (i=0; i<d->size_H; i++) d->eigr[i-d->cX_in_H] = w[i];
630: #endif
632: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_eig;
634: return(0);
635: }
640: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d)
641: {
642: PetscErrorCode ierr;
646: /* S <- H */
647: d->ldS = d->ldpX = d->size_H;
648: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
649: d->size_H, d->size_H);
651: /* S = pX' * H * pX */
652: EPSDenseHessenberg(d->size_H, 0, d->S, d->ldS, d->pX);
653: EPSDenseSchur(d->size_H, 0, d->S, d->ldS, d->pX, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
654:
656: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;
658: return(0);
659: }
663: /*
664: auxS(dgges) = size_H (beta) + 8*size_H+16 (work)
665: auxS(zgges) = size_H (beta) + 1+2*size_H (work) + 8*size_H (rwork)
666: */
667: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d)
668: {
669: #if defined(SLEPC_MISSING_LAPACK_GGES)
671: SETERRQ(((PetscObject)(d->eps))->comm,PETSC_ERR_SUP,"GGES - Lapack routine is unavailable.");
672: #else
673: PetscErrorCode ierr;
674: PetscScalar *beta = d->auxS;
675: #if !defined(PETSC_USE_COMPLEX)
676: PetscScalar *auxS = beta + d->size_H;
677: PetscBLASInt n_auxS = d->size_auxS - d->size_H;
678: #else
679: PetscReal *auxR = (PetscReal*)(beta + d->size_H);
680: PetscScalar *auxS = (PetscScalar*)(auxR+8*d->size_H);
681: PetscBLASInt n_auxS = d->size_auxS - 9*d->size_H;
682: #endif
683: PetscInt i;
684: PetscBLASInt info,n,a;
694: /* S <- H, T <- G */
695: d->ldS = d->ldT = d->ldpX = d->ldpY = d->size_H;
696: SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
697: d->size_H, d->size_H);
698: SlepcDenseCopyTriang(d->T, 0, d->size_H, d->G, d->sG, d->ldH,
699: d->size_H, d->size_H);
701: /* S = Z'*H*Q, T = Z'*G*Q */
702: n = d->size_H;
703: #if !defined(PETSC_USE_COMPLEX)
704: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
705: &a, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
706: auxS, &n_auxS, PETSC_NULL, &info);
707: #else
708: LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
709: &a, d->eigr-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
710: auxS, &n_auxS, auxR, PETSC_NULL, &info);
711: #endif
712: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack GGES %d", info);
714: /* eigr[i] <- eigr[i] / beta[i] */
715: for (i=0; i<d->size_H; i++)
716: d->eigr[i-d->cX_in_H] /= beta[i],
717: d->eigi[i-d->cX_in_H] /= beta[i];
719: d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;
721: return(0);
722: #endif
723: }
727: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n)
728: {
729: PetscErrorCode ierr;
733: EPSSortDenseHEP(d->eps, d->size_H, 0, d->eigr-d->cX_in_H, d->pX, d->ldpX);
734:
736: if (d->calcpairs_eigs_trans) {
737: d->calcpairs_eigs_trans(d);
738: }
740: return(0);
741: }
746: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n)
747: {
748: PetscErrorCode ierr;
751: if ((d->ldpX != d->size_H) ||
752: ( d->T &&
753: ((d->ldS != d->ldT) || (d->ldpX != d->ldpY) ||
754: (d->ldpX != d->size_H)) ) ) {
755: SETERRQ(PETSC_COMM_SELF,1, "Error before ordering eigenpairs");
756: }
758: if (d->T) {
759: EPSSortDenseSchurGeneralized(d->eps, d->size_H, 0, n, d->S, d->T,
760: d->ldS, d->pY, d->pX, d->eigr-d->cX_in_H,
761: d->eigi-d->cX_in_H);
762: } else {
763: EPSSortDenseSchur(d->eps, d->size_H, 0, d->S, d->ldS, d->pX,
764: d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
765: }
767: if (d->calcpairs_eigs_trans) {
768: d->calcpairs_eigs_trans(d);
769: }
771: #if defined(PETSC_USE_COMPLEX)
772: if (d->T) {
773: EPSCleanDenseSchur(d->size_H,0,d->S,d->ldS,d->T,d->ldT,d->eigi-d->cX_in_H,d->pX,d->ldpX,PETSC_TRUE);
774: }
775: #endif
777: return(0);
778: }
783: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
784: the norm associated to the Schur pair, where i = r_s..r_e
785: */
786: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
787: Vec *R)
788: {
789: PetscInt i;
790: PetscErrorCode ierr;
791: Vec *BV = d->BV?d->BV:d->V;
795: for (i=r_s; i<r_e; i++) {
796: /* nX(i) <- ||X(i)|| */
797: if (d->correctXnorm) {
798: /* R(i) <- V*pX(i) */
799: SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
800: &d->V[-d->cX_in_H], d->size_V+d->cX_in_H,
801: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
802: VecNorm(R[i-r_s], NORM_2, &d->nX[i]);
803: } else
804: d->nX[i] = 1.0;
806: /* R(i-r_s) <- AV*pX(i) */
807: SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
808: &d->AV[-d->cX_in_H], d->size_AV+d->cX_in_H,
809: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
811: /* R(i-r_s) <- R(i-r_s) - eigr(i)*BV*pX(i) */
812: SlepcUpdateVectorsZ(&R[i-r_s], 1.0, -d->eigr[i+d->cX_in_H],
813: &BV[-d->cX_in_H], d->size_V+d->cX_in_H,
814: &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
815: }
817: d->calcpairs_proj_res(d, r_s, r_e, R);
819: return(0);
820: }
824: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
825: PetscInt r_e, Vec *R)
826: {
827: PetscInt i;
828: PetscErrorCode ierr;
829: PetscBool lindep;
830: Vec *cX;
834: /* If exists the BcX, R <- orth(BcX, R), nR[i] <- ||R[i]|| */
835: if (d->BcX)
836: cX = d->BcX;
838: /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
839: else if (d->cY)
840: cX = d->cY;
842: /* If fany configurations, R <- orth(cX, R), nR[i] <- ||R[i]|| */
843: else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN)))
844: cX = d->cX;
846: /* Otherwise, nR[i] <- ||R[i]|| */
847: else
848: cX = PETSC_NULL;
850: if (cX) for (i=0; i<r_e-r_s; i++) {
851: IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX, PETSC_NULL,
852: cX, R[i], PETSC_NULL, &d->nR[r_s+i], &lindep);
853:
854: if(lindep || (d->nR[r_s+i] < PETSC_MACHINE_EPSILON)) {
855: PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %G!\n",r_s+i,d->nR[r_s+i]);
856: }
858: } else {
859: for(i=0; i<r_e-r_s; i++) {
860: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
861: }
862: for(i=0; i<r_e-r_s; i++) {
863: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
864: }
865: }
867: return(0);
868: }
872: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
873: the norm associated to the eigenpair, where i = r_s..r_e
874: R, vectors of Vec of size r_e-r_s,
875: auxV, PetscMax(r_e+cX_in_H, 2*(r_e-r_s)) vectors,
876: auxS, auxiliar vector of size (d->size_cX+r_e)^2+6(d->size_cX+r_e)+(r_e-r_s)*d->size_H
877: */
878: PetscErrorCode dvd_calcpairs_eig_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
879: {
880: PetscInt i,ldpX,size_in;
881: PetscErrorCode ierr;
882: Vec *Bx;
883: PetscScalar *pX,*pX0;
884: DvdReduction r;
885: DvdReductionChunk
886: ops[2];
887: DvdMult_copy_func
888: sr[2];
889: #if !defined(PETSC_USE_COMPLEX)
890: PetscScalar b[8];
891: Vec X[4];
892: #endif
896: /* Quick return */
897: if (!d->cS) { return(0); }
899: size_in = (d->size_cX+r_e)*(d->cX_in_AV+r_e)*(d->cT?2:1);
900: /* Check consistency */
901: if (d->size_auxV < PetscMax(2*(r_e-r_s),d->cX_in_AV+r_e) || d->size_auxS <
902: (d->size_cX+r_e)*(d->size_cX+r_e) /* pX */ + PetscMax(PetscMax(
903: (d->size_cX+r_e)*6 /* dvd_compute_eigenvectors */,
904: d->size_H*(r_e-r_s) /* pX0 */),
905: 2*size_in /* SlepcAllReduceSum */ )) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
907: /* Prepare reductions */
908: SlepcAllReduceSumBegin(ops,2,d->auxS,d->auxS+size_in,size_in,&r,((PetscObject)d->V[0])->comm);
910: /* auxV <- AV * pX(0:r_e+cX_in_H) */
911: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->AV-d->cX_in_AV,d->size_AV+d->cX_in_AV,d->pX,d->ldpX,d->size_H,d->cX_in_AV+r_e);
913: /* cS(:, size_cS:) <- cX' * auxV */
914: VecsMultS(&d->cS[d->ldcS*d->size_cS],0,d->ldcS,d->cY?d->cY:d->cX,0,d->size_cX+r_e,d->auxV,0,d->cX_in_AV+r_e,&r,&sr[0]);
916: if (d->cT) {
917: /* R <- BV * pX(0:r_e+cX_in_H) */
918: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->BV-d->cX_in_BV,d->size_BV+d->cX_in_BV,d->pX,d->ldpX,d->size_G,d->cX_in_BV+r_e);
919:
920: /* cS(:, size_cS:) <- cX' * auxV */
921: VecsMultS(&d->cT[d->ldcT*d->size_cT],0,d->ldcT,d->cY?d->cY:d->cX,0,d->size_cY+r_e,d->auxV,0,d->cX_in_BV+r_e,&r,&sr[1]);
922: }
924: /* Do reductions */
925: SlepcAllReduceSumEnd(&r);
927: /* qX <- eig(S,T) */
928: pX = d->auxS; d->auxS+= (d->size_cX+r_e)*(d->size_cX+r_e); d->size_auxS -= (d->size_cX+r_e)*(d->size_cX+r_e);
929: ldpX = d->size_cX+r_e;
930: pX0 = d->auxS;
931: EPSCleanDenseSchur(d->size_cX+r_e,0,d->cS,d->ldcS,d->cT,d->ldcT,d->ceigi,pX,ldpX,PETSC_FALSE);
932: dvd_compute_eigenvectors(d->size_cX+r_e,d->cS,d->ldcS,d->cT,d->ldcT,pX,ldpX,PETSC_NULL,0,d->auxS,d->size_auxS,PETSC_TRUE);
934: /* pX[i] <- pX[i] / ||pX[i]|| */
935: pX = &pX[ldpX*d->size_cX+r_s];
936: SlepcDenseNorm(pX,ldpX,d->size_cX+r_e,r_e-r_s,d->eigi+r_s);
938: /* pX0 <- d->pX(0:d->cX_in_AV+r_e-1) * pX(size_cX-cX_in_H:) */
939: SlepcDenseMatProd(pX0,d->size_H,0.0,1.0,d->pX,d->ldpX,d->size_H,d->cX_in_AV+r_e,PETSC_FALSE,&pX[d->size_cX-d->cX_in_H],ldpX,d->cX_in_AV+r_e,r_e-r_s,PETSC_FALSE);
941: /* auxV <- cX(0:size_cX-cX_in_AV)*pX + V*pX0 */
942: SlepcUpdateVectorsZ(d->auxV,0.0,1.0,d->cX,d->size_cX-d->cX_in_AV,pX,ldpX,d->size_cX-d->cX_in_AV,r_e-r_s);
943: 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);
944: d->auxS-= (d->size_cX+r_e)*(d->size_cX+r_e); d->size_auxS += (d->size_cX+r_e)*(d->size_cX+r_e);
945: /* R <- A*auxV */
946: for (i=0; i<r_e-r_s; i++) {
947: MatMult(d->A,d->auxV[i],R[i]);
948: }
949: /* Bx <- B*auxV */
950: if (d->B) {
951: Bx = &d->auxV[r_e-r_s];
952: for (i=0; i<r_e-r_s; i++) {
953: MatMult(d->B,d->auxV[i],Bx[i]);
954: }
955: } else {
956: Bx = d->auxV;
957: }
958: /* R <- (A - eig*B)*V*pX */
959: for(i=0; i<r_e-r_s; i++) {
960: #if !defined(PETSC_USE_COMPLEX)
961: if (d->eigi[r_s+i] != 0.0) {
962: /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [ 1 0
963: 0 1
964: -eigr_i -eigi_i
965: eigi_i -eigr_i] */
966: b[0] = b[5] = 1.0;
967: b[2] = b[7] = -d->eigr[r_s+i];
968: b[6] = -(b[3] = d->eigi[r_s+i]);
969: b[1] = b[4] = 0.0;
970: X[0] = R[i]; X[1] = R[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
971: SlepcUpdateVectorsD(X,4,1.0,b,4,4,2,d->auxS,d->size_auxS);
972: i++;
973: } else
974: #endif
975: {
976: /* R <- Ax -eig*Bx */
977: VecAXPBY(R[i], -d->eigr[r_s+i], 1.0, Bx[i]);
978: }
979: }
981: /* nR <- ||R|| */
982: for(i=0; i<r_e-r_s; i++) {
983: VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
984: }
985: for(i=0; i<r_e-r_s; i++) {
986: VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
987: }
989: return(0);
990: }
993: /**** Pattern routines ********************************************************/
995: /* BV <- BV*MTX */
998: 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)
999: {
1000: PetscErrorCode ierr;
1001: PetscInt cMT, rm, cp, tra_s, i;
1002: Vec *nBV;
1006: if (!real_BV || !*BV || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }
1008: if (d->V_tra_s > d->max_cX_in_proj && !BV_shift) {
1009: tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
1010: cMT = d->V_tra_e - tra_s;
1011: rm = d->V_tra_s - tra_s;
1012: cp = PetscMin(d->max_cX_in_proj - rm, *cX_in_proj);
1013: nBV = real_BV+d->max_cX_in_proj;
1014: /* BV(-cp-rm:-1-rm) <- BV(-cp:-1) */
1015: for (i=-cp; i<0; i++) {
1016: VecCopy((*BV)[i], nBV[i-rm]);
1017: }
1018: /* BV(-rm:) <- BV*MTX(tra_s:V_tra_e-1) */
1019: 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);
1020: *size_BV = d->V_tra_e - d->V_tra_s;
1021: *max_size_BV-= nBV - *BV;
1022: *BV = nBV;
1023: if (cX_in_proj && d->max_cX_in_proj>0) *cX_in_proj = cp+rm;
1024: } else if (d->V_tra_s <= d->max_cX_in_proj || BV_shift) {
1025: /* [BcX BV] <- [BcX BV*MTX] */
1026: 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);
1027: *BV+= d->V_tra_s-*cX_in_proj;
1028: *max_size_BV-= d->V_tra_s-*cX_in_proj;
1029: *size_BV = d->V_tra_e - d->V_tra_s;
1030: if (size_cBV && BV_shift) *size_cBV = *BV - real_BV;
1031: if (d->max_cX_in_proj>0) *cX_in_proj = PetscMin(*BV - real_BV, d->max_cX_in_proj);
1032: } else { /* !BV_shift */
1033: /* BV <- BV*MTX(V_tra_s:) */
1034: SlepcUpdateVectorsZ(*BV, 0.0, 1.0, *BV, *size_BV,
1035: &MTX[d->V_tra_s*d->ldMTX], d->ldMTX, d->size_MT, d->V_tra_e-d->V_tra_s);
1036:
1037: *size_BV = d->V_tra_e - d->V_tra_s;
1038: }
1040: return(0);
1041: }