Actual source code: dvd_calcpairs.c

  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Step: calc the best eigenpairs in the subspace V.

  6:   For that, performs these steps:
  7:     1) Update W <- A * V
  8:     2) Update H <- V' * W
  9:     3) Obtain eigenpairs of H
 10:     4) Select some eigenpairs
 11:     5) Compute the Ritz pairs of the selected ones

 13:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 14:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 15:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

 17:    This file is part of SLEPc.
 18:       
 19:    SLEPc is free software: you can redistribute it and/or modify it under  the
 20:    terms of version 3 of the GNU Lesser General Public License as published by
 21:    the Free Software Foundation.

 23:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 24:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 25:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 26:    more details.

 28:    You  should have received a copy of the GNU Lesser General  Public  License
 29:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 30:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 31: */

 33:  #include davidson.h
 34: #include <slepcblaslapack.h>

 36: PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d);
 37: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d);
 38: PetscErrorCode dvd_calcpairs_projeig_eig(dvdDashboard *d);
 39: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d);
 40: PetscErrorCode dvd_calcpairs_projeig_qz_gen(dvdDashboard *d);
 41: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n);
 42: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n);
 43: PetscErrorCode dvd_calcpairs_X(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
 44:                                Vec *X);
 45: PetscErrorCode dvd_calcpairs_Y(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
 46:                                Vec *Y);
 47: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
 48:                                    Vec *R);
 49: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
 50:                                       PetscInt r_e, Vec *R);
 51: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
 52:                                       DvdMult_copy_func **sr);
 53: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d);
 54: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
 55:                                       DvdMult_copy_func **sr);
 56: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d);
 57: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d);
 58: PetscErrorCode dvd_calcpairs_updateAV1(dvdDashboard *d, DvdReduction *r,
 59:                                        DvdMult_copy_func **sr);
 60: PetscErrorCode dvd_calcpairs_updateBV0(dvdDashboard *d);
 61: PetscErrorCode dvd_calcpairs_updateBV1(dvdDashboard *d, DvdReduction *r,
 62:                                        DvdMult_copy_func **sr);
 63: PETSC_STATIC_INLINE PetscErrorCode dvd_calcpairs_updateBV0_gen(dvdDashboard *d,Vec *real_BV,PetscInt *size_cX,Vec **BV,PetscInt *size_BV,PetscInt *max_size_BV,PetscBool BV_shift,PetscInt *cX_in_proj,PetscScalar *MTX);

 65: /**** Control routines ********************************************************/
 68: PetscErrorCode dvd_calcpairs_qz(dvdDashboard *d, dvdBlackboard *b,
 69:                                 orthoV_type_t orth, IP ipI,
 70:                                 PetscInt cX_proj, PetscBool harm)
 71: {
 72:   PetscErrorCode  ierr;
 73:   PetscInt        i,max_cS;
 74:   PetscBool       std_probl,her_probl;


 78:   std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
 79:   her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;

 81:   /* Setting configuration constrains */
 82: #if !defined(PETSC_USE_COMPLEX)
 83:   /* if the last converged eigenvalue is complex its conjugate pair is also
 84:      converged */
 85:   b->max_nev = PetscMax(b->max_nev, d->nev+(her_probl && !d->B?0:1));
 86: #else
 87:   b->max_nev = PetscMax(b->max_nev, d->nev);
 88: #endif
 89:   b->max_size_proj = PetscMax(b->max_size_proj, b->max_size_V+cX_proj);
 90:   d->size_real_V = b->max_size_V+b->max_nev;
 91:   d->W_shift = d->B?PETSC_TRUE:PETSC_FALSE;
 92:   d->size_real_W = harm?(b->max_size_V+(d->W_shift?b->max_nev:b->max_size_cP)):0;
 93:   d->size_real_AV = b->max_size_V+b->max_size_cP;
 94:   d->size_BDS = 0;
 95:   if (d->B && her_probl && (orth == DVD_ORTHOV_I || orth == DVD_ORTHOV_BOneMV)) {
 96:     d->size_real_BV = b->size_V; d->BV_shift = PETSC_TRUE;
 97:     if (orth == DVD_ORTHOV_BOneMV) d->size_BDS = d->eps->nds;
 98:   } else if (d->B) {
 99:     d->size_real_BV = b->max_size_V + b->max_size_P; d->BV_shift = PETSC_FALSE;
100:   } else {
101:     d->size_real_BV = 0; d->BV_shift = PETSC_FALSE;
102:   }
103:   b->own_vecs+= d->size_real_V + d->size_real_W + d->size_real_AV +
104:                 d->size_real_BV + d->size_BDS;
105:   b->own_scalars+= b->max_size_proj*b->max_size_proj*2*(std_probl?1:2) +
106:                                               /* H, G?, S, T? */
107:                    b->max_size_proj*b->max_size_proj*(std_probl?1:2) +
108:                                               /* pX, pY? */
109:                    b->max_nev*b->max_nev*(her_probl?0:(!d->B?1:2));
110:                                                 /* cS?, cT? */
111:   b->max_size_auxV = PetscMax(b->max_size_auxV, b->max_size_X);
112:                                                 /* updateV0 */
113:   max_cS = PetscMax(b->max_size_X,cX_proj);
114:   b->max_size_auxS = PetscMax(PetscMax(
115:     b->max_size_auxS,
116:     b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
117:       max_cS*b->max_nev*(her_probl?0:(!d->B?1:2)) + /* updateV0,W0 */
118:                                                      /* SlepcReduction: in */
119:       PetscMax(
120:         b->max_size_X*b->max_size_proj*2*(std_probl?1:2) + /* updateAV1,BV1 */
121:           max_cS*b->max_nev*(her_probl?0:(!d->B?1:2)), /* updateV0,W0 */
122:                                                     /* SlepcReduction: out */
123:         PetscMax(
124:           b->max_size_proj*b->max_size_proj, /* updateAV0,BV0 */
125:           b->max_size_proj+b->max_nev))), /* dvd_orth */
126:     std_probl?0:(b->max_size_proj*11+16) /* projeig */);
127: #if defined(PETSC_USE_COMPLEX)
128:   b->max_size_auxS = PetscMax(b->max_size_auxS, b->max_size_V);
129:                                            /* dvd_calcpairs_projeig_eig */
130: #endif

132:   /* Setup the step */
133:   if (b->state >= DVD_STATE_CONF) {
134:     d->max_cX_in_proj = cX_proj;
135:     d->max_size_P = b->max_size_P;
136:     d->real_V = b->free_vecs; b->free_vecs+= d->size_real_V;
137:     if (harm) {
138:       d->real_W = b->free_vecs; b->free_vecs+= d->size_real_W;
139:     } else {
140:       d->real_W = PETSC_NULL;
141:     }
142:     d->real_AV = d->AV = b->free_vecs; b->free_vecs+= d->size_real_AV;
143:     d->max_size_proj = b->max_size_proj;
144:     d->real_H = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
145:     d->ldH = b->max_size_proj;
146:     d->pX = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
147:     d->S = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
148:     if (!her_probl) {
149:       d->cS = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
150:       d->max_size_cS = d->ldcS = b->max_nev;
151:     } else {
152:       d->cS = PETSC_NULL;
153:       d->max_size_cS = d->ldcS = 0;
154:     }
155:     d->ipV = ipI;
156:     d->ipW = ipI;
157:     if (orth == DVD_ORTHOV_BOneMV) {
158:       d->BDS = b->free_vecs; b->free_vecs+= d->eps->nds;
159:       for (i=0; i<d->eps->nds; i++) {
160:         MatMult(d->B, d->eps->DS[i], d->BDS[i]);
161:       }
162:     } else
163:       d->BDS = PETSC_NULL;
164:     if (d->B) {
165:       d->real_BV = b->free_vecs; b->free_vecs+= d->size_real_BV;
166:     } else {
167:       d->size_real_BV = 0;
168:       d->real_BV = PETSC_NULL;
169:       d->BV_shift = PETSC_FALSE;
170:     }
171:     if (!std_probl) {
172:       d->real_G = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
173:       d->T = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
174:       d->pY = b->free_scalars; b->free_scalars+= b->max_size_proj*b->max_size_proj;
175:     } else {
176:       d->real_G = PETSC_NULL;
177:       d->T = PETSC_NULL;
178:       d->pY = PETSC_NULL;
179:     }
180:     if (d->B && !her_probl) {
181:       d->cT = b->free_scalars; b->free_scalars+= b->max_nev*b->max_nev;
182:       d->ldcT = b->max_nev;
183:     } else {
184:       d->cT = PETSC_NULL;
185:       d->ldcT = 0;
186:     }

188:     d->calcPairs = dvd_calcpairs_proj;
189:     d->calcpairs_residual = dvd_calcpairs_res_0;
190:     d->calcpairs_proj_res = dvd_calcpairs_proj_res;
191:     d->calcpairs_selectPairs = PETSC_NULL;
192:     d->ipI = ipI;
193:     DVD_FL_ADD(d->startList, dvd_calcpairs_qz_start);
194:   }

196:   return(0);
197: }


