Actual source code: dvdimprovex.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Step: improve the eigenvectors X

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

 10:    This file is part of SLEPc.

 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 <slepcblaslapack.h>

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

 31: typedef struct {
 32:   PetscInt     size_X;
 33:   KSP          ksp;                /* correction equation solver */
 34:   Vec          friends;            /* reference vector for composite vectors */
 35:   PetscScalar  theta[4],thetai[2]; /* the shifts used in the correction eq. */
 36:   PetscInt     maxits;             /* maximum number of iterations */
 37:   PetscInt     r_s,r_e;            /* the selected eigenpairs to improve */
 38:   PetscInt     ksp_max_size;       /* the ksp maximum subvectors size */
 39:   PetscReal    tol;                /* the maximum solution tolerance */
 40:   PetscReal    lastTol;            /* last tol for dynamic stopping criterion */
 41:   PetscReal    fix;                /* tolerance for using the approx. eigenvalue */
 42:   PetscBool    dynamic;            /* if the dynamic stopping criterion is applied */
 43:   dvdDashboard *d;                 /* the currect dvdDashboard reference */
 44:   PC           old_pc;             /* old pc in ksp */
 45:   BV           KZ;                 /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
 46:   BV           U;                  /* new X vectors */
 47:   PetscScalar  *XKZ;               /* X'*KZ */
 48:   PetscScalar  *iXKZ;              /* inverse of XKZ */
 49:   PetscInt     ldXKZ;              /* leading dimension of XKZ */
 50:   PetscInt     size_iXKZ;          /* size of iXKZ */
 51:   PetscInt     ldiXKZ;             /* leading dimension of iXKZ */
 52:   PetscInt     size_cX;            /* last value of d->size_cX */
 53:   PetscInt     old_size_X;         /* last number of improved vectors */
 54:   PetscBLASInt *iXKZPivots;        /* array of pivots */
 55: } dvdImprovex_jd;

 59: /*
 60:    Compute (I - KZ*iXKZ*X')*V where,
 61:    V, the vectors to apply the projector,
 62:    cV, the number of vectors in V,
 63: */
 64: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV)
 65: {
 66: #if defined(PETSC_MISSING_LAPACK_GETRS)
 68:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
 69: #else
 71:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
 72:   PetscInt       i,ldh,k,l;
 73:   PetscScalar    *h;
 74:   PetscBLASInt   cV_,n,info,ld;
 75: #if defined(PETSC_USE_COMPLEX)
 76:   PetscInt       j;
 77: #endif

 80:   if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");

 82:   /* h <- X'*V */
 83:   PetscMalloc1(data->size_iXKZ*cV,&h);
 84:   ldh = data->size_iXKZ;
 85:   BVGetActiveColumns(data->U,&l,&k);
 86:   if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
 87:   BVSetActiveColumns(data->U,0,k);
 88:   for (i=0;i<cV;i++) {
 89:     BVDotVec(data->U,V[i],&h[ldh*i]);
 90: #if defined(PETSC_USE_COMPLEX)
 91:     for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
 92: #endif
 93:   }
 94:   BVSetActiveColumns(data->U,l,k);

 96:   /* h <- iXKZ\h */
 97:   PetscBLASIntCast(cV,&cV_);
 98:   PetscBLASIntCast(data->size_iXKZ,&n);
 99:   PetscBLASIntCast(data->ldiXKZ,&ld);
101:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
102:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
103:   PetscFPTrapPop();
104:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRS %d",info);

106:   /* V <- V - KZ*h */
107:   BVSetActiveColumns(data->KZ,0,k);
108:   for (i=0;i<cV;i++) {
109:     BVMultVec(data->KZ,-1.0,1.0,V[i],&h[ldh*i]);
110:   }
111:   BVSetActiveColumns(data->KZ,l,k);
112:   PetscFree(h);
113:   return(0);
114: #endif
115: }

119: /*
120:    Compute (I - X*iXKZ*KZ')*V where,
121:    V, the vectors to apply the projector,
122:    cV, the number of vectors in V,
123: */
124: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV)
125: {
126: #if defined(PETSC_MISSING_LAPACK_GETRS)
128:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
129:   return(0);
130: #else
132:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
133:   PetscInt       i,ldh,k,l;
134:   PetscScalar    *h;
135:   PetscBLASInt   cV_, n, info, ld;
136: #if defined(PETSC_USE_COMPLEX)
137:   PetscInt       j;
138: #endif

141:   if (cV > 2) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");

143:   /* h <- KZ'*V */
144:   PetscMalloc1(data->size_iXKZ*cV,&h);
145:   ldh = data->size_iXKZ;
146:   BVGetActiveColumns(data->U,&l,&k);
147:   if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1,"Consistency broken");
148:   BVSetActiveColumns(data->KZ,0,k);
149:   for (i=0;i<cV;i++) {
150:     BVDotVec(data->KZ,V[i],&h[ldh*i]);
151: #if defined(PETSC_USE_COMPLEX)
152:     for (j=0;j<k;j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
153: #endif
154:   }
155:   BVSetActiveColumns(data->KZ,l,k);

157:   /* h <- iXKZ\h */
158:   PetscBLASIntCast(cV,&cV_);
159:   PetscBLASIntCast(data->size_iXKZ,&n);
160:   PetscBLASIntCast(data->ldiXKZ,&ld);
162:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
163:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C",&n,&cV_,data->iXKZ,&ld,data->iXKZPivots,h,&n,&info));
164:   PetscFPTrapPop();
165:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRS %d",info);

167:   /* V <- V - U*h */
168:   BVSetActiveColumns(data->U,0,k);
169:   for (i=0;i<cV;i++) {
170:     BVMultVec(data->U,-1.0,1.0,V[i],&h[ldh*i]);
171:   }
172:   BVSetActiveColumns(data->U,l,k);
173:   PetscFree(h);
174:   return(0);
175: #endif
176: }

