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 (DVD_IS(d->sEP, DVD_EP_STD)) {
322:     if (DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
323:       dvd_calcpairs_projeig_eig(d);
324:     } else {
325:       dvd_calcpairs_projeig_qz_std(d);
326:     }
327:   } else {
328:     dvd_calcpairs_projeig_qz_gen(d);
329:   }

331:   /* Check consistency */
332:   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 ||
333:       d->size_V != d->size_AV || d->cX_in_H != d->cX_in_AV ||
334:         (DVD_ISNOT(d->sEP, DVD_EP_STD) && (
335:           d->size_V+d->cX_in_G != d->size_G || d->cX_in_H != d->cX_in_G ||
336:           d->size_H != d->size_G || (d->BV && (
337:             d->size_V != d->size_BV || d->cX_in_H != d->cX_in_BV)))) ||
338:       (d->W && d->size_W != d->size_V)) {
339:     SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!");
340:   }

342:   return(0);
343: #undef MAX_OPS
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:     rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
371:     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);

373:     /* cS(:, size_cS:) <- cX' * auxV */
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->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->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_cX || (d->W && d->size_cY != d->size_cX)) { 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-rm);

438:     /* cS(:, size_cS:) <- cY' * auxV */
439:     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)++);

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-rm);

444:     /* cT(:, size_cS:) <- cY' * auxV */
445:     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)++);
446: 
447:     d->size_cS+= d->V_tra_s-rm;
448:     d->size_cT+= d->V_tra_s-rm;
449:   }

451:   return(0);
452: }


457: /* auxS: size_cX+V_new_e+1 */
458: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d)
459: {
460:   PetscErrorCode  ierr;
461:   Vec             *cY = d->cY?d->cY:d->cX;


465:   if (!d->W || d->V_new_s == d->V_new_e) { return(0); }

467:   /* Check consistency */
468:   if (d->size_W != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

470:   /* Update W */
471:   d->calcpairs_W(d);

473:   /* W <- gs([cY W(0:V_new_s-1)], W(V_new_s:V_new_e-1)) */
474:   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);
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;


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;

498:   /* Update H <- MTY(tra_s)' * (H * MTX(tra_s:)) */
499:   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);
500:   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);
501:   d->size_H = cMT;
502:   d->cX_in_H = d->cX_in_AV;

504:   return(0);
505: }


510: /* DvdMult_copy_func: 2 */
511: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
512:                                        DvdMult_copy_func **sr)
513: {
514:   PetscErrorCode  ierr;
515:   Vec             *W = d->W?d->W:d->V;


519:   if (d->V_new_s == d->V_new_e) { return(0); }

521:   /* Check consistency */
522:   if (d->size_H != d->V_new_s+d->cX_in_H || d->size_V != d->V_new_e) {
523:     SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!");
524:   }

526:   /* H = [H               W(old)'*AV(new);
527:           W(new)'*AV(old) W(new)'*AV(new) ],
528:      being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
529:   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)++);
530:   d->size_H = d->V_new_e+d->cX_in_H;

532:   return(0);
533: }

537: /* auxS: max(BcX*(size_cX+V_new_e+1), size_G*(V_tra_e-V_tra_s)) */
538: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d)
539: {
540:   PetscErrorCode  ierr;
541:   PetscScalar     *MTY = d->W?d->MTY:d->MTX;
542:   PetscInt        cMT, tra_s, i;
543:   PetscBool       lindep;
544:   PetscReal       norm;


548:   if (d->V_tra_s == 0 && d->V_tra_e == 0) { return(0); }

550:   /* BV <- BV*MTX */
551:   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);

553:   /* If BcX, BcX <- orth(BcX) */
554:   if (d->BcX) {
555:     for (i=0; i<d->V_tra_s; i++) {
556:       IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_BcX+i, PETSC_NULL,
557:                              d->BcX, d->BcX[d->size_BcX+i], PETSC_NULL,
558:                              &norm, &lindep);
559:       if(lindep) SETERRQ(((PetscObject)d->ipI)->comm,1, "Error during orth(BcX, B*cX(new))");
560:       VecScale(d->BcX[d->size_BcX+i], 1.0/norm);
561:     }
562:     d->size_BcX+= d->V_tra_s;
563:   }

565:   /* Update G <- MTY' * (G * MTX) */
566:   if (d->G) {
567:     tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
568:     cMT = d->V_tra_e - tra_s;
569:     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);
570:     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);
571:     d->size_G = cMT;
572:     d->cX_in_G = d->cX_in_V;
573:   }