202: PetscErrorCode dvd_calcpairs_qz_start(dvdDashboard *d)
203: {
204:   PetscBool her_probl;
205:   PetscInt  i;

208:   her_probl = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE;

210:   d->size_V = 0;
211:   d->V = d->real_V;
212:   d->cX = d->real_V;
213:   d->size_cX = 0;
214:   d->max_size_V = d->size_real_V;
215:   d->W = d->real_W;
216:   d->max_size_W = d->size_real_W;
217:   d->size_W = 0;
218:   d->size_AV = 0;
219:   d->AV = d->real_AV;
220:   d->max_size_AV = d->size_real_AV;
221:   d->size_H = 0;
222:   d->H = d->real_H;
223:   if (d->cS) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cS[i] = 0.0;
224:   d->size_BV = 0;
225:   d->BV = d->real_BV;
226:   d->max_size_BV = d->size_real_BV;
227:   d->size_G = 0;
228:   d->G = d->real_G;
229:   if (d->cT) for (i=0; i<d->max_size_cS*d->max_size_cS; i++) d->cT[i] = 0.0;
230:   d->cY = d->B && !her_probl ? d->W : PETSC_NULL;
231:   d->BcX = d->orthoV_type == DVD_ORTHOV_I && d->B && her_probl ? d->BcX : PETSC_NULL;
232:   d->size_cY = 0;
233:   d->size_BcX = 0;
234:   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;
235:   return(0);
236: }