180: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
181: {
183:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;

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

188:   /* Restore the pc of ksp */
189:   if (data->old_pc) {
190:     KSPSetPC(data->ksp, data->old_pc);
191:     PCDestroy(&data->old_pc);
192:   }
193:   return(0);
194: }

198: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
199: {
201:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;

204:   /* Free local data and objects */
205:   PetscFree(data->XKZ);
206:   PetscFree(data->iXKZ);
207:   PetscFree(data->iXKZPivots);
208:   BVDestroy(&data->KZ);
209:   BVDestroy(&data->U);
210:   PetscFree(data);
211:   return(0);
212: }

216: /*
217:    y <- theta[1]A*x - theta[0]*B*x
218:    auxV, two auxiliary vectors
219:  */
220: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y)
221: {
223:   PetscInt       n,i;
224:   const Vec      *Bx;
225:   Vec            *auxV;

228:   n = data->r_e - data->r_s;
229:   for (i=0;i<n;i++) {
230:     MatMult(data->d->A,x[i],y[i]);
231:   }

233:   SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
234:   for (i=0;i<n;i++) {
235: #if !defined(PETSC_USE_COMPLEX)
236:     if (data->d->eigi[data->r_s+i] != 0.0) {
237:       if (data->d->B) {
238:         MatMult(data->d->B,x[i],auxV[0]);
239:         MatMult(data->d->B,x[i+1],auxV[1]);
240:         Bx = auxV;
241:       } else Bx = &x[i];

243:       /* y_i   <- [ t_2i+1*A*x_i   - t_2i*Bx_i + ti_i*Bx_i+1;
244:          y_i+1      t_2i+1*A*x_i+1 - ti_i*Bx_i - t_2i*Bx_i+1  ] */
245:       VecAXPBYPCZ(y[i],-data->theta[2*i],data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
246:       VecAXPBYPCZ(y[i+1],-data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
247:       i++;
248:     } else
249: #endif
250:     {
251:       if (data->d->B) {
252:         MatMult(data->d->B,x[i],auxV[0]);
253:         Bx = auxV;
254:       } else Bx = &x[i];
255:       VecAXPBY(y[i],-data->theta[i*2],data->theta[i*2+1],Bx[0]);
256:     }
257:   }
258:   SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
259:   return(0);
260: }

264: /*
265:    y <- theta[1]'*A'*x - theta[0]'*B'*x
266:  */
267: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y)
268: {
270:   PetscInt       n,i;
271:   const Vec      *Bx;
272:   Vec            *auxV;

275:   n = data->r_e - data->r_s;
276:   for (i=0;i<n;i++) {
277:     MatMultTranspose(data->d->A,x[i],y[i]);
278:   }

280:   SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
281:   for (i=0;i<n;i++) {
282: #if !defined(PETSC_USE_COMPLEX)
283:     if (data->d->eigi[data->r_s+i] != 0.0) {
284:       if (data->d->B) {
285:         MatMultTranspose(data->d->B,x[i],auxV[0]);
286:         MatMultTranspose(data->d->B,x[i+1],auxV[1]);
287:         Bx = auxV;
288:       } else Bx = &x[i];

290:       /* y_i   <- [ t_2i+1*A*x_i   - t_2i*Bx_i - ti_i*Bx_i+1;
291:          y_i+1      t_2i+1*A*x_i+1 + ti_i*Bx_i - t_2i*Bx_i+1  ] */
292:       VecAXPBYPCZ(y[i],-data->theta[2*i],-data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
293:       VecAXPBYPCZ(y[i+1],data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
294:       i++;
295:     } else
296: #endif
297:     {
298:       if (data->d->B) {
299:         MatMultTranspose(data->d->B,x[i],auxV[0]);
300:         Bx = auxV;
301:       } else Bx = &x[i];
302:       VecAXPBY(y[i],PetscConj(-data->theta[i*2]),PetscConj(data->theta[i*2+1]),Bx[0]);
303:     }
304:   }
305:   SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
306:   return(0);
307: }

311: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w)
312: {
314:   dvdImprovex_jd *data;
315:   PetscInt       n,i;
316:   const Vec      *inx,*outx,*wx;
317:   Vec            *auxV;
318:   Mat            A;

321:   PCGetOperators(pc,&A,NULL);
322:   MatShellGetContext(A,(void**)&data);
323:   VecCompGetSubVecs(in,NULL,&inx);
324:   VecCompGetSubVecs(out,NULL,&outx);
325:   VecCompGetSubVecs(w,NULL,&wx);
326:   n = data->r_e - data->r_s;
327:   SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
328:   switch (side) {
329:   case PC_LEFT:
330:     /* aux <- theta[1]A*in - theta[0]*B*in */
331:     dvd_aux_matmult(data,inx,auxV);

333:     /* out <- K * aux */
334:     for (i=0;i<n;i++) {
335:       data->d->improvex_precond(data->d,data->r_s+i,auxV[i],outx[i]);
336:     }
337:     break;
338:   case PC_RIGHT:
339:     /* aux <- K * in */
340:     for (i=0;i<n;i++) {
341:       data->d->improvex_precond(data->d,data->r_s+i,inx[i],auxV[i]);
342:     }

344:     /* out <- theta[1]A*auxV - theta[0]*B*auxV */
345:     dvd_aux_matmult(data,auxV,outx);
346:     break;
347:   case PC_SYMMETRIC:
348:     /* aux <- K^{1/2} * in */
349:     for (i=0;i<n;i++) {
350:       PCApplySymmetricRight(data->old_pc,inx[i],auxV[i]);
351:     }

353:     /* wx <- theta[1]A*auxV - theta[0]*B*auxV */
354:     dvd_aux_matmult(data,auxV,wx);

356:     /* aux <- K^{1/2} * in */
357:     for (i=0;i<n;i++) {
358:       PCApplySymmetricLeft(data->old_pc,wx[i],outx[i]);
359:     }
360:     break;
361:   default:
362:     SETERRQ(PETSC_COMM_SELF,1,"Unsupported KSP side");
363:   }
364:   /* out <- out - v*(u'*out) */
365:   dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
366:   SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
367:   return(0);
368: }

372: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out)
373: {
375:   dvdImprovex_jd *data;
376:   PetscInt       n,i;
377:   const Vec      *inx, *outx;
378:   Mat            A;

381:   PCGetOperators(pc,&A,NULL);
382:   MatShellGetContext(A,(void**)&data);
383:   VecCompGetSubVecs(in,NULL,&inx);
384:   VecCompGetSubVecs(out,NULL,&outx);
385:   n = data->r_e - data->r_s;
386:   /* out <- K * in */
387:   for (i=0;i<n;i++) {
388:     data->d->improvex_precond(data->d,data->r_s+i,inx[i],outx[i]);
389:   }
390:   /* out <- out - v*(u'*out) */
391:   dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
392:   return(0);
393: }

397: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out)
398: {
400:   dvdImprovex_jd *data;
401:   PetscInt       n,i;
402:   const Vec      *inx, *outx;
403:   Vec            *auxV;
404:   Mat            A;

407:   PCGetOperators(pc,&A,NULL);
408:   MatShellGetContext(A,(void**)&data);
409:   VecCompGetSubVecs(in,NULL,&inx);
410:   VecCompGetSubVecs(out,NULL,&outx);
411:   n = data->r_e - data->r_s;
412:   SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
413:   /* auxV <- in */
414:   for (i=0;i<n;i++) {
415:     VecCopy(inx[i],auxV[i]);
416:   }
417:   /* auxV <- auxV - u*(v'*auxV) */
418:   dvd_improvex_applytrans_proj(data->d,auxV,n);
419:   /* out <- K' * aux */
420:   for (i=0;i<n;i++) {
421:     PCApplyTranspose(data->old_pc,auxV[i],outx[i]);
422:   }
423:   SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
424:   return(0);
425: }