575:   return(0);
576: }


581: /* DvdMult_copy_func: 2 */
582: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
583:                                        DvdMult_copy_func **sr)
584: {
585:   PetscErrorCode  ierr;
586:   Vec             *W = d->W?d->W:d->V, *BV = d->BV?d->BV:d->V;


590:   if (!d->G || d->V_new_s == d->V_new_e) { return(0); }

592:   /* G = [G               W(old)'*BV(new);
593:           W(new)'*BV(old) W(new)'*BV(new) ],
594:      being old=0:V_new_s-1, new=V_new_s:V_new_e-1 */
595:   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)++);
596:   d->size_G = d->V_new_e+d->cX_in_G;

598:   return(0);
599: }

601: /* in complex, d->size_H real auxiliar values are needed */
604: PetscErrorCode dvd_calcpairs_projeig_eig(dvdDashboard *d)
605: {
606:   PetscErrorCode  ierr;
607:   PetscReal       *w;
608: #if defined(PETSC_USE_COMPLEX)
609:   PetscInt        i;
610: #endif


614:   /* S <- H */
615:   d->ldS = d->ldpX = d->size_H;
616:   SlepcDenseCopyTriang(d->S, DVD_MAT_LTRIANG, d->size_H, d->H, d->sH, d->ldH,
617:                               d->size_H, d->size_H);

619:   /* S = pX' * L * pX */
620: #if !defined(PETSC_USE_COMPLEX)
621:   w = d->eigr-d->cX_in_H;
622: #else
623:   w = (PetscReal*)d->auxS;
624: #endif
625:   EPSDenseHEP(d->size_H, d->S, d->ldS, w, d->pX);
626: #if defined(PETSC_USE_COMPLEX)
627:   for (i=0; i<d->size_H; i++) d->eigr[i-d->cX_in_H] = w[i];
628: #endif

630:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_eig;

632:   return(0);
633: }


638: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d)
639: {
640:   PetscErrorCode  ierr;


644:   /* S <- H */
645:   d->ldS = d->ldpX = d->size_H;
646:   SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
647:                               d->size_H, d->size_H);

649:   /* S = pX' * H * pX */
650:   EPSDenseHessenberg(d->size_H, 0, d->S, d->ldS, d->pX);
651:   EPSDenseSchur(d->size_H, 0, d->S, d->ldS, d->pX, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
652: 

654:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;

656:   return(0);
657: }

661: /*
662:   auxS(dgges) = size_H (beta) + 8*size_H+16 (work)
663:   auxS(zgges) = size_H (beta) + 1+2*size_H (work) + 8*size_H (rwork)
664: */
665: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d)
666: {
667: #if defined(SLEPC_MISSING_LAPACK_GGES)
669:   SETERRQ(((PetscObject)(d->eps))->comm,PETSC_ERR_SUP,"GGES - Lapack routine is unavailable.");
670: #else
671:   PetscErrorCode  ierr;
672:   PetscScalar     *beta = d->auxS;
673: #if !defined(PETSC_USE_COMPLEX)
674:   PetscScalar     *auxS = beta + d->size_H;
675:   PetscBLASInt    n_auxS = d->size_auxS - d->size_H;
676: #else
677:   PetscReal       *auxR = (PetscReal*)(beta + d->size_H);
678:   PetscScalar     *auxS = (PetscScalar*)(auxR+8*d->size_H);
679:   PetscBLASInt    n_auxS = d->size_auxS - 9*d->size_H;
680: #endif
681:   PetscInt        i;
682:   PetscBLASInt    info,n,a;


692:   /* S <- H, T <- G */
693:   d->ldS = d->ldT = d->ldpX = d->ldpY = d->size_H;
694:   SlepcDenseCopyTriang(d->S, 0, d->size_H, d->H, d->sH, d->ldH,
695:                               d->size_H, d->size_H);
696:   SlepcDenseCopyTriang(d->T, 0, d->size_H, d->G, d->sG, d->ldH,
697:                               d->size_H, d->size_H);

699:   /* S = Z'*H*Q, T = Z'*G*Q */
700:   n = d->size_H;
701: #if !defined(PETSC_USE_COMPLEX)
702:   LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
703:               &a, d->eigr-d->cX_in_H, d->eigi-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
704:               auxS, &n_auxS, PETSC_NULL, &info);
705: #else
706:   LAPACKgges_(d->pY?"V":"N", "V", "N", PETSC_NULL, &n, d->S, &n, d->T, &n,
707:               &a, d->eigr-d->cX_in_H, beta, d->pY, &n, d->pX, &n,
708:               auxS, &n_auxS, auxR, PETSC_NULL, &info);
709: #endif
710:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack GGES %d", info);