241: PetscErrorCode dvd_calcpairs_proj(dvdDashboard *d)
242: {
243:   PetscErrorCode  ierr;
244:   DvdReduction    r;
245:   static const PetscInt MAX_OPS = 7;
246:   DvdReductionChunk
247:                   ops[MAX_OPS];
248:   DvdMult_copy_func
249:                   sr[MAX_OPS], *sr0 = sr;
250:   PetscInt        size_in, i;
251:   PetscScalar     *in = d->auxS, *out;
252:   PetscBool       stdp;


256:   stdp = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
257:   size_in =
258:     (d->size_cX+d->V_tra_s-d->cX_in_H)*d->V_tra_s*(d->cT?2:(d->cS?1:0)) + /* updateV0,W0 */
259:     (d->size_H*(d->V_new_e-d->V_new_s)*2+
260:       (d->V_new_e-d->V_new_s)*(d->V_new_e-d->V_new_s))*(!stdp?2:1); /* updateAV1,BV1 */
261: 
262:   out = in+size_in;

264:    /* Check consistency */
265:    if (2*size_in > d->size_auxS) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

267:   /* Prepare reductions */
268:   SlepcAllReduceSumBegin(ops, MAX_OPS, in, out, size_in, &r,
269:                                 ((PetscObject)d->V[0])->comm);
270:   /* Allocate size_in */
271:   d->auxS+= size_in;
272:   d->size_auxS-= size_in;

274:   /* Update AV, BV, W and the projected matrices */
275:   /* 1. S <- S*MT */
276:   dvd_calcpairs_updateV0(d, &r, &sr0);
277:   dvd_calcpairs_updateW0(d, &r, &sr0);
278:   dvd_calcpairs_updateAV0(d);
279:   dvd_calcpairs_updateBV0(d);
280:   /* 2. V <- orth(V, V_new) */
281:   dvd_calcpairs_updateV1(d);
282:   /* 3. AV <- [AV A * V(V_new_s:V_new_e-1)] */
283:   /* Check consistency */
284:   if (d->size_AV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
285:   for (i=d->V_new_s; i<d->V_new_e; i++) {
286:     MatMult(d->A, d->V[i], d->AV[i]);
287:   }
288:   d->size_AV = d->V_new_e;
289:   /* 4. BV <- [BV B * V(V_new_s:V_new_e-1)] */
290:   if (d->B && d->orthoV_type != DVD_ORTHOV_BOneMV) {
291:     /* Check consistency */
292:     if (d->size_BV != d->V_new_s) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
293:     for (i=d->V_new_s; i<d->V_new_e; i++) {
294:       MatMult(d->B, d->V[i], d->BV[i]);
295:     }
296:     d->size_BV = d->V_new_e;
297:   }
298:   /* 5 <- W <- [W f(AV,BV)] */
299:   dvd_calcpairs_updateW1(d);
300:   dvd_calcpairs_updateAV1(d, &r, &sr0);
301:   dvd_calcpairs_updateBV1(d, &r, &sr0);

303:   /* Deallocate size_in */
304:   d->auxS-= size_in;
305:   d->size_auxS+= size_in;

307:   /* Do reductions */
308:   SlepcAllReduceSumEnd(&r);

310:   /* Perform the transformation on the projected problem */
311:   if (d->calcpairs_proj_trans) {
312:     d->calcpairs_proj_trans(d);
313:   }

315:   d->V_tra_s = d->V_tra_e = 0;
316:   d->V_new_s = d->V_new_e;

318:   /* Solve the projected problem */
319:   if (DVD_IS(d->sEP, DVD_EP_STD)) {
320:     if (DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
321:       dvd_calcpairs_projeig_eig(d);
322:     } else {
323:       dvd_calcpairs_projeig_qz_std(d);
324:     }
325:   } else {
326:     dvd_calcpairs_projeig_qz_gen(d);
327:   }

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

340:   return(0);
341: }

343: /**** Basic routines **********************************************************/

347: /* auxV: V_tra_s, DvdMult_copy_func: 1 */
348: PetscErrorCode dvd_calcpairs_updateV0(dvdDashboard *d, DvdReduction *r,
349:                                       DvdMult_copy_func **sr)
350: {
351:   PetscErrorCode  ierr;
352:   PetscInt        rm;


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

358:   /* cX <- [cX V*MT(0:V_tra_s-1)], V <- V*MT(V_tra_s:V_tra_e) */
359:   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);

361:   /* Udpate cS for standard problems */
362:   if (d->cS && !d->cT && !d->cY && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
363:     /* Check consistency */
364:     if (d->size_cS+d->V_tra_s != d->size_cX) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

366:     /* auxV <- AV * MTX(0:V_tra_e-1) */
367:     rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
368:     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);

370:     /* cS(:, size_cS:) <- cX' * auxV */
371:     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)++);
372:     d->size_cS+= d->V_tra_s-rm;
373:   }
374: 
375:   return(0);
376: }