429: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out)
430: {
432:   dvdImprovex_jd *data;
433:   PetscInt       n;
434:   const Vec      *inx, *outx;
435:   PCSide         side;

438:   MatShellGetContext(A,(void**)&data);
439:   VecCompGetSubVecs(in,NULL,&inx);
440:   VecCompGetSubVecs(out,NULL,&outx);
441:   n = data->r_e - data->r_s;
442:   /* out <- theta[1]A*in - theta[0]*B*in */
443:   dvd_aux_matmult(data,inx,outx);
444:   KSPGetPCSide(data->ksp,&side);
445:   if (side == PC_RIGHT) {
446:     /* out <- out - v*(u'*out) */
447:     dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
448:   }
449:   return(0);
450: }

454: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out)
455: {
457:   dvdImprovex_jd *data;
458:   PetscInt       n,i;
459:   const Vec      *inx,*outx,*r;
460:   Vec            *auxV;
461:   PCSide         side;

464:   MatShellGetContext(A,(void**)&data);
465:   VecCompGetSubVecs(in,NULL,&inx);
466:   VecCompGetSubVecs(out,NULL,&outx);
467:   n = data->r_e - data->r_s;
468:   KSPGetPCSide(data->ksp,&side);
469:   if (side == PC_RIGHT) {
470:     /* auxV <- in */
471:     SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
472:     for (i=0;i<n;i++) {
473:       VecCopy(inx[i],auxV[i]);
474:     }
475:     /* auxV <- auxV - v*(u'*auxV) */
476:     dvd_improvex_applytrans_proj(data->d,auxV,n);
477:     r = auxV;
478:   } else r = inx;
479:   /* out <- theta[1]A*r - theta[0]*B*r */
480:   dvd_aux_matmulttrans(data,r,outx);
481:   if (side == PC_RIGHT) {
482:     SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
483:   }
484:   return(0);
485: }

489: static PetscErrorCode MatCreateVecs_dvd_jd(Mat A,Vec *right,Vec *left)
490: {
492:   Vec            *r,*l;
493:   dvdImprovex_jd *data;
494:   PetscInt       n,i;

497:   MatShellGetContext(A,(void**)&data);
498:   n = data->ksp_max_size;
499:   if (right) {
500:     PetscMalloc1(n,&r);
501:   }
502:   if (left) {
503:     PetscMalloc1(n,&l);
504:   }
505:   for (i=0;i<n;i++) {
506:     MatCreateVecs(data->d->A,right?&r[i]:NULL,left?&l[i]:NULL);
507:   }
508:   if (right) {
509:     VecCreateCompWithVecs(r,n,data->friends,right);
510:     for (i=0;i<n;i++) {
511:       VecDestroy(&r[i]);
512:     }
513:   }
514:   if (left) {
515:     VecCreateCompWithVecs(l,n,data->friends,left);
516:     for (i=0;i<n;i++) {
517:       VecDestroy(&l[i]);
518:     }
519:   }

521:   if (right) {
522:     PetscFree(r);
523:   }
524:   if (left) {
525:     PetscFree(l);
526:   }
527:   return(0);
528: }