712:   /* eigr[i] <- eigr[i] / beta[i] */
713:   for (i=0; i<d->size_H; i++)
714:     d->eigr[i-d->cX_in_H] /= beta[i],
715:     d->eigi[i-d->cX_in_H] /= beta[i];

717:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;

719:   return(0);
720: #endif
721: }

725: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n)
726: {
727:   PetscErrorCode  ierr;


731:   EPSSortDenseHEP(d->eps, d->size_H, 0, d->eigr-d->cX_in_H, d->pX, d->ldpX);
732: 

734:   if (d->calcpairs_eigs_trans) {
735:     d->calcpairs_eigs_trans(d);
736:   }

738:   return(0);
739: }


744: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n)
745: {
746:   PetscErrorCode  ierr;

749:   if ((d->ldpX != d->size_H) ||
750:       ( d->T &&
751:         ((d->ldS != d->ldT) || (d->ldpX != d->ldpY) ||
752:          (d->ldpX != d->size_H)) ) ) {
753:      SETERRQ(PETSC_COMM_SELF,1, "Error before ordering eigenpairs");
754:   }

756:   if (d->T) {
757:     EPSSortDenseSchurGeneralized(d->eps, d->size_H, 0, n, d->S, d->T,
758:                                         d->ldS, d->pY, d->pX, d->eigr-d->cX_in_H,
759:                                         d->eigi-d->cX_in_H);
760:   } else {
761:     EPSSortDenseSchur(d->eps, d->size_H, 0, d->S, d->ldS, d->pX,
762:                              d->eigr-d->cX_in_H, d->eigi-d->cX_in_H);
763:   }

765:   if (d->calcpairs_eigs_trans) {
766:     d->calcpairs_eigs_trans(d);
767:   }

769: #if defined(PETSC_USE_COMPLEX)
770:   if (d->T) {
771:     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);
772:   }
773: #endif

775:   return(0);
776: }


781: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
782:    the norm associated to the Schur pair, where i = r_s..r_e
783: */
784: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
785:                              Vec *R)
786: {
787:   PetscInt        i;
788:   PetscErrorCode  ierr;
789:   Vec             *BV = d->BV?d->BV:d->V;


793:   for (i=r_s; i<r_e; i++) {
794:     /* nX(i) <- ||X(i)|| */
795:     if (d->correctXnorm) {
796:       /* R(i) <- V*pX(i) */
797:       SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
798:         &d->V[-d->cX_in_H], d->size_V+d->cX_in_H,
799:         &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
800:       VecNorm(R[i-r_s], NORM_2, &d->nX[i]);
801:     } else
802:       d->nX[i] = 1.0;

804:     /* R(i-r_s) <- AV*pX(i) */
805:     SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
806:       &d->AV[-d->cX_in_H], d->size_AV+d->cX_in_H,
807:       &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);

809:     /* R(i-r_s) <- R(i-r_s) - eigr(i)*BV*pX(i) */
810:     SlepcUpdateVectorsZ(&R[i-r_s], 1.0, -d->eigr[i+d->cX_in_H],
811:       &BV[-d->cX_in_H], d->size_V+d->cX_in_H,
812:       &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
813:   }

815:   d->calcpairs_proj_res(d, r_s, r_e, R);

817:   return(0);
818: }