381: /* auxS: size_cX+V_new_e+1 */
382: PetscErrorCode dvd_calcpairs_updateV1(dvdDashboard *d)
383: {
384:   PetscErrorCode  ierr;
385:   Vec             *cX = d->BcX? d->BcX : ( (d->cY && !d->W)? d->cY : d->cX );


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

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

394:   /* V <- gs([cX V(0:V_new_s-1)], V(V_new_s:V_new_e-1)) */
395:   if (d->orthoV_type == DVD_ORTHOV_BOneMV) {
396:     dvd_BorthV(d->ipV, d->eps->DS, d->BDS, d->eps->nds, d->cX, d->real_BV,
397:                       d->size_cX, d->V, d->BV, d->V_new_s, d->V_new_e,
398:                       d->auxS, d->eps->rand);
399:     d->size_BV = d->V_new_e;
400:   } else {
401:     dvd_orthV(d->ipV, d->eps->DS, d->eps->nds, cX, d->size_cX, d->V,
402:                    d->V_new_s, d->V_new_e, d->auxS, d->eps->rand);
403: 
404:   }
405:   d->size_V = d->V_new_e;

407:   return(0);
408: }

412: /* auxV: V_tra_s, DvdMult_copy_func: 2 */
413: PetscErrorCode dvd_calcpairs_updateW0(dvdDashboard *d, DvdReduction *r,
414:                                       DvdMult_copy_func **sr)
415: {
416:   PetscErrorCode  ierr;
417:   PetscInt        rm;


421:   if (!d->W || (d->V_tra_s == 0 && d->V_tra_e == 0)) { return(0); }

423:   /* cY <- [cY W*MTY(0:V_tra_s-1)], W <- W*MTY(V_tra_s:V_tra_e) */
424:   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);