532: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
533: {
535:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
536:   PetscInt       rA, cA, rlA, clA;
537:   Mat            A;
538:   PetscBool      t;
539:   PC             pc;
540:   Vec            v0[2];

543:   data->size_cX = data->old_size_X = 0;
544:   data->lastTol = data->dynamic?0.5:0.0;

546:   /* Setup the ksp */
547:   if (data->ksp) {
548:     /* Create the reference vector */
549:     BVGetColumn(d->eps->V,0,&v0[0]);
550:     v0[1] = v0[0];
551:     VecCreateCompWithVecs(v0,data->ksp_max_size,NULL,&data->friends);
552:     BVRestoreColumn(d->eps->V,0,&v0[0]);
553:     PetscLogObjectParent((PetscObject)d->eps,(PetscObject)data->friends);

555:     /* Save the current pc and set a PCNONE */
556:     KSPGetPC(data->ksp, &data->old_pc);
557:     PetscObjectTypeCompare((PetscObject)data->old_pc,PCNONE,&t);
558:     data->lastTol = 0.5;
559:     if (t) data->old_pc = 0;
560:     else {
561:       PetscObjectReference((PetscObject)data->old_pc);
562:       PCCreate(PetscObjectComm((PetscObject)d->eps),&pc);
563:       PCSetType(pc,PCSHELL);
564:       PCSetOperators(pc,d->A,d->A);
565:       PCSetReusePreconditioner(pc,PETSC_TRUE);
566:       PCShellSetApply(pc,PCApply_dvd);
567:       PCShellSetApplyBA(pc,PCApplyBA_dvd);
568:       PCShellSetApplyTranspose(pc,PCApplyTranspose_dvd);
569:       KSPSetPC(data->ksp,pc);
570:       PCDestroy(&pc);
571:     }

573:     /* Create the (I-v*u')*K*(A-s*B) matrix */
574:     MatGetSize(d->A,&rA,&cA);
575:     MatGetLocalSize(d->A,&rlA,&clA);
576:     MatCreateShell(PetscObjectComm((PetscObject)d->A),rlA*data->ksp_max_size,clA*data->ksp_max_size,rA*data->ksp_max_size,cA*data->ksp_max_size,data,&A);
577:     MatShellSetOperation(A,MATOP_MULT,(void(*)(void))MatMult_dvd_jd);
578:     MatShellSetOperation(A,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_dvd_jd);
579:     MatShellSetOperation(A,MATOP_GET_VECS,(void(*)(void))MatCreateVecs_dvd_jd);

581:     /* Try to avoid KSPReset */
582:     KSPGetOperatorsSet(data->ksp,&t,NULL);
583:     if (t) {
584:       Mat      M;
585:       PetscInt rM;
586:       KSPGetOperators(data->ksp,&M,NULL);
587:       MatGetSize(M,&rM,NULL);
588:       if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
589:     }
590:     KSPSetOperators(data->ksp,A,A);
591:     KSPSetReusePreconditioner(data->ksp,PETSC_TRUE);
592:     KSPSetUp(data->ksp);
593:     MatDestroy(&A);
594:   } else {
595:     data->old_pc = 0;
596:     data->friends = NULL;
597:   }
598:   BVSetActiveColumns(data->KZ,0,0);
599:   BVSetActiveColumns(data->U,0,0);
600:   return(0);
601: }

605: /*
606:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
607:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
608:   where
609:   pX,pY, the right and left eigenvectors of the projected system
610:   ld, the leading dimension of pX and pY
611: */
612: static PetscErrorCode dvd_improvex_jd_proj_cuv(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
613: {
614: #if defined(PETSC_MISSING_LAPACK_GETRF)
616:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable");
617: #else
619:   PetscInt       n=i_e-i_s,size_KZ,V_new,rm,i,lv,kv,lKZ,kKZ;
620:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
621:   Vec            u[2],v[2];
622:   PetscBLASInt   s,ldXKZ,info;

625:   /* Check consistency */
626:   BVGetActiveColumns(d->eps->V,&lv,&kv);
627:   V_new = lv - data->size_cX;
628:   if (V_new > data->old_size_X) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
629:   data->old_size_X = n;
630:   data->size_cX = lv;

632:   /* KZ <- KZ(rm:rm+max_cX-1) */
633:   BVGetActiveColumns(data->KZ,&lKZ,&kKZ);
634:   rm = PetscMax(V_new+lKZ-d->max_cX_in_impr,0);
635:   if (rm > 0) {
636:     for (i=0;i<lKZ;i++) {
637:       BVCopyColumn(data->KZ,i+rm,i);
638:       BVCopyColumn(data->U,i+rm,i);
639:     }
640:   }

642:   /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
643:   if (rm > 0) {
644:     for (i=0; i<lKZ; i++) {
645:       PetscMemcpy(&data->XKZ[i*data->ldXKZ+i],&data->XKZ[(i+rm)*data->ldXKZ+i+rm],sizeof(PetscScalar)*lKZ);
646:     }
647:   }
648:   lKZ = PetscMin(d->max_cX_in_impr,lKZ+V_new);
649:   BVSetActiveColumns(data->KZ,lKZ,lKZ+n);
650:   BVSetActiveColumns(data->U,lKZ,lKZ+n);

652:   /* Compute X, KZ and KR */
653:   BVGetColumn(data->U,lKZ,u);
654:   if (n>1) { BVGetColumn(data->U,lKZ+1,&u[1]); }
655:   BVGetColumn(data->KZ,lKZ,v);
656:   if (n>1) { BVGetColumn(data->KZ,lKZ+1,&v[1]); }
657:   d->improvex_jd_proj_uv(d,i_s,i_e,u,v,kr,theta,thetai,pX,pY,ld);
658:   BVRestoreColumn(data->U,lKZ,u);
659:   if (n>1) { BVRestoreColumn(data->U,lKZ+1,&u[1]); }
660:   BVRestoreColumn(data->KZ,lKZ,v);
661:   if (n>1) { BVRestoreColumn(data->KZ,lKZ+1,&v[1]); }

663:   /* XKZ <- U'*KZ */
664:   BVMultS(data->KZ,data->U,data->XKZ,data->ldXKZ);

666:   /* iXKZ <- inv(XKZ) */
667:   size_KZ = lKZ+n;
668:   PetscBLASIntCast(lKZ+n,&s);
669:   data->ldiXKZ = data->size_iXKZ = size_KZ;
670:   for (i=0;i<size_KZ;i++) {
671:     PetscMemcpy(&data->iXKZ[data->ldiXKZ*i],&data->XKZ[data->ldXKZ*i],sizeof(PetscScalar)*size_KZ);
672:   }
673:   PetscBLASIntCast(data->ldiXKZ,&ldXKZ);
674:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
675:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&s,&s,data->iXKZ,&ldXKZ,data->iXKZPivots,&info));
676:   PetscFPTrapPop();
677:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack XGETRF %d",info);
678:   return(0);
679: #endif
680: }

