Actual source code: dvd_improvex.c

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

  4:   Step: improve the eigenvectors X

  6:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  7:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  8:    Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain

 10:    This file is part of SLEPc.
 11:       
 12:    SLEPc is free software: you can redistribute it and/or modify it under  the
 13:    terms of version 3 of the GNU Lesser General Public License as published by
 14:    the Free Software Foundation.

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

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

 26:  #include davidson.h
 27: #include <slepcvec.h>
 28: #include <slepcblaslapack.h>

 30: PetscErrorCode dvd_improvex_PfuncV(dvdDashboard *d, void *funcV, Vec *D,
 31:   PetscInt max_size_D, PetscInt r_s, PetscInt r_e, Vec *auxV,
 32:   PetscScalar *auxS);
 33: PetscErrorCode dvd_matmult_jd(Mat A, Vec in, Vec out);
 34: PetscErrorCode dvd_matgetvecs_jd(Mat A, Vec *right, Vec *left);
 35: PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d);
 36: PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d);
 37: PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d);
 38: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
 39:                                    PetscInt max_size_D, PetscInt r_s,
 40:                                    PetscInt r_e, PetscInt *size_D);
 41: PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d, PetscInt i_s,
 42:   PetscInt i_e, Vec **u, Vec **v, Vec *kr, Vec **auxV, PetscScalar **auxS,
 43:   PetscScalar *theta, PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY,
 44:   PetscInt ld);
 45: PetscErrorCode dvd_improvex_jd_proj_uv_KXX(dvdDashboard *d, PetscInt i_s,
 46:   PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
 47:   PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
 48: PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d, PetscInt i_s,
 49:   PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
 50:   PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld);
 51: PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d, PetscInt i,
 52:   PetscScalar* theta, PetscScalar* thetai,
 53:   PetscInt *maxits, PetscReal *tol);
 54: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
 55:   PetscScalar *pY, PetscInt ld_,
 56:   PetscScalar *auxS, PetscInt size_auxS);
 57: PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d, Vec *V, PetscInt cV, PetscScalar *auxS);

 59: #define size_Z (64*4)

 61: /**** JD update step (I - Kfg'/(g'Kf)) K(A - sB) (I - Kfg'/(g'Kf)) t = (I - Kfg'/(g'Kf))r  ****/

 63: typedef struct {
 64:   PetscInt size_X;
 65:   void
 66:     *old_improveX_data;   /* old improveX_data */
 67:   improveX_type
 68:     old_improveX;         /* old improveX */
 69:   KSP ksp;                /* correction equation solver */
 70:   Vec
 71:     friends,              /* reference vector for composite vectors */
 72:     *auxV;                /* auxiliar vectors */
 73:   PetscScalar *auxS,      /* auxiliar scalars */
 74:     *theta,
 75:     *thetai;              /* the shifts used in the correction eq. */
 76:   PetscInt maxits,        /* maximum number of iterations */
 77:     r_s, r_e,             /* the selected eigenpairs to improve */
 78:     ksp_max_size;         /* the ksp maximum subvectors size */
 79:   PetscReal tol,          /* the maximum solution tolerance */
 80:     lastTol,              /* last tol for dynamic stopping criterion */
 81:     fix;                  /* tolerance for using the approx. eigenvalue */
 82:   PetscBool
 83:     dynamic;              /* if the dynamic stopping criterion is applied */
 84:   dvdDashboard
 85:     *d;                   /* the currect dvdDashboard reference */
 86:   PC old_pc;              /* old pc in ksp */
 87:   Vec
 88:     *u,                   /* new X vectors */
 89:     *real_KZ,             /* original KZ */
 90:     *KZ;                  /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
 91:   PetscScalar
 92:    *XKZ,                  /* X'*KZ */
 93:    *iXKZ;                 /* inverse of XKZ */
 94:   PetscInt
 95:     ldXKZ,                /* leading dimension of XKZ */
 96:     size_iXKZ,            /* size of iXKZ */
 97:     ldiXKZ,               /* leading dimension of iXKZ */
 98:     size_KZ,              /* size of converged KZ */
 99:     size_real_KZ,         /* original size of KZ */
100:     size_cX,              /* last value of d->size_cX */
101:     old_size_X;           /* last number of improved vectors */
102:   PetscBLASInt
103:     *iXKZPivots;          /* array of pivots */
104: } dvdImprovex_jd;

106: #define _Ceil(A,B) ((A)/(B)+((A)%(B)==0?0:1))
107: #define FromIntToScalar(S) ((PetscInt)_Ceil((S)*sizeof(PetscBLASInt),sizeof(PetscScalar)))