822: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
823:                                       PetscInt r_e, Vec *R)
824: {
825:   PetscInt        i;
826:   PetscErrorCode  ierr;
827:   PetscBool       lindep;
828:   Vec             *cX;


832:   /* If exists the BcX, R <- orth(BcX, R), nR[i] <- ||R[i]|| */
833:   if (d->BcX)
834:     cX = d->BcX;

836:   /* If exists left subspace, R <- orth(cY, R), nR[i] <- ||R[i]|| */
837:   else if (d->cY)
838:     cX = d->cY;

840:   /* If fany configurations, R <- orth(cX, R), nR[i] <- ||R[i]|| */
841:   else if (!(DVD_IS(d->sEP, DVD_EP_STD) && DVD_IS(d->sEP, DVD_EP_HERMITIAN)))
842:     cX = d->cX;

844:   /* Otherwise, nR[i] <- ||R[i]|| */
845:   else
846:     cX = PETSC_NULL;

848:   if (cX) for (i=0; i<r_e-r_s; i++) {
849:     IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX, PETSC_NULL,
850:                            cX, R[i], PETSC_NULL, &d->nR[r_s+i], &lindep);
851: 
852:     if(lindep || (d->nR[r_s+i] < PETSC_MACHINE_EPSILON)) {
853:       PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %G!\n",r_s+i,d->nR[r_s+i]);
854:     }

856:   } else {
857:     for(i=0; i<r_e-r_s; i++) {
858:       VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
859:     }
860:     for(i=0; i<r_e-r_s; i++) {
861:       VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
862:     }
863:   }

865:   return(0);
866: }

870: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
871:    the norm associated to the eigenpair, where i = r_s..r_e
872:    R, vectors of Vec of size r_e-r_s,
873:    auxV, PetscMax(r_e+cX_in_H, 2*(r_e-r_s)) vectors,
874:    auxS, auxiliar vector of size (d->size_cX+r_e)^2+6(d->size_cX+r_e)+(r_e-r_s)*d->size_H
875: */
876: PetscErrorCode dvd_calcpairs_eig_res_0(dvdDashboard *d,PetscInt r_s,PetscInt r_e,Vec *R)
877: {
878:   PetscInt        i,ldpX,size_in;
879:   PetscErrorCode  ierr;
880:   Vec             *Bx;
881:   PetscScalar     *pX,*pX0;
882:   DvdReduction    r;
883:   DvdReductionChunk
884:                   ops[2];
885:   DvdMult_copy_func
886:                   sr[2];
887: #if !defined(PETSC_USE_COMPLEX)
888:   PetscScalar     b[8];
889:   Vec             X[4];
890: #endif


894:   /* Quick return */
895:   if (!d->cS) { return(0); }

897:   size_in = (d->size_cX+r_e)*(d->cX_in_AV+r_e)*(d->cT?2:1);
898:   /* Check consistency */
899:   if (d->size_auxV < PetscMax(2*(r_e-r_s),d->cX_in_AV+r_e) || d->size_auxS <
900:       (d->size_cX+r_e)*(d->size_cX+r_e) /* pX */ + PetscMax(PetscMax(
901:         (d->size_cX+r_e)*6 /* dvd_compute_eigenvectors */,
902:         d->size_H*(r_e-r_s) /* pX0 */),
903:         2*size_in /* SlepcAllReduceSum */ )) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

905:   /* Prepare reductions */
906:   SlepcAllReduceSumBegin(ops,2,d->auxS,d->auxS+size_in,size_in,&r,((PetscObject)d->V[0])->comm);

908:   /* auxV <- AV * pX(0:r_e+cX_in_H) */
909:   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);

911:   /* cS(:, size_cS:) <- cX' * auxV */
912:   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]);

914:   if (d->cT) {
915:     /* R <- BV * pX(0:r_e+cX_in_H) */
916:     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);
917: 
918:     /* cS(:, size_cS:) <- cX' * auxV */
919:     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]);
920:   }

922:   /* Do reductions */
923:   SlepcAllReduceSumEnd(&r);

925:   /* qX <- eig(S,T) */
926:   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);
927:   ldpX = d->size_cX+r_e;
928:   pX0 = d->auxS;
929:   EPSCleanDenseSchur(d->size_cX+r_e,0,d->cS,d->ldcS,d->cT,d->ldcT,d->ceigi,pX,ldpX,PETSC_FALSE);
930:   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);

932:   /* pX[i] <- pX[i] / ||pX[i]|| */
933:   pX = &pX[ldpX*d->size_cX+r_s];
934:   SlepcDenseNorm(pX,ldpX,d->size_cX+r_e,r_e-r_s,d->eigi+r_s);

936:   /* pX0 <- d->pX(0:d->cX_in_AV+r_e-1) * pX(size_cX-cX_in_H:) */
937:   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);