684: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
685: {
686:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
688:   PetscInt       i,j,n,maxits,maxits0,lits,s,ld,k,max_size_D,lV,kV;
689:   PetscScalar    *pX,*pY;
690:   PetscReal      tol,tol0;
691:   Vec            *kr,kr_comp,D_comp,D[2],kr0[2];
692:   PetscBool      odd_situation = PETSC_FALSE;

695:   BVGetActiveColumns(d->eps->V,&lV,&kV);
696:   max_size_D = d->eps->ncv-kV;
697:   /* Quick exit */
698:   if ((max_size_D == 0) || r_e-r_s <= 0) {
699:    *size_D = 0;
700:     return(0);
701:   }

703:   n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
704:   if (n == 0) SETERRQ(PETSC_COMM_SELF,1,"n == 0");
705:   if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1,"size_X < r_e-r_s");

707:   DSGetLeadingDimension(d->eps->ds,&ld);

709:   /* Restart lastTol if a new pair converged */
710:   if (data->dynamic && data->size_cX < lV)
711:     data->lastTol = 0.5;

713:   for (i=0,s=0;i<n;i+=s) {
714:     /* If the selected eigenvalue is complex, but the arithmetic is real... */
715: #if !defined(PETSC_USE_COMPLEX)
716:     if (d->eigi[i] != 0.0) {
717:       if (i+2 <= max_size_D) s=2;
718:       else break;
719:     } else
720: #endif
721:       s=1;

723:     data->r_s = r_s+i;
724:     data->r_e = r_s+i+s;
725:     SlepcVecPoolGetVecs(d->auxV,s,&kr);

727:     /* Compute theta, maximum iterations and tolerance */
728:     maxits = 0;
729:     tol = 1;
730:     for (j=0;j<s;j++) {
731:       d->improvex_jd_lit(d,r_s+i+j,&data->theta[2*j],&data->thetai[j],&maxits0,&tol0);
732:       maxits += maxits0;
733:       tol *= tol0;
734:     }
735:     maxits/= s;
736:     tol = data->dynamic?data->lastTol:PetscExpReal(PetscLogReal(tol)/s);

738:     /* Compute u, v and kr */
739:     k = r_s+i;
740:     DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
741:     DSNormalize(d->eps->ds,DS_MAT_X,r_s+i);
742:     k = r_s+i;
743:     DSVectors(d->eps->ds,DS_MAT_Y,&k,NULL);
744:     DSNormalize(d->eps->ds,DS_MAT_Y,r_s+i);
745:     DSGetArray(d->eps->ds,DS_MAT_X,&pX);
746:     DSGetArray(d->eps->ds,DS_MAT_Y,&pY);
747:     dvd_improvex_jd_proj_cuv(d,r_s+i,r_s+i+s,kr,data->theta,data->thetai,pX,pY,ld);
748:     DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
749:     DSRestoreArray(d->eps->ds,DS_MAT_Y,&pY);

751:     /* Check if the first eigenpairs are converged */
752:     if (i == 0) {
753:       d->preTestConv(d,0,s,s,&d->npreconv);
754:       if (d->npreconv > 0) break;
755:     }

757:     /* Test the odd situation of solving Ax=b with A=I */
758: #if !defined(PETSC_USE_COMPLEX)
759:     odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && data->thetai[0] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
760: #else
761:     odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && d->B == NULL)? PETSC_TRUE: PETSC_FALSE;
762: #endif
763:     /* If JD */
764:     if (data->ksp && !odd_situation) {
765:       /* kr <- -kr */
766:       for (j=0;j<s;j++) {
767:         VecScale(kr[j],-1.0);
768:       }

770:       /* Compose kr and D */
771:       kr0[0] = kr[0];
772:       kr0[1] = (s==2 ? kr[1] : NULL);
773:       VecCreateCompWithVecs(kr0,data->ksp_max_size,data->friends,&kr_comp);
774:       BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
775:       if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
776:       else D[1] = NULL;
777:       VecCreateCompWithVecs(D,data->ksp_max_size,data->friends,&D_comp);
778:       VecCompSetSubVecs(data->friends,s,NULL);

780:       /* Solve the correction equation */
781:       KSPSetTolerances(data->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,maxits);
782:       KSPSolve(data->ksp,kr_comp,D_comp);
783:       KSPGetIterationNumber(data->ksp,&lits);

785:       /* Destroy the composed ks and D */
786:       VecDestroy(&kr_comp);
787:       VecDestroy(&D_comp);
788:       BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
789:       if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }

791:     /* If GD */
792:     } else {
793:       BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
794:       if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
795:       for (j=0;j<s;j++) {
796:         d->improvex_precond(d,r_s+i+j,kr[j],D[j]);
797:       }
798:       dvd_improvex_apply_proj(d,D,s);
799:       BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
800:       if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
801:     }
802:     /* Prevent that short vectors are discarded in the orthogonalization */
803:     if (i == 0 && d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
804:       for (j=0;j<s;j++) {
805:         BVScaleColumn(d->eps->V,kV+r_s+i+j,1.0/d->eps->errest[d->nconv+r_s]);
806:       }
807:     }
808:     SlepcVecPoolRestoreVecs(d->auxV,s,&kr);
809:   }
810:   *size_D = i;
811:   if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);
812:   return(0);
813: }