426:   /* Udpate cS and cT */
427:   if (d->cY && d->cT && (d->V_tra_s > d->max_cX_in_proj || d->size_cX >= d->nev)) {
428:     /* Check consistency */
429:     if (d->size_cS+d->V_tra_s != d->size_cY) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

431:     /* auxV <- AV * MTX(0:V_tra_e-1) */
432:     rm = d->size_cX>=d->nev?0:d->max_cX_in_proj;
433:     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);

435:     /* cS(:, size_cS:) <- cY' * auxV */
436:     VecsMultS(&d->cS[d->ldcS*d->size_cS], 0, d->ldcS, d->cY, 0, d->size_cY-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);

438:     /* auxV <- BV * MTX(0:V_tra_e-1) */
439:     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);

441:     /* cT(:, size_cS:) <- cY' * auxV */
442:     VecsMultS(&d->cT[d->ldcS*d->size_cS], 0, d->ldcS, d->cY, 0, d->size_cY-rm, d->auxV, 0, d->V_tra_s-rm, r, (*sr)++);
443: 
444:     d->size_cS+= d->V_tra_s-rm;
445:     d->size_cT+= d->V_tra_s-rm;
446:   }

448:   return(0);
449: }


454: /* auxS: size_cX+V_new_e+1 */
455: PetscErrorCode dvd_calcpairs_updateW1(dvdDashboard *d)
456: {
457:   PetscErrorCode  ierr;


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

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

466:   /* Update W */
467:   d->calcpairs_W(d);

469:   /* W <- gs([cY W(0:V_new_s-1)], W(V_new_s:V_new_e-1)) */
470:   dvd_orthV(d->ipW, PETSC_NULL, 0, d->cY, d->size_cY, d->W, d->V_new_s,
471:                    d->V_new_e, d->auxS, d->eps->rand);
472: 
473:   d->size_W = d->V_new_e;

475:   return(0);
476: }