111: PetscErrorCode dvd_improvex_jd(dvdDashboard *d, dvdBlackboard *b, KSP ksp,
112:                                PetscInt max_bs, PetscInt cX_impr, PetscBool dynamic)
113: {
114:   PetscErrorCode  ierr;
115:   dvdImprovex_jd  *data;
116:   PetscBool       useGD, herm = DVD_IS(d->sEP, DVD_EP_HERMITIAN)?PETSC_TRUE:PETSC_FALSE, std_probl = DVD_IS(d->sEP, DVD_EP_STD)?PETSC_TRUE:PETSC_FALSE;
117:   PC              pc;
118:   PetscInt        size_P,s=1;


122:   /* Setting configuration constrains */
123:   PetscTypeCompare((PetscObject)ksp, KSPPREONLY, &useGD);

125:   /* If the arithmetic is real and the problem is not Hermitian, then
126:      the block size is incremented in one */
127: #if !defined(PETSC_USE_COMPLEX)
128:   if (!herm) {
129:     max_bs++;
130:     b->max_size_P = PetscMax(b->max_size_P, 2);
131:     s = 2;
132:   } else
133: #endif
134:     b->max_size_P = PetscMax(b->max_size_P, 1);
135:   b->max_size_X = PetscMax(b->max_size_X, max_bs);
136:   size_P = b->max_size_P+cX_impr;
137:   b->max_size_auxV = PetscMax(b->max_size_auxV,
138:      b->max_size_X*(useGD?2:3)+ /* u, kr?, auxV */
139:      ((herm || !d->eps->trueres)?1:PetscMax(s*2,b->max_size_cX_proj+b->max_size_X))); /* testConv */
140: 
141:   b->own_scalars+= size_P*size_P; /* XKZ */
142:   b->max_size_auxS = PetscMax(b->max_size_auxS,
143:     (herm?0:1)*2*b->max_size_proj*b->max_size_proj + /* pX, pY */
144:     b->max_size_X*3 + /* theta, thetai */
145:     size_P*size_P + /* iXKZ */
146:     FromIntToScalar(size_P) + /* iXkZPivots */
147:     PetscMax(PetscMax(PetscMax(
148:       3*b->max_size_proj*b->max_size_X, /* dvd_improvex_apply_proj */
149:       8*cX_impr*b->max_size_X), /* dvd_improvex_jd_proj_cuv_KZX */
150:       (herm?0:1)*6*b->max_size_proj), /* dvd_improvex_get_eigenvectors */
151:       (herm || !d->eps->trueres)?0:b->max_nev*b->max_nev+PetscMax(b->max_nev*6,(b->max_nev+b->max_size_proj)*s+b->max_nev*(b->max_size_X+b->max_size_cX_proj)*(std_probl?2:4)+64))); /* preTestConv */
152:   b->own_vecs+= size_P; /* KZ */

154:   /* Setup the preconditioner */
155:   if (ksp) {
156:     KSPGetPC(ksp, &pc);
157:     dvd_static_precond_PC(d, b, pc);
158:   } else {
159:     dvd_static_precond_PC(d, b, 0);
160:   }

162:   /* Setup the step */
163:   if (b->state >= DVD_STATE_CONF) {
164:     PetscMalloc(sizeof(dvdImprovex_jd), &data);
165:     data->dynamic = dynamic;
166:     data->size_real_KZ = size_P;
167:     data->real_KZ = b->free_vecs; b->free_vecs+= data->size_real_KZ;
168:     d->max_cX_in_impr = cX_impr;
169:     data->XKZ = b->free_scalars; b->free_scalars+= size_P*size_P;
170:     data->ldXKZ = size_P;
171:     data->size_X = b->max_size_X;
172:     data->old_improveX_data = d->improveX_data;
173:     d->improveX_data = data;
174:     data->old_improveX = d->improveX;
175:     data->ksp = useGD?PETSC_NULL:ksp;
176:     data->d = d;
177:     d->improveX = dvd_improvex_jd_gen;
178:     data->ksp_max_size = max_bs;

180:     DVD_FL_ADD(d->startList, dvd_improvex_jd_start);
181:     DVD_FL_ADD(d->endList, dvd_improvex_jd_end);
182:     DVD_FL_ADD(d->destroyList, dvd_improvex_jd_d);
183:   }

185:   return(0);
186: }


191: PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
192: {
193:   PetscErrorCode  ierr;
194:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
195:   PetscInt        rA, cA, rlA, clA;
196:   Mat             A;
197:   PetscBool       t;
198:   PC              pc;


202:   data->KZ = data->real_KZ;
203:   data->size_KZ = data->size_cX = data->old_size_X = 0;
204:   data->lastTol = data->dynamic?0.5:0.0;

206:   /* Setup the ksp */
207:   if(data->ksp) {
208:     /* Create the reference vector */
209:     VecCreateCompWithVecs(d->V, data->ksp_max_size, PETSC_NULL,
210:                                  &data->friends);

212:     /* Save the current pc and set a PCNONE */
213:     KSPGetPC(data->ksp, &data->old_pc);
214:     PetscTypeCompare((PetscObject)data->old_pc, PCNONE, &t);
215: 
216:     data->lastTol = 0.5;
217:     if (t) {
218:       data->old_pc = 0;
219:     } else {
220:       PetscObjectReference((PetscObject)data->old_pc);
221:       PCCreate(((PetscObject)d->eps)->comm, &pc);
222:       PCSetType(pc, PCNONE);
223:       PCSetOperators(pc, d->A, d->A, SAME_PRECONDITIONER);
224:       KSPSetPC(data->ksp, pc);
225:       PCDestroy(&pc);
226:     }

228:     /* Create the (I-v*u')*K*(A-s*B) matrix */
229:     MatGetSize(d->A, &rA, &cA);
230:     MatGetLocalSize(d->A, &rlA, &clA);
231:     MatCreateShell(((PetscObject)d->A)->comm, rlA*data->ksp_max_size,
232:                           clA*data->ksp_max_size, rA*data->ksp_max_size,
233:                           cA*data->ksp_max_size, data, &A);
234:     MatShellSetOperation(A, MATOP_MULT,
235:                                 (void(*)(void))dvd_matmult_jd);
236:     MatShellSetOperation(A, MATOP_GET_VECS,
237:                                 (void(*)(void))dvd_matgetvecs_jd);
238: 

240:     /* Try to avoid KSPReset */
241:     KSPGetOperatorsSet(data->ksp,&t,PETSC_NULL);
242:     if (t) {
243:       Mat M;
244:       PetscInt rM;
245:       KSPGetOperators(data->ksp,&M,PETSC_NULL,PETSC_NULL);
246:       MatGetSize(M,&rM,PETSC_NULL);
247:       if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
248:     }
249:     KSPSetOperators(data->ksp, A, A, SAME_PRECONDITIONER);
250: 
251:     KSPSetUp(data->ksp);
252:     MatDestroy(&A);
253:   } else {
254:     data->old_pc = 0;
255:     data->friends = PETSC_NULL;
256:   }
257: 
258:   return(0);
259: }


264: PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
265: {
266:   PetscErrorCode  ierr;
267:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;


271:   if (data->friends) { VecDestroy(&data->friends);  }

273:   /* Restore the pc of ksp */
274:   if (data->old_pc) {
275:     KSPSetPC(data->ksp, data->old_pc);
276:     PCDestroy(&data->old_pc);
277:   }

279:   return(0);
280: }