939:   /* auxV <- cX(0:size_cX-cX_in_AV)*pX + V*pX0 */
940:   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);
941:   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);
942:   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);
943:   /* R <- A*auxV */
944:   for (i=0; i<r_e-r_s; i++) {
945:     MatMult(d->A,d->auxV[i],R[i]);
946:   }
947:   /* Bx <- B*auxV */
948:   if (d->B) {
949:     Bx = &d->auxV[r_e-r_s];
950:     for (i=0; i<r_e-r_s; i++) {
951:       MatMult(d->B,d->auxV[i],Bx[i]);
952:     }
953:   } else {
954:     Bx = d->auxV;
955:   }
956:   /* R <- (A - eig*B)*V*pX */
957:   for(i=0; i<r_e-r_s; i++) {
958: #if !defined(PETSC_USE_COMPLEX)
959:     if (d->eigi[r_s+i] != 0.0) {
960:       /* [Ax_i Ax_i+1 Bx_i Bx_i+1]*= [   1        0 
961:                                          0        1
962:                                       -eigr_i -eigi_i
963:                                        eigi_i -eigr_i] */
964:       b[0] = b[5] = 1.0;
965:       b[2] = b[7] = -d->eigr[r_s+i];
966:       b[6] = -(b[3] = d->eigi[r_s+i]);
967:       b[1] = b[4] = 0.0;
968:       X[0] = R[i]; X[1] = R[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
969:       SlepcUpdateVectorsD(X,4,1.0,b,4,4,2,d->auxS,d->size_auxS);
970:       i++;
971:     } else
972: #endif
973:     {
974:       /* R <- Ax -eig*Bx */
975:       VecAXPBY(R[i], -d->eigr[r_s+i], 1.0, Bx[i]);
976:     }
977:   }

979:   /* nR <- ||R|| */
980:   for(i=0; i<r_e-r_s; i++) {
981:     VecNormBegin(R[i],NORM_2,&d->nR[r_s+i]);
982:   }
983:   for(i=0; i<r_e-r_s; i++) {
984:     VecNormEnd(R[i],NORM_2,&d->nR[r_s+i]);
985:   }

987:   return(0);
988: }


991: /**** Pattern routines ********************************************************/

993: /* BV <- BV*MTX */
996: 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)
997: {
998:   PetscErrorCode  ierr;
999:   PetscInt        cMT, rm, cp, tra_s, i;
1000:   Vec             *nBV;


1004:   if (!real_BV || !*BV || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }

1006:   if (d->V_tra_s > d->max_cX_in_proj && !BV_shift) {
1007:     tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
1008:     cMT = d->V_tra_e - tra_s;
1009:     rm = d->V_tra_s - tra_s;
1010:     cp = PetscMin(d->max_cX_in_proj - rm, *cX_in_proj);
1011:     nBV = real_BV+d->max_cX_in_proj;
1012:     /* BV(-cp-rm:-1-rm) <- BV(-cp:-1) */
1013:     for (i=-cp; i<0; i++) {
1014:       VecCopy((*BV)[i], nBV[i-rm]);
1015:     }
1016:     /* BV(-rm:) <- BV*MTX(tra_s:V_tra_e-1) */
1017:     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);
1018:     *size_BV = d->V_tra_e  - d->V_tra_s;
1019:     *max_size_BV-= nBV - *BV;
1020:     *BV = nBV;
1021:     if (cX_in_proj && d->max_cX_in_proj>0) *cX_in_proj = cp+rm;
1022:   } else if (d->V_tra_s <= d->max_cX_in_proj || BV_shift) {
1023:     /* [BcX BV] <- [BcX BV*MTX] */
1024:     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);
1025:     *BV+= d->V_tra_s-*cX_in_proj;
1026:     *max_size_BV-= d->V_tra_s-*cX_in_proj;
1027:     *size_BV = d->V_tra_e  - d->V_tra_s;
1028:     if (size_cBV && BV_shift) *size_cBV = *BV - real_BV;
1029:     if (d->max_cX_in_proj>0) *cX_in_proj = PetscMin(*BV - real_BV, d->max_cX_in_proj);
1030:   } else { /* !BV_shift */
1031:     /* BV <- BV*MTX(V_tra_s:) */
1032:     SlepcUpdateVectorsZ(*BV, 0.0, 1.0, *BV, *size_BV,
1033:       &MTX[d->V_tra_s*d->ldMTX], d->ldMTX, d->size_MT, d->V_tra_e-d->V_tra_s);
1034: 
1035:     *size_BV = d->V_tra_e - d->V_tra_s;
1036:   }

1038:   return(0);
1039: }