817: PetscErrorCode dvd_improvex_jd(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs,PetscInt cX_impr,PetscBool dynamic)
818: {
820:   dvdImprovex_jd *data;
821:   PetscBool      useGD;
822:   PC             pc;
823:   PetscInt       size_P;

826:   /* Setting configuration constrains */
827:   PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&useGD);

829:   /* If the arithmetic is real and the problem is not Hermitian, then
830:      the block size is incremented in one */
831: #if !defined(PETSC_USE_COMPLEX)
832:   if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) {
833:     max_bs++;
834:     b->max_size_P = PetscMax(b->max_size_P,2);
835:   } else
836: #endif
837:   {
838:     b->max_size_P = PetscMax(b->max_size_P,1);
839:   }
840:   b->max_size_X = PetscMax(b->max_size_X,max_bs);
841:   size_P = b->max_size_P+cX_impr;

843:   /* Setup the preconditioner */
844:   if (ksp) {
845:     KSPGetPC(ksp,&pc);
846:     dvd_static_precond_PC(d,b,pc);
847:   } else {
848:     dvd_static_precond_PC(d,b,0);
849:   }

851:   /* Setup the step */
852:   if (b->state >= DVD_STATE_CONF) {
853:     PetscNewLog(d->eps,&data);
854:     data->dynamic = dynamic;
855:     d->max_cX_in_impr = cX_impr;
856:     PetscMalloc1(size_P*size_P,&data->XKZ);
857:     PetscMalloc1(size_P*size_P,&data->iXKZ);
858:     PetscMalloc1(size_P,&data->iXKZPivots);
859:     data->ldXKZ = size_P;
860:     data->size_X = b->max_size_X;
861:     d->improveX_data = data;
862:     data->ksp = useGD? NULL: ksp;
863:     data->d = d;
864:     d->improveX = dvd_improvex_jd_gen;
865: #if !defined(PETSC_USE_COMPLEX)
866:     if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) data->ksp_max_size = 2;
867:     else
868: #endif
869:       data->ksp_max_size = 1;
870:     /* Create various vector basis */
871:     BVDuplicateResize(d->eps->V,size_P,&data->KZ);
872:     BVSetMatrix(data->KZ,NULL,PETSC_FALSE);
873:     BVDuplicate(data->KZ,&data->U);

875:     EPSDavidsonFLAdd(&d->startList,dvd_improvex_jd_start);
876:     EPSDavidsonFLAdd(&d->endList,dvd_improvex_jd_end);
877:     EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_jd_d);
878:   }
879:   return(0);
880: }

882: #if !defined(PETSC_USE_COMPLEX)
885: PETSC_STATIC_INLINE PetscErrorCode dvd_complex_rayleigh_quotient(Vec ur,Vec ui,Vec Axr,Vec Axi,Vec Bxr,Vec Bxi,PetscScalar *eigr,PetscScalar *eigi)
886: { 
888:   PetscScalar    rAr,iAr,rAi,iAi,rBr,iBr,rBi,iBi,b0,b2,b4,b6,b7;

891:   /* eigr = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k
892:      eigi = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k
893:      k    =  (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr)    */
894:   VecDotBegin(Axr,ur,&rAr); /* r*A*r */
895:   VecDotBegin(Axr,ui,&iAr); /* i*A*r */
896:   VecDotBegin(Axi,ur,&rAi); /* r*A*i */
897:   VecDotBegin(Axi,ui,&iAi); /* i*A*i */
898:   VecDotBegin(Bxr,ur,&rBr); /* r*B*r */
899:   VecDotBegin(Bxr,ui,&iBr); /* i*B*r */
900:   VecDotBegin(Bxi,ur,&rBi); /* r*B*i */
901:   VecDotBegin(Bxi,ui,&iBi); /* i*B*i */
902:   VecDotEnd(Axr,ur,&rAr); /* r*A*r */
903:   VecDotEnd(Axr,ui,&iAr); /* i*A*r */
904:   VecDotEnd(Axi,ur,&rAi); /* r*A*i */
905:   VecDotEnd(Axi,ui,&iAi); /* i*A*i */
906:   VecDotEnd(Bxr,ur,&rBr); /* r*B*r */
907:   VecDotEnd(Bxr,ui,&iBr); /* i*B*r */
908:   VecDotEnd(Bxi,ur,&rBi); /* r*B*i */
909:   VecDotEnd(Bxi,ui,&iBi); /* i*B*i */
910:   b0 = rAr+iAi; /* rAr+iAi */
911:   b2 = rAi-iAr; /* rAi-iAr */
912:   b4 = rBr+iBi; /* rBr+iBi */
913:   b6 = rBi-iBr; /* rBi-iBr */
914:   b7 = b4*b4 + b6*b6; /* k */
915:   *eigr = (b0*b4 + b2*b6) / b7; /* eig_r */
916:   *eigi = (b2*b4 - b0*b6) / b7; /* eig_i */
917:   return(0);
918: }
919: #endif