285: PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
286: {
287:   PetscErrorCode  ierr;
288:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;

291: 
292:   /* Restore changes in dvdDashboard */
293:   d->improveX_data = data->old_improveX_data;

295:   /* Free local data and objects */
296:   PetscFree(data);

298:   return(0);
299: }


304: PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d, Vec *D,
305:                              PetscInt max_size_D, PetscInt r_s,
306:                              PetscInt r_e, PetscInt *size_D)
307: {
308:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
309:   PetscErrorCode  ierr;
310:   PetscInt        i, j, n, maxits, maxits0, lits, s, ldpX;
311:   PetscScalar     *pX, *pY, *auxS = d->auxS, *auxS0;
312:   PetscReal       tol, tol0;
313:   Vec             *u, *v, *kr, kr_comp, D_comp;


317:   /* Quick exit */
318:   if ((max_size_D == 0) || r_e-r_s <= 0) {
319:    /* Callback old improveX */
320:     if (data->old_improveX) {
321:       d->improveX_data = data->old_improveX_data;
322:       data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
323:       d->improveX_data = data;
324:     }
325:     return(0);
326:   }
327: 
328:   n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
329:   if (n == 0) SETERRQ(PETSC_COMM_SELF,1, "n == 0!\n");
330:   if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1, "size_X < r_e-r_s!\n");

332:   /* Compute the eigenvectors of the selected pairs */
333:   if (DVD_IS(d->sEP, DVD_EP_HERMITIAN)) {
334:     pX = d->pX;
335:     pY = d->pY?d->pY:d->pX;
336:     ldpX = d->ldpX;
337:   } else {
338:     pX = auxS; auxS+= d->size_H*d->size_H;
339:     pY = auxS; auxS+= d->size_H*d->size_H;
340:     dvd_improvex_get_eigenvectors(d, pX, pY, d->size_H, auxS,
341:                                          d->size_auxS-(auxS-d->auxS));
342: 
343:     ldpX = d->size_H;
344:   }

346:   /* Restart lastTol if a new pair converged */
347:   if (data->dynamic && data->size_cX < d->size_cX)
348:     data->lastTol = 0.5;

350:   for(i=0, s=0, auxS0=auxS; i<n; i+=s) {
351:     /* If the selected eigenvalue is complex, but the arithmetic is real... */
352: #if !defined(PETSC_USE_COMPLEX)
353:     if (PetscAbsScalar(d->eigi[i] != 0.0)) {
354:       if (i+2 <= max_size_D) s=2; else break;
355:     } else
356: #endif
357:       s=1;

359:     data->auxV = d->auxV;
360:     data->r_s = r_s+i; data->r_e = r_s+i+s;
361:     auxS = auxS0;
362:     data->theta = auxS; auxS+= 2*s;
363:     data->thetai = auxS; auxS+= s;
364:     /* If GD, kr = D */
365:     if (!data->ksp) {
366:       kr = &D[i];

368:     /* If JD, kr = auxV */
369:     } else {
370:       kr = data->auxV; data->auxV+= s;
371:     }

373:     /* Compute theta, maximum iterations and tolerance */
374:     maxits = 0; tol = 1;
375:     for(j=0; j<s; j++) {
376:       d->improvex_jd_lit(d, r_s+i+j, &data->theta[2*j],
377:                                 &data->thetai[j], &maxits0, &tol0);
378: 
379:       maxits+= maxits0; tol*= tol0;
380:     }
381:     maxits/= s; tol = data->dynamic?data->lastTol:exp(log(tol)/s);

383:     /* Compute u, v and kr */
384:     dvd_improvex_jd_proj_cuv(d, r_s+i, r_s+i+s, &u, &v, kr,
385:              &data->auxV, &auxS, data->theta, data->thetai,
386:              &pX[d->size_H*(r_s+i+d->cX_in_H)], &pY[d->size_H*(r_s+i+d->cX_in_H)], ldpX);
387: 
388:     data->u = u;

390:     /* Check if the first eigenpairs are converged */
391:     if (i == 0) {
392:       PetscInt n_auxV = data->auxV-d->auxV+s, n_auxS = auxS - d->auxS;
393:       d->auxV+= n_auxV; d->size_auxV-= n_auxV;
394:       d->auxS+= n_auxS; d->size_auxS-= n_auxS;
395:       d->preTestConv(d,0,s,s,d->auxV-s,PETSC_NULL,&d->npreconv);
396:       d->auxV-= n_auxV; d->size_auxV+= n_auxV;
397:       d->auxS-= n_auxS; d->size_auxS+= n_auxS;
398:       if (d->npreconv > 0) break;
399:     }

401:     /* Compute kr <- kr - v*(u'*kr) */
402:     dvd_improvex_apply_proj(d, kr, s, auxS);

404:     /* If JD */
405:     if (data->ksp) {
406:       data->auxS = auxS;

408:       /* kr <- -kr */
409:       for(j=0; j<s; j++) {
410:         VecScale(kr[j], -1.0);
411:       }

413:       /* Compouse kr and D */
414:       VecCreateCompWithVecs(kr, data->ksp_max_size, data->friends,
415:                                    &kr_comp);
416:       VecCreateCompWithVecs(&D[i], data->ksp_max_size, data->friends,
417:                                    &D_comp);
418:       VecCompSetSubVecs(data->friends,s,PETSC_NULL);
419: 
420:       /* Solve the correction equation */
421:       KSPSetTolerances(data->ksp, tol, PETSC_DEFAULT, PETSC_DEFAULT,
422:                               maxits);
423:       KSPSolve(data->ksp, kr_comp, D_comp);
424:       KSPGetIterationNumber(data->ksp, &lits);
425:       d->eps->OP->lineariterations+= lits;
426: 
427:       /* Destroy the composed ks and D */
428:       VecDestroy(&kr_comp);
429:       VecDestroy(&D_comp);
430:     }
431:   }
432:   *size_D = i;
433:   if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
434: 
435:   /* Callback old improveX */
436:   if (data->old_improveX) {
437:     d->improveX_data = data->old_improveX_data;
438:     data->old_improveX(d, PETSC_NULL, 0, 0, 0, PETSC_NULL);
439:     d->improveX_data = data;
440:   }