480: /* auxS: size_H*(V_tra_e-V_tra_s) */
481: PetscErrorCode dvd_calcpairs_updateAV0(dvdDashboard *d)
482: {
483:   PetscErrorCode  ierr;
484:   PetscScalar     *MTY = d->W?d->MTY:d->MTX;
485:   PetscInt        cMT, tra_s;


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

491:   /* AV(V_tra_s-cp-1:) = cAV*MTX(V_tra_s:) */
492:   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);
493:   tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
494:   cMT = d->V_tra_e - tra_s;

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

502:   return(0);
503: }


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


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

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

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

530:   return(0);
531: }

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


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

548:   /* BV <- BV*MTX */
549:   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);

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

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

573:   return(0);
574: }


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


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

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

596:   return(0);
597: }

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


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

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

628:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_eig;

630:   return(0);
631: }


636: PetscErrorCode dvd_calcpairs_projeig_qz_std(dvdDashboard *d)
637: {
638:   PetscErrorCode  ierr;


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

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

652:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;

654:   return(0);
655: }

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


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

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

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

715:   d->calcpairs_selectPairs = dvd_calcpairs_selectPairs_qz;

717:   return(0);
718: #endif
719: }

723: PetscErrorCode dvd_calcpairs_selectPairs_eig(dvdDashboard *d, PetscInt n)
724: {
725:   PetscErrorCode  ierr;


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

732:   if (d->calcpairs_eigs_trans) {
733:     d->calcpairs_eigs_trans(d);
734:   }

736:   return(0);
737: }


742: PetscErrorCode dvd_calcpairs_selectPairs_qz(dvdDashboard *d, PetscInt n)
743: {
744:   PetscErrorCode  ierr;
745: #if defined(PETSC_USE_COMPLEX)
746:   PetscScalar     s;
747:   PetscInt        i, j;
748: #endif

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:   /* Some functions need the diagonal elements in T be real */
772: #if defined(PETSC_USE_COMPLEX)
773:   if (d->T) for(i=0; i<d->size_H; i++)
774:     if (PetscImaginaryPart(d->T[d->ldT*i+i]) != 0.0) {
775:       s = PetscConj(d->T[d->ldT*i+i])/PetscAbsScalar(d->T[d->ldT*i+i]);
776:       for(j=0; j<=i; j++)
777:         d->T[d->ldT*i+j] = PetscRealPart(d->T[d->ldT*i+j]*s),
778:         d->S[d->ldS*i+j]*= s;
779:       for(j=0; j<d->size_H; j++) d->pX[d->ldpX*i+j]*= s;
780:     }
781: #endif

783:   return(0);
784: }


789: /* Compute the residual vectors R(i) <- (AV - BV*eigr(i))*pX(i), and also
790:    the norm, where i = r_s..r_e
791: */
792: PetscErrorCode dvd_calcpairs_res_0(dvdDashboard *d, PetscInt r_s, PetscInt r_e,
793:                              Vec *R)
794: {
795:   PetscInt        i;
796:   PetscErrorCode  ierr;
797:   Vec             *BV = d->BV?d->BV:d->V;


801:   for (i=r_s; i<r_e; i++) {
802:     /* nX(i) <- ||X(i)|| */
803:     if (d->correctXnorm) {
804:       /* R(i) <- V*pX(i) */
805:       SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
806:         &d->V[-d->cX_in_H], d->size_V+d->cX_in_H,
807:         &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
808:       VecNorm(R[i-r_s], NORM_2, &d->nX[i]);
809:     } else
810:       d->nX[i] = 1.0;

812:     /* R(i-r_s) <- AV*pX(i) */
813:     SlepcUpdateVectorsZ(&R[i-r_s], 0.0, 1.0,
814:       &d->AV[-d->cX_in_H], d->size_AV+d->cX_in_H,
815:       &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);

817:     /* R(i-r_s) <- R(i-r_s) - eigr(i)*BV*pX(i) */
818:     SlepcUpdateVectorsZ(&R[i-r_s], 1.0, -d->eigr[i+d->cX_in_H],
819:       &BV[-d->cX_in_H], d->size_V+d->cX_in_H,
820:       &d->pX[d->ldpX*(i+d->cX_in_H)], d->ldpX, d->size_H, 1);
821:   }

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

825:   return(0);
826: }