923: PETSC_STATIC_INLINE PetscErrorCode dvd_compute_n_rr(PetscInt i_s,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,Vec *u,Vec *Ax,Vec *Bx)
924: {
926:   PetscInt       i;
927:   PetscScalar    b0,b1;

930:   for (i=0; i<n; i++) {
931: #if !defined(PETSC_USE_COMPLEX)
932:     if (eigi[i_s+i] != 0.0) {
933:       PetscScalar eigr0=0.0,eigi0=0.0;
934:       dvd_complex_rayleigh_quotient(u[i],u[i+1],Ax[i],Ax[i+1],Bx[i],Bx[i+1],&eigr0,&eigi0);
935:       if (PetscAbsScalar(eigr[i_s+i]-eigr0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10 || PetscAbsScalar(eigi[i_s+i]-eigi0)/PetscAbsScalar(eigi[i_s+i]) > 1e-10) {
936:         PetscInfo4(u[0],"The eigenvalue %g%+gi is far from its Rayleigh quotient value %g%+gi\n",(double)eigr[i_s+i],(double)eigi[i_s+i],(double)eigr0,(double)eigi0);
937:       }
938:       i++;
939:     } else
940: #endif
941:     {
942:       VecDotBegin(Ax[i],u[i],&b0);
943:       VecDotBegin(Bx[i],u[i],&b1);
944:       VecDotEnd(Ax[i],u[i],&b0);
945:       VecDotEnd(Bx[i],u[i],&b1);
946:       b0 = b0/b1;
947:       if (PetscAbsScalar(eigr[i_s+i]-b0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10) {
948:         PetscInfo4(u[0],"The eigenvalue %g+%g is far from its Rayleigh quotient value %g+%g\n",(double)PetscRealPart(eigr[i_s+i]),(double)PetscImaginaryPart(eigr[i_s+i]),(double)PetscRealPart(b0),(double)PetscImaginaryPart(b0));
949:       }
950:     }
951:   }
952:   return(0);
953: }

957: /*
958:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
959:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
960:   where
961:   pX,pY, the right and left eigenvectors of the projected system
962:   ld, the leading dimension of pX and pY
963: */
964: static PetscErrorCode dvd_improvex_jd_proj_uv_KZX(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,Vec *v,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
965: {
967:   PetscInt       n = i_e-i_s,i;
968:   PetscScalar    *b;
969:   Vec            *Ax,*Bx,*r;
970:   Mat            M;
971:   BV             X;

974:   BVDuplicateResize(d->eps->V,4,&X);
975:   MatCreateSeqDense(PETSC_COMM_SELF,4,4,NULL,&M);
976:   /* u <- X(i) */
977:   dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);

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

981:   /* Bx <- B*X(i) */
982:   Bx = kr;
983:   if (d->BX) {
984:     for (i=i_s; i<i_e; ++i) {
985:       BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
986:     }
987:   } else {
988:     for (i=0;i<n;i++) {
989:       if (d->B) {
990:         MatMult(d->B, u[i], Bx[i]);
991:       } else {
992:         VecCopy(u[i], Bx[i]);
993:       }
994:     }
995:   }

997:   /* Ax <- A*X(i) */
998:   SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
999:   Ax = r;
1000:   for (i=i_s; i<i_e; ++i) {
1001:     BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
1002:   }

1004:   /* v <- Y(i) */
1005:   for (i=i_s; i<i_e; ++i) {
1006:     BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,v[i-i_s],&pY[ld*i]);
1007:   }

1009:   /* Recompute the eigenvalue */
1010:   dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,v,Ax,Bx);

1012:   for (i=0;i<n;i++) {
1013: #if !defined(PETSC_USE_COMPLEX)
1014:     if (d->eigi[i_s+i] != 0.0) {
1015:       /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i'    0            1        0
1016:                                        0         theta_2i'     0        1
1017:                                      theta_2i+1 -thetai_i   -eigr_i -eigi_i
1018:                                      thetai_i    theta_2i+1  eigi_i -eigr_i ] */
1019:       MatDenseGetArray(M,&b);
1020:       b[0] = b[5] = PetscConj(theta[2*i]);
1021:       b[2] = b[7] = -theta[2*i+1];
1022:       b[6] = -(b[3] = thetai[i]);
1023:       b[1] = b[4] = 0.0;
1024:       b[8] = b[13] = 1.0/d->nX[i_s+i];
1025:       b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
1026:       b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
1027:       b[9] = b[12] = 0.0;
1028:       MatDenseRestoreArray(M,&b);
1029:       BVInsertVec(X,0,Ax[i]);
1030:       BVInsertVec(X,1,Ax[i+1]);
1031:       BVInsertVec(X,2,Bx[i]);
1032:       BVInsertVec(X,3,Bx[i+1]);
1033:       BVSetActiveColumns(X,0,4);
1034:       BVMultInPlace(X,M,0,4);
1035:       BVCopyVec(X,0,Ax[i]);
1036:       BVCopyVec(X,1,Ax[i+1]);
1037:       BVCopyVec(X,2,Bx[i]);
1038:       BVCopyVec(X,3,Bx[i+1]);
1039:       i++;
1040:     } else
1041: #endif
1042:     {
1043:       /* [Ax_i Bx_i]*= [ theta_2i'    1/nX_i
1044:                         theta_2i+1  -eig_i/nX_i ] */
1045:       MatDenseGetArray(M,&b);
1046:       b[0] = PetscConj(theta[i*2]);
1047:       b[1] = theta[i*2+1];
1048:       b[4] = 1.0/d->nX[i_s+i];
1049:       b[5] = -d->eigr[i_s+i]/d->nX[i_s+i];
1050:       MatDenseRestoreArray(M,&b);
1051:       BVInsertVec(X,0,Ax[i]);
1052:       BVInsertVec(X,1,Bx[i]);
1053:       BVSetActiveColumns(X,0,2);
1054:       BVMultInPlace(X,M,0,2);
1055:       BVCopyVec(X,0,Ax[i]);
1056:       BVCopyVec(X,1,Bx[i]);
1057:     }
1058:   }
1059:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

1061:   /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
1062:   for (i=0;i<n;i++) {
1063:     d->improvex_precond(d,i_s+i,r[i],v[i]);
1064:   }
1065:   SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);

1067:   /* kr <- P*(Ax - eig_i*Bx) */
1068:   d->calcpairs_proj_res(d,i_s,i_e,kr);
1069:   BVDestroy(&X);
1070:   MatDestroy(&M);
1071:   return(0);
1072: }