442:   return(0);
443: }


448: PetscErrorCode dvd_matmult_jd(Mat A, Vec in, Vec out)
449: {
450:   PetscErrorCode  ierr;
451:   dvdImprovex_jd  *data;
452:   PetscInt        n, i;
453:   const Vec       *inx, *outx, *Bx;


457:   MatShellGetContext(A, (void**)&data);
458:   VecCompGetSubVecs(in,PETSC_NULL,&inx);
459:   VecCompGetSubVecs(out,PETSC_NULL,&outx);
460:   n = data->r_e - data->r_s;

462:   /* aux <- theta[1]A*in - theta[0]*B*in */
463:   for(i=0; i<n; i++) {
464:     MatMult(data->d->A, inx[i], data->auxV[i]);
465:   }
466:   if (data->d->B) {
467:     for(i=0; i<n; i++) {
468:       MatMult(data->d->B, inx[i], outx[i]);
469:     }
470:     Bx = outx;
471:   } else
472:     Bx = inx;

474:   for(i=0; i<n; i++) {
475: #if !defined(PETSC_USE_COMPLEX)
476:     if(data->d->eigi[data->r_s+i] != 0.0) {
477:       /* aux_i   <- [ t_2i+1*A*inx_i   - t_2i*Bx_i + ti_i*Bx_i+1;
478:          aux_i+1      t_2i+1*A*inx_i+1 - ti_i*Bx_i - t_2i*Bx_i+1  ] */
479:       VecAXPBYPCZ(data->auxV[i], -data->theta[2*i], data->thetai[i],
480:                          data->theta[2*i+1], Bx[i], Bx[i+1]);
481: 
482:       VecAXPBYPCZ(data->auxV[i+1], -data->thetai[i],
483:                          -data->theta[2*i], data->theta[2*i+1], Bx[i],
484:                          Bx[i+1]);
485:       i++;
486:     } else
487: #endif
488:     {
489:       VecAXPBY(data->auxV[i], -data->theta[i*2], data->theta[i*2+1],
490:                       Bx[i]);
491:     }
492:   }

494:   /* out <- K * aux */
495:   for(i=0; i<n; i++) {
496:     data->d->improvex_precond(data->d, data->r_s+i, data->auxV[i],
497:                                      outx[i]);
498:   }

500:   /* out <- out - v*(u'*out) */
501:   dvd_improvex_apply_proj(data->d, (Vec*)outx, n, data->auxS);

503:   return(0);
504: }


509: PetscErrorCode dvd_matgetvecs_jd(Mat A, Vec *right, Vec *left)
510: {
511:   PetscErrorCode  ierr;
512:   Vec             *r, *l;
513:   dvdImprovex_jd  *data;
514:   PetscInt        n, i;


518:   MatShellGetContext(A, (void**)&data);
519:   n = data->ksp_max_size;
520:   if (right) {
521:     PetscMalloc(sizeof(Vec)*n, &r);
522:   }
523:   if (left) {
524:     PetscMalloc(sizeof(Vec)*n, &l);
525:   }
526:   for (i=0; i<n; i++) {
527:     MatGetVecs(data->d->A, right?&r[i]:PETSC_NULL,
528:                       left?&l[i]:PETSC_NULL);
529:   }
530:   if(right) {
531:     VecCreateCompWithVecs(r, n, data->friends, right);
532:     for (i=0; i<n; i++) {
533:       VecDestroy(&r[i]);
534:     }
535:   }
536:   if(left) {
537:     VecCreateCompWithVecs(l, n, data->friends, left);
538:     for (i=0; i<n; i++) {
539:       VecDestroy(&l[i]);
540:     }
541:   }

543:   if (right) {
544:     PetscFree(r);
545:   }
546:   if (left) {
547:     PetscFree(l);
548:   }

550:   return(0);
551: }


556: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d, dvdBlackboard *b,
557:                                        ProjType_t p)
558: {

561:   /* Setup the step */
562:   if (b->state >= DVD_STATE_CONF) {
563:     switch(p) {
564:     case DVD_PROJ_KXX:
565:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KXX; break;
566:     case DVD_PROJ_KZX:
567:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX; break;
568:     }
569:   }

571:   return(0);
572: }

576: /* 
577:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
578:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
579:   where
580:   auxV, 4*(i_e-i_s) auxiliar global vectors
581:   pX,pY, the right and left eigenvectors of the projected system
582:   ld, the leading dimension of pX and pY
583:   auxS: max(8*bs*max_cX_in_proj,size_V*size_V)
584: */
585: PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d, PetscInt i_s,
586:   PetscInt i_e, Vec **u, Vec **v, Vec *kr, Vec **auxV, PetscScalar **auxS,
587:   PetscScalar *theta, PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY,
588:   PetscInt ld)
589: {
590:   PetscErrorCode  ierr;
591:   PetscInt        n = i_e - i_s, size_KZ, V_new, rm, i, size_in;
592:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
593:   PetscBLASInt    s, ldXKZ, info;
594:   DvdReduction    r;
595:   DvdReductionChunk
596:                   ops[2];
597:   DvdMult_copy_func
598:                   sr[2];


602:   /* Check consistency */
603:   V_new = d->size_cX - data->size_cX;
604:   if (V_new > data->old_size_X) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }
605:   data->old_size_X = n;

607:   /* KZ <- KZ(rm:rm+max_cX-1) */
608:   rm = PetscMax(V_new+data->size_KZ-d->max_cX_in_impr, 0);
609:   for (i=0; i<d->max_cX_in_impr; i++) {
610:     VecCopy(data->KZ[i+rm], data->KZ[i]);
611:   }
612:   data->size_cX = d->size_cX;