830: PetscErrorCode dvd_calcpairs_proj_res(dvdDashboard *d, PetscInt r_s,
831:                                       PetscInt r_e, Vec *R)
832: {
833:   PetscInt        i;
834:   PetscErrorCode  ierr;
835:   PetscBool       lindep;
836:   Vec             *cX;


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

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

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

852:   /* Otherwise, nR[i] <- ||R[i]|| */
853:   else
854:     cX = PETSC_NULL;

856:   if (cX) for (i=0; i<r_e-r_s; i++) {
857:     IPOrthogonalize(d->ipI, 0, PETSC_NULL, d->size_cX, PETSC_NULL,
858:                            cX, R[i], PETSC_NULL, &d->nR[r_s+i], &lindep);
859: 
860:     if(lindep || (d->nR[r_s+i] < PETSC_MACHINE_EPSILON)) {
861:       PetscInfo2(d->eps,"The computed eigenvector residual %D is too low, %G!\n",r_s+i,d->nR[r_s+i]);
862:     }

864:   } else for(i=0; i<r_e-r_s; i++) {
865:     VecNorm(R[i], NORM_2, &d->nR[r_s+i]);
866:   }

868:   return(0);
869: }

871: /**** Pattern routines ********************************************************/

873: /* BV <- BV*MTX */
876: 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)
877: {
878:   PetscErrorCode  ierr;
879:   PetscInt        cMT, rm, cp, tra_s, i;
880:   Vec             *nBV;


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

886:   if (d->V_tra_s > d->max_cX_in_proj && !BV_shift) {
887:     tra_s = PetscMax(d->V_tra_s-d->max_cX_in_proj, 0);
888:     cMT = d->V_tra_e - tra_s;
889:     rm = d->V_tra_s - tra_s;
890:     cp = PetscMin(d->max_cX_in_proj - rm, *cX_in_proj);
891:     nBV = real_BV+d->max_cX_in_proj;
892:     /* BV(-cp-rm:-1-rm) <- BV(-cp:-1) */
893:     for (i=-cp; i<0; i++) {
894:       VecCopy((*BV)[i], nBV[i-rm]);
895:     }
896:     /* BV(-rm:) <- BV*MTX(tra_s:V_tra_e-1) */
897:     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);
898:     *size_BV = d->V_tra_e  - d->V_tra_s;
899:     *max_size_BV-= nBV - *BV;
900:     *BV = nBV;
901:     if (cX_in_proj && d->max_cX_in_proj>0) *cX_in_proj = cp+rm;
902:   } else if (d->V_tra_s <= d->max_cX_in_proj || BV_shift) {
903:     /* [BcX BV] <- [BcX BV*MTX] */
904:     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);
905:     *BV+= d->V_tra_s-*cX_in_proj;
906:     *max_size_BV-= d->V_tra_s-*cX_in_proj;
907:     *size_BV = d->V_tra_e  - d->V_tra_s;
908:     if (size_cBV && BV_shift) *size_cBV = *BV - real_BV;
909:     if (d->max_cX_in_proj>0) *cX_in_proj = PetscMin(*BV - real_BV, d->max_cX_in_proj);
910:   } else { /* !BV_shift */
911:     /* BV <- BV*MTX(V_tra_s:) */
912:     SlepcUpdateVectorsZ(*BV, 0.0, 1.0, *BV, *size_BV,
913:       &MTX[d->V_tra_s*d->ldMTX], d->ldMTX, d->size_MT, d->V_tra_e-d->V_tra_s);
914: 
915:     *size_BV = d->V_tra_e - d->V_tra_s;
916:   }

918:   return(0);
919: }