1076: /*
1077:   Compute: u <- K^{-1}*X, v <- X,
1078:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
1079:   where
1080:   pX,pY, the right and left eigenvectors of the projected system
1081:   ld, the leading dimension of pX and pY
1082: */
1083: static PetscErrorCode dvd_improvex_jd_proj_uv_KXX(dvdDashboard *d,PetscInt i_s,PetscInt i_e,Vec *u,Vec *v,Vec *kr,PetscScalar *theta,PetscScalar *thetai,PetscScalar *pX,PetscScalar *pY,PetscInt ld)
1084: {
1086:   PetscInt       n = i_e - i_s,i;
1087:   PetscScalar    *b;
1088:   Vec            *Ax,*Bx,*r;
1089:   Mat            M;
1090:   BV             X;

1093:   BVDuplicateResize(d->eps->V,4,&X);
1094:   MatCreateSeqDense(PETSC_COMM_SELF,4,2,NULL,&M);
1095:   /* [v u] <- X(i) Y(i) */
1096:   dvd_improvex_compute_X(d,i_s,i_e,v,pX,ld);
1097:   for (i=i_s; i<i_e; ++i) {
1098:     BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,u[i-i_s],&pY[ld*i]);
1099:   }

1101:   /* Bx <- B*X(i) */
1102:   SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
1103:   Bx = r;
1104:   if (d->BX) {
1105:     for (i=i_s; i<i_e; ++i) {
1106:       BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
1107:     }
1108:   } else {
1109:     if (d->B) {
1110:       for (i=0;i<n;i++) {
1111:         MatMult(d->B,v[i],Bx[i]);
1112:       }
1113:     } else Bx = v;
1114:   }

1116:   /* Ax <- A*X(i) */
1117:   Ax = kr;
1118:   for (i=i_s; i<i_e; ++i) {
1119:     BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
1120:   }

1122:   /* Recompute the eigenvalue */
1123:   dvd_compute_n_rr(i_s,n,d->eigr,d->eigi,u,Ax,Bx);

1125:   for (i=0;i<n;i++) {
1126:     if (d->eigi[i_s+i] == 0.0) {
1127:       /* kr <- Ax -eig*Bx */
1128:       VecAXPBY(kr[i],-d->eigr[i_s+i]/d->nX[i_s+i],1.0/d->nX[i_s+i],Bx[i]);
1129:     } else {
1130:       /* [kr_i kr_i+1 r_i r_i+1]*= [   1        0
1131:                                        0        1
1132:                                     -eigr_i -eigi_i
1133:                                      eigi_i -eigr_i] */
1134:       MatDenseGetArray(M,&b);
1135:       b[0] = b[5] = 1.0/d->nX[i_s+i];
1136:       b[2] = b[7] = -d->eigr[i_s+i]/d->nX[i_s+i];
1137:       b[6] = -(b[3] = d->eigi[i_s+i]/d->nX[i_s+i]);
1138:       b[1] = b[4] = 0.0;
1139:       MatDenseRestoreArray(M,&b);
1140:       BVInsertVec(X,0,kr[i]);
1141:       BVInsertVec(X,1,kr[i+1]);
1142:       BVInsertVec(X,2,r[i]);
1143:       BVInsertVec(X,3,r[i+1]);
1144:       BVSetActiveColumns(X,0,4);
1145:       BVMultInPlace(X,M,0,2);
1146:       BVCopyVec(X,0,kr[i]);
1147:       BVCopyVec(X,1,kr[i+1]);
1148:       i++;
1149:     }
1150:   }
1151:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

1153:   /* kr <- P*kr */
1154:   d->calcpairs_proj_res(d,i_s,i_e,r);
1155:   SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);

1157:   /* u <- K^{-1} X(i) */
1158:   for (i=0;i<n;i++) {
1159:     d->improvex_precond(d,i_s+i,v[i],u[i]);
1160:   }
1161:   return(0);
1162: }

1166: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol)
1167: {
1168:   dvdImprovex_jd *data = (dvdImprovex_jd*)d->improveX_data;
1169:   PetscReal      a;

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

1174:   if (d->nR[i]/a < data->fix) {
1175:     theta[0] = d->eigr[i];
1176:     theta[1] = 1.0;
1177: #if !defined(PETSC_USE_COMPLEX)
1178:     *thetai = d->eigi[i];
1179: #endif
1180:   } else {
1181:     theta[0] = d->target[0];
1182:     theta[1] = d->target[1];
1183: #if !defined(PETSC_USE_COMPLEX)
1184:     *thetai = 0.0;
1185: #endif
1186: }

1188: #if defined(PETSC_USE_COMPLEX)
1189:   if (thetai) *thetai = 0.0;
1190: #endif
1191:   *maxits = data->maxits;
1192:   *tol = data->tol;
1193:   return(0);
1194: }

1198: PetscErrorCode dvd_improvex_jd_lit_const(dvdDashboard *d,dvdBlackboard *b,PetscInt maxits,PetscReal tol,PetscReal fix)
1199: {
1200:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;

1203:   /* Setup the step */
1204:   if (b->state >= DVD_STATE_CONF) {
1205:     data->maxits = maxits;
1206:     data->tol = tol;
1207:     data->fix = fix;
1208:     d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
1209:   }
1210:   return(0);
1211: }

1215: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d,dvdBlackboard *b,ProjType_t p)
1216: {
1218:   /* Setup the step */
1219:   if (b->state >= DVD_STATE_CONF) {
1220:     switch (p) {
1221:     case DVD_PROJ_KXX:
1222:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KXX;
1223:       break;
1224:     case DVD_PROJ_KZX:
1225:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX;
1226:       break;
1227:     }
1228:   }
1229:   return(0);
1230: }