614:   /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
615:   for (i=0; i<d->max_cX_in_impr; i++) {
616:     SlepcDenseCopy(&data->XKZ[i*data->ldXKZ+i], data->ldXKZ, &data->XKZ[(i+rm)*data->ldXKZ+i+rm], data->ldXKZ, data->size_KZ, 1);
617:   }
618:   data->size_KZ = PetscMin(d->max_cX_in_impr, data->size_KZ+V_new);

620:   /* Compute X, KZ and KR */
621:   *u = *auxV; *auxV+= n;
622:   *v = &data->KZ[data->size_KZ];
623:   d->improvex_jd_proj_uv(d, i_s, i_e, *u, *v, kr, *auxV, theta, thetai,
624:                                 pX, pY, ld);

626:   /* XKZ <- X'*KZ */
627:   size_KZ = data->size_KZ+n;
628:   size_in = 2*n*data->size_KZ+n*n;
629:   SlepcAllReduceSumBegin(ops,2,*auxS,*auxS+size_in,size_in,&r,((PetscObject)d->V[0])->comm);
630:   VecsMultS(data->XKZ,0,data->ldXKZ,d->V-data->size_KZ,0,data->size_KZ,data->KZ,data->size_KZ,size_KZ,&r,&sr[0]);
631:   VecsMultS(&data->XKZ[data->size_KZ],0,data->ldXKZ,*u,0,n,data->KZ,0,size_KZ,&r,&sr[1]);
632:   SlepcAllReduceSumEnd(&r);

634:   /* iXKZ <- inv(XKZ) */
635:   s = PetscBLASIntCast(size_KZ);
636:   data->ldiXKZ = data->size_iXKZ = size_KZ;
637:   data->iXKZ = *auxS; *auxS+= size_KZ*size_KZ;
638:   data->iXKZPivots = (PetscBLASInt*)*auxS;
639:   *auxS += FromIntToScalar(size_KZ);
640:   SlepcDenseCopy(data->iXKZ,data->ldiXKZ,data->XKZ,data->ldXKZ,size_KZ,size_KZ);
641:   ldXKZ = PetscBLASIntCast(data->ldiXKZ);
642:   LAPACKgetrf_(&s, &s, data->iXKZ, &ldXKZ, data->iXKZPivots, &info);
643:   if (info) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_LIB, "Error in Lapack XGETRF %d", info);

645:   return(0);
646: }

648: #define DVD_COMPLEX_RAYLEIGH_QUOTIENT(ur,ui,Axr,Axi,Bxr,Bxi,eigr,eigi,b,ierr)\
649: { \
650:   VecDot((Axr), (ur), &(b)[0]);  /* r*A*r */ \
651:   VecDot((Axr), (ui), &(b)[1]);  /* i*A*r */ \
652:   VecDot((Axi), (ur), &(b)[2]);  /* r*A*i */ \
653:   VecDot((Axi), (ui), &(b)[3]);  /* i*A*i */ \
654:   VecDot((Bxr), (ur), &(b)[4]);  /* r*B*r */ \
655:   VecDot((Bxr), (ui), &(b)[5]);  /* i*B*r */ \
656:   VecDot((Bxi), (ur), &(b)[6]);  /* r*B*i */ \
657:   VecDot((Bxi), (ui), &(b)[7]);  /* i*B*i */ \
658:   (b)[0]  = (b)[0]+(b)[3]; /* rAr+iAi */ \
659:   (b)[2] =  (b)[2]-(b)[1]; /* rAi-iAr */ \
660:   (b)[4] = (b)[4]+(b)[7]; /* rBr+iBi */ \
661:   (b)[6] = (b)[6]-(b)[5]; /* rBi-iBr */ \
662:   (b)[7] = (b)[4]*(b)[4] + (b)[6]*(b)[6]; /* k */ \
663:   *(eigr) = ((b)[0]*(b)[4] + (b)[2]*(b)[6]) / (b)[7]; /* eig_r */ \
664:   *(eigi) = ((b)[2]*(b)[4] - (b)[0]*(b)[6]) / (b)[7]; /* eig_i */ \
665: }

667: #if !defined(PETSC_USE_COMPLEX)
668: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
669:   for((i)=0; (i)<(n); (i)++) { \
670:     if ((eigi)[(i_s)+(i)] != 0.0) { \
671:       /* eig_r = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k \
672:          eig_i = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k \
673:          k     =  (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr)    */ \
674:       DVD_COMPLEX_RAYLEIGH_QUOTIENT((u)[(i)], (u)[(i)+1], (Ax)[(i)], \
675:         (Ax)[(i)+1], (Bx)[(i)], (Bx)[(i)+1], &(b)[8], &(b)[9], (b), (ierr)); \
676:       if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[8])/ \
677:             PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10    || \
678:           PetscAbsScalar((eigi)[(i_s)+(i)] - (b)[9])/ \
679:             PetscAbsScalar((eigi)[(i_s)+(i)]) > 1e-10         ) { \
680:         (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its "\
681:                             "Rayleigh quotient value %G+%G\n", \
682:                             (eigr)[(i_s)+(i)], \
683:                             (eigi)[(i_s)+1], (b)[8], (b)[9]); \
684:       } \
685:       (i)++; \
686:     } \
687:   }
688: #else
689: #define DVD_COMPUTE_N_RR(eps,i,i_s,n,eigr,eigi,u,Ax,Bx,b,ierr) \
690:   for((i)=0; (i)<(n); (i)++) { \
691:       (ierr) = VecDot((Ax)[(i)], (u)[(i)], &(b)[0]);  \
692:       (ierr) = VecDot((Bx)[(i)], (u)[(i)], &(b)[1]);  \
693:       (b)[0] = (b)[0]/(b)[1]; \
694:       if (PetscAbsScalar((eigr)[(i_s)+(i)] - (b)[0])/ \
695:             PetscAbsScalar((eigr)[(i_s)+(i)]) > 1e-10     ) { \
696:         (ierr) = PetscInfo4((eps), "The eigenvalue %G+%G is far from its " \
697:                "Rayleigh quotient value %G+%G\n", \
698:                PetscRealPart((eigr)[(i_s)+(i)]), \
699:                PetscImaginaryPart((eigr)[(i_s)+(i)]), PetscRealPart((b)[0]), \
700:                PetscImaginaryPart((b)[0])); \
701:       } \
702:     }
703: #endif

707: PETSC_STATIC_INLINE PetscErrorCode dvd_improvex_compute_X(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,PetscScalar *pX,PetscInt ld)
708: {
709:   PetscErrorCode  ierr;
710:   PetscInt        n = i_e - i_s, i;

713:   SlepcUpdateVectorsZ(u, 0.0, 1.0, d->V-d->cX_in_H, d->size_V+d->cX_in_H, pX, ld, d->size_H, n);
714:   /* nX(i) <- ||X(i)|| */
715:   if (d->correctXnorm) {
716:     for (i=0; i<n; i++) {
717:       VecNormBegin(u[i], NORM_2, &d->nX[i_s+i]);
718:     }
719:     for (i=0; i<n; i++) {
720:       VecNormEnd(u[i], NORM_2, &d->nX[i_s+i]);
721:     }
722: #if !defined(PETSC_USE_COMPLEX)
723:     for(i=0; i<n; i++)
724:       if(d->eigi[i_s+i] != 0.0)
725:         d->nX[i_s+i] = d->nX[i_s+i+1] = PetscSqrtScalar(d->nX[i_s+i]*d->nX[i_s+i]+d->nX[i_s+i+1]*d->nX[i_s+i+1]);
726: #endif
727:   } else {
728:     for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;
729:   }
730:   return(0);
731: }

735: /* 
736:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
737:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
738:   where
739:   auxV, 4*(i_e-i_s) auxiliar global vectors
740:   pX,pY, the right and left eigenvectors of the projected system
741:   ld, the leading dimension of pX and pY
742: */
743: PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d, PetscInt i_s,
744:   PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
745:   PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
746: {
747:   PetscErrorCode  ierr;
748:   PetscInt        n = i_e - i_s, i;
749:   PetscScalar     b[16];
750:   Vec             *Ax, *Bx, *r=auxV, X[4];
751:   /* The memory manager doen't allow to call a subroutines */
752:   PetscScalar     Z[size_Z];


756:   /* u <- X(i) */
757:   dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);

759:   /* v <- theta[0]A*u + theta[1]*B*u */

761:   /* Bx <- B*X(i) */
762:   Bx = kr;
763:   if (d->BV) {
764:     SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV-d->cX_in_H, d->size_BV+d->cX_in_H, pX, ld, d->size_H, n);
765:   } else {
766:     for(i=0; i<n; i++) {
767:       if (d->B) {
768:         MatMult(d->B, u[i], Bx[i]);
769:       } else {
770:         VecCopy(u[i], Bx[i]);
771:       }
772:     }
773:   }

775:   /* Ax <- A*X(i) */
776:   Ax = r;
777:   SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV-d->cX_in_H, d->size_AV+d->cX_in_H, pX, ld, d->size_H, n);

779:   /* v <- Y(i) */
780:   SlepcUpdateVectorsZ(v, 0.0, 1.0, (d->W?d->W:d->V)-d->cX_in_H, d->size_V+d->cX_in_H, pY, ld, d->size_H, n);

782:   /* Recompute the eigenvalue */
783:   DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, v, Ax, Bx, b, ierr);

785:   for(i=0; i<n; i++) {
786: #if !defined(PETSC_USE_COMPLEX)
787:     if(d->eigi[i_s+i] != 0.0) {
788:       /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i'    0            1        0   
789:                                        0         theta_2i'     0        1   
790:                                      theta_2i+1 -thetai_i   -eigr_i -eigi_i 
791:                                      thetai_i    theta_2i+1  eigi_i -eigr_i ] */
792:       b[0] = b[5] = PetscConj(theta[2*i]);
793:       b[2] = b[7] = -theta[2*i+1];
794:       b[6] = -(b[3] = thetai[i]);
795:       b[1] = b[4] = 0.0;
796:       b[8] = b[13] = 1.0/d->nX[i_s+i];
797:       b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
798:       b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
799:       b[9] = b[12] = 0.0;
800:       X[0] = Ax[i]; X[1] = Ax[i+1]; X[2] = Bx[i]; X[3] = Bx[i+1];
801:       SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 4, Z, size_Z);
802: 
803:       i++;
804:     } else
805: #endif
806:     {
807:       /* [Ax_i Bx_i]*= [ theta_2i'    1/nX_i
808:                         theta_2i+1  -eig_i/nX_i ] */
809:       b[0] = PetscConj(theta[i*2]);
810:       b[1] = theta[i*2+1];
811:       b[2] = 1.0/d->nX[i_s+i];
812:       b[3] = -d->eigr[i_s+i]/d->nX[i_s+i];
813:       X[0] = Ax[i]; X[1] = Bx[i];
814:       SlepcUpdateVectorsD(X, 2, 1.0, b, 2, 2, 2, Z, size_Z);
815: 
816:     }
817:   }
818:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

820:   /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
821:   for(i=0; i<n; i++) {
822:     d->improvex_precond(d, i_s+i, r[i], v[i]);
823:   }

825:   /* kr <- K^{-1}*kr = K^{-1}*(Ax - eig_i*Bx) */
826:   d->calcpairs_proj_res(d, i_s, i_e, Bx);
827:   for(i=0; i<n; i++) {
828:     VecCopy(Bx[i], r[i]);
829:     d->improvex_precond(d, i_s+i, r[i], kr[i]);
830:   }

832:   return(0);
833: }

837: /* 
838:   Compute: u <- K^{-1}*X, v <- X,
839:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
840:   where
841:   auxV, 4*(i_e-i_s) auxiliar global vectors
842:   pX,pY, the right and left eigenvectors of the projected system
843:   ld, the leading dimension of pX and pY
844: */
845: PetscErrorCode dvd_improvex_jd_proj_uv_KXX(dvdDashboard *d, PetscInt i_s,
846:   PetscInt i_e, Vec *u, Vec *v, Vec *kr, Vec *auxV, PetscScalar *theta,
847:   PetscScalar *thetai, PetscScalar *pX, PetscScalar *pY, PetscInt ld)
848: {
849:   PetscErrorCode  ierr;
850:   PetscInt        n = i_e - i_s, i;
851:   PetscScalar     b[16];
852:   Vec             *Ax, *Bx, *r = auxV, X[4];
853:   /* The memory manager doen't allow to call a subroutines */
854:   PetscScalar     Z[size_Z];


858:   /* [v u] <- X(i) Y(i) */
859:   dvd_improvex_compute_X(d,i_s,i_e,v,pX,ld);
860:   SlepcUpdateVectorsZ(u, 0.0, 1.0, (d->W?d->W:d->V)-d->cX_in_H, d->size_V+d->cX_in_H, pY, ld, d->size_H, n);

862:   /* Bx <- B*X(i) */
863:   Bx = kr;
864:   if (d->BV) {
865:     SlepcUpdateVectorsZ(Bx, 0.0, 1.0, d->BV-d->cX_in_H, d->size_BV+d->cX_in_H, pX, ld, d->size_H, n);
866:   } else {
867:     if (d->B) {
868:       for(i=0; i<n; i++) {
869:         MatMult(d->B, v[i], Bx[i]);
870:       }
871:     } else
872:       Bx = v;
873:   }

875:   /* Ax <- A*X(i) */
876:   Ax = r;
877:   SlepcUpdateVectorsZ(Ax, 0.0, 1.0, d->AV-d->cX_in_H, d->size_AV+d->cX_in_H, pX, ld, d->size_H, n);

879:   /* Recompute the eigenvalue */
880:   DVD_COMPUTE_N_RR(d->eps, i, i_s, n, d->eigr, d->eigi, u, Ax, Bx, b, ierr);

882:   for(i=0; i<n; i++) {
883:     if (d->eigi[i_s+i] == 0.0) {
884:       /* r <- Ax -eig*Bx */
885:       VecAXPBY(r[i], -d->eigr[i_s+i]/d->nX[i_s+i], 1.0/d->nX[i_s+i], Bx[i]);
886:     } else {
887:       /* [r_i r_i+1 kr_i kr_i+1]*= [   1        0 
888:                                        0        1
889:                                     -eigr_i -eigi_i
890:                                      eigi_i -eigr_i] */
891:       b[0] = b[5] = 1.0/d->nX[i_s+i];
892:       b[2] = b[7] = -d->eigr[i_s+i]/d->nX[i_s+i];
893:       b[6] = -(b[3] = d->eigi[i_s+i]/d->nX[i_s+i]);
894:       b[1] = b[4] = 0.0;
895:       X[0] = r[i]; X[1] = r[i+1]; X[2] = kr[i]; X[3] = kr[i+1];
896:       SlepcUpdateVectorsD(X, 4, 1.0, b, 4, 4, 2, Z, size_Z);
897: 
898:       i++;
899:     }
900:   }
901:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

903:   /* kr <- K^{-1}*r */
904:   d->calcpairs_proj_res(d, i_s, i_e, r);
905:   for(i=0; i<n; i++) {
906:     d->improvex_precond(d, i_s+i, r[i], kr[i]);
907:   }

909:   /* u <- K^{-1} X(i) */
910:   for(i=0; i<n; i++) {
911:     d->improvex_precond(d, i_s+i, v[i], u[i]);
912:   }

914:   return(0);
915: }


920: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d, dvdBlackboard *b,
921:                                          PetscInt maxits, PetscReal tol,
922:                                          PetscReal fix)
923: {
924:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;


928:   /* Setup the step */
929:   if (b->state >= DVD_STATE_CONF) {
930:     data->maxits = maxits;
931:     data->tol = tol;
932:     data->fix = fix;
933:     d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
934:   }

936:   return(0);
937: }


942: PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d, PetscInt i,
943:   PetscScalar* theta, PetscScalar* thetai, PetscInt *maxits, PetscReal *tol)
944: {
945:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
946:   PetscReal       a;

949:   a = SlepcAbsEigenvalue(d->eigr[i],d->eigi[i]);

951:   if (d->nR[i]/a < data->fix) {
952:     theta[0] = d->eigr[i];
953:     theta[1] = 1.0;
954: #if !defined(PETSC_USE_COMPLEX)
955:     *thetai = d->eigi[i];
956: #endif
957:   } else {
958:     theta[0] = d->target[0];
959:     theta[1] = d->target[1];
960: #if !defined(PETSC_USE_COMPLEX)
961:     *thetai = 0.0;
962: #endif
963: }

965: #if defined(PETSC_USE_COMPLEX)
966:   if(thetai) *thetai = 0.0;
967: #endif
968:   *maxits = data->maxits;
969:   *tol = data->tol;

971:   return(0);
972: }


975: /**** Patterns implementation *************************************************/

977: typedef PetscInt (*funcV0_t)(dvdDashboard*, PetscInt, PetscInt, Vec*);
978: typedef PetscInt (*funcV1_t)(dvdDashboard*, PetscInt, PetscInt, Vec*,
979:                              PetscScalar*, Vec);

983: /* Compute D <- K^{-1} * funcV[r_s..r_e] */
984: PetscErrorCode dvd_improvex_PfuncV(dvdDashboard *d, void *funcV, Vec *D,
985:   PetscInt max_size_D, PetscInt r_s, PetscInt r_e, Vec *auxV,
986:   PetscScalar *auxS)
987: {
988:   PetscErrorCode  ierr;
989:   PetscInt        i;


993:   if (max_size_D >= r_e-r_s+1) {
994:     /* The optimized version needs one vector extra of D */
995:     /* D(1:r.size) = R(r_s:r_e-1) */
996:     if (auxS) ((funcV1_t)funcV)(d, r_s, r_e, D+1, auxS, auxV[0]);
997:     else      ((funcV0_t)funcV)(d, r_s, r_e, D+1);

999:     /* D = K^{-1} * R */
1000:     for (i=0; i<r_e-r_s; i++) {
1001:       d->improvex_precond(d, i+r_s, D[i+1], D[i]);
1002:     }
1003:   } else if (max_size_D == r_e-r_s) {
1004:     /* Non-optimized version */
1005:     /* auxV <- R[r_e-1] */
1006:     if (auxS) ((funcV1_t)funcV)(d, r_e-1, r_e, auxV, auxS, auxV[1]);
1007:     else      ((funcV0_t)funcV)(d, r_e-1, r_e, auxV);

1009:     /* D(1:r.size-1) = R(r_s:r_e-2) */
1010:     if (auxS) ((funcV1_t)funcV)(d, r_s, r_e-1, D+1, auxS, auxV[1]);
1011:     else      ((funcV0_t)funcV)(d, r_s, r_e-1, D+1);

1013:     /* D = K^{-1} * R */
1014:     for (i=0; i<r_e-r_s-1; i++) {
1015:       d->improvex_precond(d, i+r_s, D[i+1], D[i]);
1016:     }
1017:     d->improvex_precond(d, r_e-1, auxV[0], D[r_e-r_s-1]);
1018:   } else {
1019:     SETERRQ(PETSC_COMM_SELF,1, "Problem: r_e-r_s > max_size_D!");
1020:   }

1022:   return(0);
1023: }


1028: /* Compute the left and right projected eigenvectors where,
1029:    pX, the returned right eigenvectors
1030:    pY, the returned left eigenvectors,
1031:    ld_, the leading dimension of pX and pY,
1032:    auxS, auxiliar vector of size length 6*d->size_H
1033: */
1034: PetscErrorCode dvd_improvex_get_eigenvectors(dvdDashboard *d, PetscScalar *pX,
1035:   PetscScalar *pY, PetscInt ld, PetscScalar *auxS, PetscInt size_auxS)
1036: {
1037:   PetscErrorCode  ierr;
1038: 

1041:   SlepcDenseCopy(pY, ld, d->T?d->pY:d->pX, d->ldpX, d->size_H,
1042:                         d->size_H);
1043:   SlepcDenseCopy(pX, ld, d->pX, d->ldpX, d->size_H, d->size_H);
1044: 
1045: 
1046:   /* [qX, qY] <- eig(S, T); pX <- d->pX * qX; pY <- d->pY * qY */
1047:   dvd_compute_eigenvectors(d->size_H, d->S, d->ldS, d->T, d->ldT, pX,
1048:                                   ld, pY, ld, auxS, size_auxS, PETSC_TRUE);
1049: 

1051:   /* 2-Normalize the columns of pX an pY */
1052:   SlepcDenseNorm(pX, ld, d->size_H, d->size_H, d->eigi-d->cX_in_H);
1053:   SlepcDenseNorm(pY, ld, d->size_H, d->size_H, d->eigi-d->cX_in_H);

1055:   return(0);
1056: }

1060: /* Compute (I - KZ*iXKZ*X')*V where,
1061:    V, the vectors to apply the projector,
1062:    cV, the number of vectors in V,
1063:    auxS, auxiliar vector of size length 3*size_iXKZ*cV
1064: */
1065: PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d, Vec *V, PetscInt cV, PetscScalar *auxS)
1066: {
1067:   PetscErrorCode  ierr;
1068:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
1069:   PetscInt        size_in = data->size_iXKZ*cV, i, ldh;
1070:   PetscScalar     *h, *in, *out;
1071:   PetscBLASInt    cV_, n, info, ld;
1072:   DvdReduction    r;
1073:   DvdReductionChunk
1074:                   ops[4];
1075:   DvdMult_copy_func
1076:                   sr[4];
1077: 

1080:   /* Check consistency */
1081:   if (cV > 2) { SETERRQ(PETSC_COMM_SELF,1, "Consistency broken!"); }

1083:   /* h <- X'*V */
1084:   h = auxS; in = h+size_in; out = in+size_in; ldh = data->size_iXKZ;
1085:   SlepcAllReduceSumBegin(ops, 4, in, out, size_in, &r,
1086:                                 ((PetscObject)d->V[0])->comm);
1087:   for (i=0; i<cV; i++) {
1088:     VecsMultS(&h[i*ldh],0,ldh,d->V-data->size_KZ,0,data->size_KZ,V+i,0,1,&r,&sr[i*2]);
1089:     VecsMultS(&h[i*ldh+data->size_KZ],0,ldh,data->u,0,data->size_iXKZ-data->size_KZ,V+i,0,1,&r,&sr[i*2+1]);
1090:   }
1091:   SlepcAllReduceSumEnd(&r);

1093:   /* h <- iXKZ\h */
1094:   cV_ = PetscBLASIntCast(cV);
1095:   n = PetscBLASIntCast(data->size_iXKZ);
1096:   ld = PetscBLASIntCast(data->ldiXKZ);
1098:   LAPACKgetrs_("N", &n, &cV_, data->iXKZ, &ld, data->iXKZPivots, h, &n, &info);
1099:   if (info) SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_LIB, "Error in Lapack XGETRS %d", info);

1101:   /* V <- V - KZ*h */
1102:   for (i=0; i<cV; i++) {
1103:     SlepcUpdateVectorsZ(V+i,1.0,-1.0,data->KZ,data->size_iXKZ,&h[ldh*i],ldh,data->size_iXKZ,1);
1104:   }

1106:   return(0);
1107: }