Actual source code: dvd_improvex.c

slepc-3.5.2 2014-10-10
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-2014, 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: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w);
 30: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out);
 31: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out);
 32: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out);
 33: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out);
 34: static PetscErrorCode MatGetVecs_dvd_jd(Mat A,Vec *right,Vec *left);
 35: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d);
 36: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d);
 37: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d);
 38: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D);
 39: 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);
 40: 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);
 41: 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);
 42: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol);
 43: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV);
 44: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV);

 46: #define size_Z (64*4)

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

 50: typedef struct {
 51:   PetscInt size_X;
 52:   KSP ksp;                /* correction equation solver */
 53:   Vec friends;            /* reference vector for composite vectors */
 54:   PetscScalar theta[4], thetai[2];
 55:                           /* the shifts used in the correction eq. */
 56:   PetscInt maxits,        /* maximum number of iterations */
 57:     r_s, r_e,             /* the selected eigenpairs to improve */
 58:     ksp_max_size;         /* the ksp maximum subvectors size */
 59:   PetscReal tol,          /* the maximum solution tolerance */
 60:     lastTol,              /* last tol for dynamic stopping criterion */
 61:     fix;                  /* tolerance for using the approx. eigenvalue */
 62:   PetscBool
 63:     dynamic;              /* if the dynamic stopping criterion is applied */
 64:   dvdDashboard
 65:     *d;                   /* the currect dvdDashboard reference */
 66:   PC old_pc;              /* old pc in ksp */
 67:   BV KZ;                  /* KZ vecs for the projector KZ*inv(X'*KZ)*X' */
 68:   BV U;                   /* new X vectors */
 69:   PetscScalar
 70:    *XKZ,                  /* X'*KZ */
 71:    *iXKZ;                 /* inverse of XKZ */
 72:   PetscInt
 73:     ldXKZ,                /* leading dimension of XKZ */
 74:     size_iXKZ,            /* size of iXKZ */
 75:     ldiXKZ,               /* leading dimension of iXKZ */
 76:     size_cX,              /* last value of d->size_cX */
 77:     old_size_X;           /* last number of improved vectors */
 78:   PetscBLASInt
 79:     *iXKZPivots;          /* array of pivots */
 80: } dvdImprovex_jd;

 82: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y);
 83: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y);

 87: PetscErrorCode dvd_improvex_jd(dvdDashboard *d,dvdBlackboard *b,KSP ksp,PetscInt max_bs,PetscInt cX_impr,PetscBool dynamic)
 88: {
 89:   PetscErrorCode  ierr;
 90:   dvdImprovex_jd  *data;
 91:   PetscBool       useGD;
 92:   PC              pc;
 93:   PetscInt        size_P;


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

100:   /* If the arithmetic is real and the problem is not Hermitian, then
101:      the block size is incremented in one */
102: #if !defined(PETSC_USE_COMPLEX)
103:   if (!DVD_IS(d->sEP,DVD_EP_HERMITIAN)) {
104:     max_bs++;
105:     b->max_size_P = PetscMax(b->max_size_P,2);
106:   } else
107: #endif
108:   {
109:     b->max_size_P = PetscMax(b->max_size_P,1);
110:   }
111:   b->max_size_X = PetscMax(b->max_size_X,max_bs);
112:   size_P = b->max_size_P+cX_impr;

114:   /* Setup the preconditioner */
115:   if (ksp) {
116:     KSPGetPC(ksp,&pc);
117:     dvd_static_precond_PC(d,b,pc);
118:   } else {
119:     dvd_static_precond_PC(d,b,0);
120:   }

122:   /* Setup the step */
123:   if (b->state >= DVD_STATE_CONF) {
124:     PetscMalloc(sizeof(dvdImprovex_jd),&data);
125:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdImprovex_jd));
126:     data->dynamic = dynamic;
127:     d->max_cX_in_impr = cX_impr;
128:     PetscMalloc1(size_P*size_P,&data->XKZ);
129:     PetscMalloc1(size_P*size_P,&data->iXKZ);
130:     PetscMalloc1(size_P,&data->iXKZPivots);
131:     data->ldXKZ = size_P;
132:     data->size_X = b->max_size_X;
133:     d->improveX_data = data;
134:     data->ksp = useGD?NULL:ksp;
135:     data->d = d;
136:     d->improveX = dvd_improvex_jd_gen;
137:     data->ksp_max_size = max_bs;

139:     /* Create various vector basis */
140:     BVDuplicateResize(d->eps->V,size_P,&data->KZ);
141:     BVSetMatrix(data->KZ,NULL,PETSC_FALSE);
142:     BVDuplicate(data->KZ,&data->U);

144:     EPSDavidsonFLAdd(&d->startList,dvd_improvex_jd_start);
145:     EPSDavidsonFLAdd(&d->endList,dvd_improvex_jd_end);
146:     EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_jd_d);
147:   }
148:   return(0);
149: }

153: static PetscErrorCode dvd_improvex_jd_start(dvdDashboard *d)
154: {
155:   PetscErrorCode  ierr;
156:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
157:   PetscInt        rA, cA, rlA, clA;
158:   Mat             A;
159:   PetscBool       t;
160:   PC              pc;
161:   Vec             v0[2];

164:   data->size_cX = data->old_size_X = 0;
165:   data->lastTol = data->dynamic?0.5:0.0;

167:   /* Setup the ksp */
168:   if (data->ksp) {
169:     /* Create the reference vector */
170:     BVGetColumn(d->eps->V,0,&v0[0]);
171:     v0[1] = v0[0];
172:     VecCreateCompWithVecs(v0,data->ksp_max_size,NULL,&data->friends);
173:     BVRestoreColumn(d->eps->V,0,&v0[0]);
174:     PetscLogObjectParent((PetscObject)d->eps,(PetscObject)data->friends);

176:     /* Save the current pc and set a PCNONE */
177:     KSPGetPC(data->ksp, &data->old_pc);
178:     PetscObjectTypeCompare((PetscObject)data->old_pc,PCNONE,&t);
179:     data->lastTol = 0.5;
180:     if (t) {
181:       data->old_pc = 0;
182:     } else {
183:       PetscObjectReference((PetscObject)data->old_pc);
184:       PCCreate(PetscObjectComm((PetscObject)d->eps),&pc);
185:       PCSetType(pc,PCSHELL);
186:       PCSetOperators(pc,d->A,d->A);
187:       PCSetReusePreconditioner(pc,PETSC_TRUE);
188:       PCShellSetApply(pc,PCApply_dvd);
189:       PCShellSetApplyBA(pc,PCApplyBA_dvd);
190:       PCShellSetApplyTranspose(pc,PCApplyTranspose_dvd);
191:       KSPSetPC(data->ksp,pc);
192:       PCDestroy(&pc);
193:     }

195:     /* Create the (I-v*u')*K*(A-s*B) matrix */
196:     MatGetSize(d->A, &rA, &cA);
197:     MatGetLocalSize(d->A, &rlA, &clA);
198:     MatCreateShell(PetscObjectComm((PetscObject)d->A), rlA*data->ksp_max_size,
199:                           clA*data->ksp_max_size, rA*data->ksp_max_size,
200:                           cA*data->ksp_max_size, data, &A);
201:     MatShellSetOperation(A, MATOP_MULT,(void(*)(void))MatMult_dvd_jd);
202:     MatShellSetOperation(A, MATOP_MULT_TRANSPOSE,(void(*)(void))MatMultTranspose_dvd_jd);
203:     MatShellSetOperation(A, MATOP_GET_VECS,(void(*)(void))MatGetVecs_dvd_jd);

205:     /* Try to avoid KSPReset */
206:     KSPGetOperatorsSet(data->ksp,&t,NULL);
207:     if (t) {
208:       Mat M;
209:       PetscInt rM;
210:       KSPGetOperators(data->ksp,&M,NULL);
211:       MatGetSize(M,&rM,NULL);
212:       if (rM != rA*data->ksp_max_size) { KSPReset(data->ksp); }
213:     }
214:     KSPSetOperators(data->ksp,A,A);
215:     KSPSetReusePreconditioner(data->ksp,PETSC_TRUE);
216:     KSPSetUp(data->ksp);
217:     MatDestroy(&A);
218:   } else {
219:     data->old_pc = 0;
220:     data->friends = NULL;
221:   }
222:   BVSetActiveColumns(data->KZ,0,0);
223:   BVSetActiveColumns(data->U,0,0);
224:   return(0);
225: }

229: static PetscErrorCode dvd_improvex_jd_end(dvdDashboard *d)
230: {
231:   PetscErrorCode  ierr;
232:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;

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

237:   /* Restore the pc of ksp */
238:   if (data->old_pc) {
239:     KSPSetPC(data->ksp, data->old_pc);
240:     PCDestroy(&data->old_pc);
241:   }
242:   return(0);
243: }

247: static PetscErrorCode dvd_improvex_jd_d(dvdDashboard *d)
248: {
249:   PetscErrorCode  ierr;
250:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;

253:   /* Free local data and objects */
254:   PetscFree(data->XKZ);
255:   PetscFree(data->iXKZ);
256:   PetscFree(data->iXKZPivots);
257:   BVDestroy(&data->KZ);
258:   BVDestroy(&data->U);
259:   PetscFree(data);
260:   return(0);
261: }

265: static PetscErrorCode dvd_improvex_jd_gen(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
266: {
267:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
268:   PetscErrorCode  ierr;
269:   PetscInt        i,j,n,maxits,maxits0,lits,s,ld,k,max_size_D,lV,kV;
270:   PetscScalar     *pX,*pY;
271:   PetscReal       tol,tol0;
272:   Vec             *kr,kr_comp,D_comp,D[2],kr0[2];
273:   PetscBool       odd_situation = PETSC_FALSE;

276:   BVGetActiveColumns(d->eps->V,&lV,&kV);
277:   max_size_D = d->eps->ncv-kV;
278:   /* Quick exit */
279:   if ((max_size_D == 0) || r_e-r_s <= 0) {
280:    *size_D = 0;
281:     return(0);
282:   }

284:   n = PetscMin(PetscMin(data->size_X, max_size_D), r_e-r_s);
285:   if (n == 0) SETERRQ(PETSC_COMM_SELF,1, "n == 0");
286:   if (data->size_X < r_e-r_s) SETERRQ(PETSC_COMM_SELF,1, "size_X < r_e-r_s");

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

290:   /* Restart lastTol if a new pair converged */
291:   if (data->dynamic && data->size_cX < lV)
292:     data->lastTol = 0.5;

294:   for (i=0,s=0;i<n;i+=s) {
295:     /* If the selected eigenvalue is complex, but the arithmetic is real... */
296: #if !defined(PETSC_USE_COMPLEX)
297:     if (d->eigi[i] != 0.0) {
298:       if (i+2 <= max_size_D) s=2;
299:       else break;
300:     } else
301: #endif
302:       s=1;

304:     data->r_s = r_s+i; data->r_e = r_s+i+s;
305:     SlepcVecPoolGetVecs(d->auxV,s,&kr);

307:     /* Compute theta, maximum iterations and tolerance */
308:     maxits = 0; tol = 1;
309:     for (j=0;j<s;j++) {
310:       d->improvex_jd_lit(d,r_s+i+j,&data->theta[2*j],&data->thetai[j],&maxits0,&tol0);
311:       maxits+= maxits0; tol*= tol0;
312:     }
313:     maxits/= s; tol = data->dynamic?data->lastTol:PetscExpReal(PetscLogReal(tol)/s);

315:     /* Compute u, v and kr */
316:     k = r_s+i;
317:     DSVectors(d->eps->ds,DS_MAT_X,&k,NULL);
318:     DSNormalize(d->eps->ds,DS_MAT_X,r_s+i);
319:     k = r_s+i;
320:     DSVectors(d->eps->ds,DS_MAT_Y,&k,NULL);
321:     DSNormalize(d->eps->ds,DS_MAT_Y,r_s+i);
322:     DSGetArray(d->eps->ds,DS_MAT_X,&pX);
323:     DSGetArray(d->eps->ds,DS_MAT_Y,&pY);
324:     dvd_improvex_jd_proj_cuv(d,r_s+i,r_s+i+s,kr,data->theta,data->thetai,pX,pY,ld);
325:     DSRestoreArray(d->eps->ds,DS_MAT_X,&pX);
326:     DSRestoreArray(d->eps->ds,DS_MAT_Y,&pY);

328:     /* Check if the first eigenpairs are converged */
329:     if (i == 0) {
330:       d->preTestConv(d,0,s,s,&d->npreconv);
331:       if (d->npreconv > 0) break;
332:     }

334:     /* Test the odd situation of solving Ax=b with A=I */
335: #if !defined(PETSC_USE_COMPLEX)
336:     odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && data->thetai[0] == 0. && d->B == NULL)?PETSC_TRUE:PETSC_FALSE;
337: #else
338:     odd_situation = (data->ksp && data->theta[0] == 1. && data->theta[1] == 0. && d->B == NULL)?PETSC_TRUE:PETSC_FALSE;
339: #endif
340:     /* If JD */
341:     if (data->ksp && !odd_situation) {
342:       /* kr <- -kr */
343:       for (j=0;j<s;j++) {
344:         VecScale(kr[j],-1.0);
345:       }

347:       /* Compose kr and D */
348:       kr0[0] = kr[0];
349:       kr0[1] = (s==2 ? kr[1] : NULL);
350:       VecCreateCompWithVecs(kr0,data->ksp_max_size,data->friends,&kr_comp);
351:       BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
352:       if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
353:       else D[1] = NULL;
354:       VecCreateCompWithVecs(D,data->ksp_max_size,data->friends,&D_comp);
355:       VecCompSetSubVecs(data->friends,s,NULL);

357:       /* Solve the correction equation */
358:       KSPSetTolerances(data->ksp,tol,PETSC_DEFAULT,PETSC_DEFAULT,maxits);
359:       KSPSolve(data->ksp,kr_comp,D_comp);
360:       KSPGetIterationNumber(data->ksp,&lits);

362:       /* Destroy the composed ks and D */
363:       VecDestroy(&kr_comp);
364:       VecDestroy(&D_comp);
365:       BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
366:       if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }

368:     /* If GD */
369:     } else {
370:       BVGetColumn(d->eps->V,kV+r_s+i,&D[0]);
371:       if (s==2) { BVGetColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
372:       for (j=0;j<s;j++) {
373:         d->improvex_precond(d,r_s+i+j,kr[j],D[j]);
374:       }
375:       dvd_improvex_apply_proj(d,D,s);
376:       BVRestoreColumn(d->eps->V,kV+r_s+i,&D[0]);
377:       if (s==2) { BVRestoreColumn(d->eps->V,kV+r_s+i+1,&D[1]); }
378:     }
379:     /* Prevent that short vectors are discarded in the orthogonalization */
380:     if (i == 0 && d->eps->errest[d->nconv+r_s] > PETSC_MACHINE_EPSILON && d->eps->errest[d->nconv+r_s] < PETSC_MAX_REAL) {
381:       for (j=0;j<s;j++) {
382:         BVScaleColumn(d->eps->V,kV+r_s+i+j,1.0/d->eps->errest[d->nconv+r_s]);
383:       }
384:     }
385:     SlepcVecPoolRestoreVecs(d->auxV,s,&kr);
386:   }
387:   *size_D = i;
388:   if (data->dynamic) data->lastTol = PetscMax(data->lastTol/2.0,PETSC_MACHINE_EPSILON*10.0);

390:   return(0);
391: }

393: /* y <- theta[1]A*x - theta[0]*B*x
394:    auxV, two auxiliary vectors */
397: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmult(dvdImprovex_jd *data,const Vec *x,const Vec *y)
398: {
399:   PetscErrorCode  ierr;
400:   PetscInt        n,i;
401:   const Vec       *Bx;
402:   Vec             *auxV;

405:   n = data->r_e - data->r_s;
406:   for (i=0;i<n;i++) {
407:     MatMult(data->d->A,x[i],y[i]);
408:   }

410:   SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
411:   for (i=0;i<n;i++) {
412: #if !defined(PETSC_USE_COMPLEX)
413:     if (data->d->eigi[data->r_s+i] != 0.0) {
414:       if (data->d->B) {
415:         MatMult(data->d->B,x[i],auxV[0]);
416:         MatMult(data->d->B,x[i+1],auxV[1]);
417:         Bx = auxV;
418:       } else Bx = &x[i];

420:       /* y_i   <- [ t_2i+1*A*x_i   - t_2i*Bx_i + ti_i*Bx_i+1;
421:          y_i+1      t_2i+1*A*x_i+1 - ti_i*Bx_i - t_2i*Bx_i+1  ] */
422:       VecAXPBYPCZ(y[i],-data->theta[2*i],data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
423:       VecAXPBYPCZ(y[i+1],-data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
424:       i++;
425:     } else
426: #endif
427:     {
428:       if (data->d->B) {
429:         MatMult(data->d->B,x[i],auxV[0]);
430:         Bx = auxV;
431:       } else Bx = &x[i];
432:       VecAXPBY(y[i],-data->theta[i*2],data->theta[i*2+1],Bx[0]);
433:     }
434:   }
435:   SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
436:   return(0);
437: }

439: /* y <- theta[1]'*A'*x - theta[0]'*B'*x */
442: PETSC_STATIC_INLINE PetscErrorCode dvd_aux_matmulttrans(dvdImprovex_jd *data,const Vec *x,const Vec *y)
443: {
444:   PetscErrorCode  ierr;
445:   PetscInt        n,i;
446:   const Vec       *Bx;
447:   Vec             *auxV;

450:   n = data->r_e - data->r_s;
451:   for (i=0;i<n;i++) {
452:     MatMultTranspose(data->d->A,x[i],y[i]);
453:   }

455:   SlepcVecPoolGetVecs(data->d->auxV,2,&auxV);
456:   for (i=0;i<n;i++) {
457: #if !defined(PETSC_USE_COMPLEX)
458:     if (data->d->eigi[data->r_s+i] != 0.0) {
459:       if (data->d->B) {
460:         MatMultTranspose(data->d->B,x[i],auxV[0]);
461:         MatMultTranspose(data->d->B,x[i+1],auxV[1]);
462:         Bx = auxV;
463:       } else Bx = &x[i];

465:       /* y_i   <- [ t_2i+1*A*x_i   - t_2i*Bx_i - ti_i*Bx_i+1;
466:          y_i+1      t_2i+1*A*x_i+1 + ti_i*Bx_i - t_2i*Bx_i+1  ] */
467:       VecAXPBYPCZ(y[i],-data->theta[2*i],-data->thetai[i],data->theta[2*i+1],Bx[0],Bx[1]);
468:       VecAXPBYPCZ(y[i+1],data->thetai[i],-data->theta[2*i],data->theta[2*i+1],Bx[0],Bx[1]);
469:       i++;
470:     } else
471: #endif
472:     {
473:       if (data->d->B) {
474:         MatMultTranspose(data->d->B,x[i],auxV[0]);
475:         Bx = auxV;
476:       } else Bx = &x[i];
477:       VecAXPBY(y[i],PetscConj(-data->theta[i*2]),PetscConj(data->theta[i*2+1]),Bx[0]);
478:     }
479:   }
480:   SlepcVecPoolRestoreVecs(data->d->auxV,2,&auxV);
481:   return(0);
482: }

486: static PetscErrorCode PCApplyBA_dvd(PC pc,PCSide side,Vec in,Vec out,Vec w)
487: {
488:   PetscErrorCode  ierr;
489:   dvdImprovex_jd  *data;
490:   PetscInt        n,i;
491:   const Vec       *inx,*outx,*wx;
492:   Vec             *auxV;
493:   Mat             A;

496:   PCGetOperators(pc,&A,NULL);
497:   MatShellGetContext(A,(void**)&data);
498:   VecCompGetSubVecs(in,NULL,&inx);
499:   VecCompGetSubVecs(out,NULL,&outx);
500:   VecCompGetSubVecs(w,NULL,&wx);
501:   n = data->r_e - data->r_s;
502:   SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
503:   switch (side) {
504:   case PC_LEFT:
505:     /* aux <- theta[1]A*in - theta[0]*B*in */
506:     dvd_aux_matmult(data,inx,auxV);

508:     /* out <- K * aux */
509:     for (i=0;i<n;i++) {
510:       data->d->improvex_precond(data->d,data->r_s+i,auxV[i],outx[i]);
511:     }
512:     break;
513:   case PC_RIGHT:
514:     /* aux <- K * in */
515:     for (i=0;i<n;i++) {
516:       data->d->improvex_precond(data->d,data->r_s+i,inx[i],auxV[i]);
517:     }

519:     /* out <- theta[1]A*auxV - theta[0]*B*auxV */
520:     dvd_aux_matmult(data,auxV,outx);
521:     break;
522:   case PC_SYMMETRIC:
523:     /* aux <- K^{1/2} * in */
524:     for (i=0;i<n;i++) {
525:       PCApplySymmetricRight(data->old_pc,inx[i],auxV[i]);
526:     }

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

531:     /* aux <- K^{1/2} * in */
532:     for (i=0;i<n;i++) {
533:       PCApplySymmetricLeft(data->old_pc,wx[i],outx[i]);
534:     }
535:     break;
536:   default:
537:     SETERRQ(PETSC_COMM_SELF,1, "Unsupported KSP side");
538:   }
539:   /* out <- out - v*(u'*out) */
540:   dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
541:   SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
542:   return(0);
543: }

547: static PetscErrorCode PCApply_dvd(PC pc,Vec in,Vec out)
548: {
549:   PetscErrorCode  ierr;
550:   dvdImprovex_jd  *data;
551:   PetscInt        n,i;
552:   const Vec       *inx, *outx;
553:   Mat             A;

556:   PCGetOperators(pc,&A,NULL);
557:   MatShellGetContext(A,(void**)&data);
558:   VecCompGetSubVecs(in,NULL,&inx);
559:   VecCompGetSubVecs(out,NULL,&outx);
560:   n = data->r_e - data->r_s;
561:   /* out <- K * in */
562:   for (i=0;i<n;i++) {
563:     data->d->improvex_precond(data->d,data->r_s+i,inx[i],outx[i]);
564:   }
565:   /* out <- out - v*(u'*out) */
566:   dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
567:   return(0);
568: }

572: static PetscErrorCode PCApplyTranspose_dvd(PC pc,Vec in,Vec out)
573: {
574:   PetscErrorCode  ierr;
575:   dvdImprovex_jd  *data;
576:   PetscInt        n,i;
577:   const Vec       *inx, *outx;
578:   Vec             *auxV;
579:   Mat             A;

582:   PCGetOperators(pc,&A,NULL);
583:   MatShellGetContext(A,(void**)&data);
584:   VecCompGetSubVecs(in,NULL,&inx);
585:   VecCompGetSubVecs(out,NULL,&outx);
586:   n = data->r_e - data->r_s;
587:   SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
588:   /* auxV <- in */
589:   for (i=0;i<n;i++) {
590:     VecCopy(inx[i],auxV[i]);
591:   }
592:   /* auxV <- auxV - u*(v'*auxV) */
593:   dvd_improvex_applytrans_proj(data->d,auxV,n);
594:   /* out <- K' * aux */
595:   for (i=0;i<n;i++) {
596:     PCApplyTranspose(data->old_pc,auxV[i],outx[i]);
597:   }
598:   SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
599:   return(0);
600: }

604: static PetscErrorCode MatMult_dvd_jd(Mat A,Vec in,Vec out)
605: {
606:   PetscErrorCode  ierr;
607:   dvdImprovex_jd  *data;
608:   PetscInt        n;
609:   const Vec       *inx, *outx;
610:   PCSide          side;

613:   MatShellGetContext(A,(void**)&data);
614:   VecCompGetSubVecs(in,NULL,&inx);
615:   VecCompGetSubVecs(out,NULL,&outx);
616:   n = data->r_e - data->r_s;
617:   /* out <- theta[1]A*in - theta[0]*B*in */
618:   dvd_aux_matmult(data,inx,outx);
619:   KSPGetPCSide(data->ksp,&side);
620:   if (side == PC_RIGHT) {
621:     /* out <- out - v*(u'*out) */
622:     dvd_improvex_apply_proj(data->d,(Vec*)outx,n);
623:   }
624:   return(0);
625: }

629: static PetscErrorCode MatMultTranspose_dvd_jd(Mat A,Vec in,Vec out)
630: {
631:   PetscErrorCode  ierr;
632:   dvdImprovex_jd  *data;
633:   PetscInt        n,i;
634:   const Vec       *inx,*outx,*r;
635:   Vec             *auxV;
636:   PCSide          side;

639:   MatShellGetContext(A,(void**)&data);
640:   VecCompGetSubVecs(in,NULL,&inx);
641:   VecCompGetSubVecs(out,NULL,&outx);
642:   n = data->r_e - data->r_s;
643:   KSPGetPCSide(data->ksp,&side);
644:   if (side == PC_RIGHT) {
645:     /* auxV <- in */
646:     SlepcVecPoolGetVecs(data->d->auxV,n,&auxV);
647:     for (i=0;i<n;i++) {
648:       VecCopy(inx[i],auxV[i]);
649:     }
650:     /* auxV <- auxV - v*(u'*auxV) */
651:     dvd_improvex_applytrans_proj(data->d,auxV,n);
652:     r = auxV;
653:   } else {
654:     r = inx;
655:   }
656:   /* out <- theta[1]A*r - theta[0]*B*r */
657:   dvd_aux_matmulttrans(data,r,outx);
658:   if (side == PC_RIGHT) {
659:     SlepcVecPoolRestoreVecs(data->d->auxV,n,&auxV);
660:   }
661:   return(0);
662: }

666: static PetscErrorCode MatGetVecs_dvd_jd(Mat A,Vec *right,Vec *left)
667: {
668:   PetscErrorCode  ierr;
669:   Vec             *r, *l;
670:   dvdImprovex_jd  *data;
671:   PetscInt        n, i;

674:   MatShellGetContext(A,(void**)&data);
675:   n = data->ksp_max_size;
676:   if (right) {
677:     PetscMalloc(sizeof(Vec)*n,&r);
678:   }
679:   if (left) {
680:     PetscMalloc(sizeof(Vec)*n,&l);
681:   }
682:   for (i=0; i<n; i++) {
683:     MatGetVecs(data->d->A,right?&r[i]:NULL,left?&l[i]:NULL);
684:   }
685:   if (right) {
686:     VecCreateCompWithVecs(r,n,data->friends,right);
687:     for (i=0; i<n; i++) {
688:       VecDestroy(&r[i]);
689:     }
690:   }
691:   if (left) {
692:     VecCreateCompWithVecs(l,n,data->friends,left);
693:     for (i=0; i<n; i++) {
694:       VecDestroy(&l[i]);
695:     }
696:   }

698:   if (right) {
699:     PetscFree(r);
700:   }
701:   if (left) {
702:     PetscFree(l);
703:   }
704:   return(0);
705: }

709: PetscErrorCode dvd_improvex_jd_proj_uv(dvdDashboard *d,dvdBlackboard *b,ProjType_t p)
710: {
712:   /* Setup the step */
713:   if (b->state >= DVD_STATE_CONF) {
714:     switch (p) {
715:     case DVD_PROJ_KXX:
716:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KXX; break;
717:     case DVD_PROJ_KZX:
718:       d->improvex_jd_proj_uv = dvd_improvex_jd_proj_uv_KZX; break;
719:     }
720:   }
721:   return(0);
722: }

726: /*
727:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
728:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
729:   where
730:   pX,pY, the right and left eigenvectors of the projected system
731:   ld, the leading dimension of pX and pY
732: */
733: 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)
734: {
735: #if defined(PETSC_MISSING_LAPACK_GETRF)
737:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRF - Lapack routine is unavailable");
738: #else
739:   PetscErrorCode    ierr;
740:   PetscInt          n=i_e-i_s,size_KZ,V_new,rm,i,lv,kv,lKZ,kKZ;
741:   dvdImprovex_jd    *data = (dvdImprovex_jd*)d->improveX_data;
742:   Vec               u[2],v[2];
743:   PetscBLASInt      s, ldXKZ, info;

746:   /* Check consistency */
747:   BVGetActiveColumns(d->eps->V,&lv,&kv);
748:   V_new = lv - data->size_cX;
749:   if (V_new > data->old_size_X) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
750:   data->old_size_X = n;
751:   data->size_cX = lv;

753:   /* KZ <- KZ(rm:rm+max_cX-1) */
754:   BVGetActiveColumns(data->KZ,&lKZ,&kKZ);
755:   rm = PetscMax(V_new+lKZ-d->max_cX_in_impr,0);
756:   if (rm > 0) for (i=0; i<lKZ; i++) {
757:       BVCopyColumn(data->KZ,i+rm,i);
758:       BVCopyColumn(data->U,i+rm,i);
759:   }

761:   /* XKZ <- XKZ(rm:rm+max_cX-1,rm:rm+max_cX-1) */
762:   if (rm > 0) for (i=0; i<lKZ; i++) {
763:     PetscMemcpy(&data->XKZ[i*data->ldXKZ+i],&data->XKZ[(i+rm)*data->ldXKZ+i+rm],sizeof(PetscScalar)*lKZ);
764:   }
765:   lKZ = PetscMin(d->max_cX_in_impr,lKZ+V_new);
766:   BVSetActiveColumns(data->KZ,lKZ,lKZ+n);
767:   BVSetActiveColumns(data->U,lKZ,lKZ+n);

769:   /* Compute X, KZ and KR */
770:   BVGetColumn(data->U,lKZ,u);
771:   if (n>1) {BVGetColumn(data->U,lKZ+1,&u[1]);}
772:   BVGetColumn(data->KZ,lKZ,v);
773:   if (n>1) {BVGetColumn(data->KZ,lKZ+1,&v[1]);}
774:   d->improvex_jd_proj_uv(d,i_s,i_e,u,v,kr,theta,thetai,pX,pY,ld);
775:   BVRestoreColumn(data->U,lKZ,u);
776:   if (n>1) {BVRestoreColumn(data->U,lKZ+1,&u[1]);}
777:   BVRestoreColumn(data->KZ,lKZ,v);
778:   if (n>1) {BVRestoreColumn(data->KZ,lKZ+1,&v[1]);}

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

783:   /* iXKZ <- inv(XKZ) */
784:   size_KZ = lKZ+n;
785:   PetscBLASIntCast(lKZ+n,&s);
786:   data->ldiXKZ = data->size_iXKZ = size_KZ;
787:   for (i=0;i<size_KZ;i++) {
788:     PetscMemcpy(&data->iXKZ[data->ldiXKZ*i],&data->XKZ[data->ldXKZ*i],sizeof(PetscScalar)*size_KZ);
789:   }
790:   PetscBLASIntCast(data->ldiXKZ,&ldXKZ);
791:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
792:   PetscStackCallBLAS("LAPACKgetrf",LAPACKgetrf_(&s, &s, data->iXKZ, &ldXKZ, data->iXKZPivots, &info));
793:   PetscFPTrapPop();
794:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack XGETRF %d", info);
795:   return(0);
796: #endif
797: }

801: PETSC_STATIC_INLINE PetscErrorCode dvd_complex_rayleigh_quotient(Vec ur,Vec ui,Vec Axr,Vec Axi,Vec Bxr,Vec Bxi,PetscScalar *eigr,PetscScalar *eigi)
802: { 
803:   PetscErrorCode  ierr;
804:   PetscScalar     rAr,iAr,rAi,iAi,rBr,iBr,rBi,iBi,b0,b2,b4,b6,b7;

807:   /* eigr = [(rAr+iAi)*(rBr+iBi) + (rAi-iAr)*(rBi-iBr)]/k
808:      eigi = [(rAi-iAr)*(rBr+iBi) - (rAr+iAi)*(rBi-iBr)]/k
809:      k    =  (rBr+iBi)*(rBr+iBi) + (rBi-iBr)*(rBi-iBr)    */
810:   VecDotBegin(Axr,ur,&rAr); /* r*A*r */
811:   VecDotBegin(Axr,ui,&iAr); /* i*A*r */
812:   VecDotBegin(Axi,ur,&rAi); /* r*A*i */
813:   VecDotBegin(Axi,ui,&iAi); /* i*A*i */
814:   VecDotBegin(Bxr,ur,&rBr); /* r*B*r */
815:   VecDotBegin(Bxr,ui,&iBr); /* i*B*r */
816:   VecDotBegin(Bxi,ur,&rBi); /* r*B*i */
817:   VecDotBegin(Bxi,ui,&iBi); /* i*B*i */
818:   VecDotEnd(Axr,ur,&rAr); /* r*A*r */
819:   VecDotEnd(Axr,ui,&iAr); /* i*A*r */
820:   VecDotEnd(Axi,ur,&rAi); /* r*A*i */
821:   VecDotEnd(Axi,ui,&iAi); /* i*A*i */
822:   VecDotEnd(Bxr,ur,&rBr); /* r*B*r */
823:   VecDotEnd(Bxr,ui,&iBr); /* i*B*r */
824:   VecDotEnd(Bxi,ur,&rBi); /* r*B*i */
825:   VecDotEnd(Bxi,ui,&iBi); /* i*B*i */
826:   b0 = rAr+iAi; /* rAr+iAi */
827:   b2 = rAi-iAr; /* rAi-iAr */
828:   b4 = rBr+iBi; /* rBr+iBi */
829:   b6 = rBi-iBr; /* rBi-iBr */
830:   b7 = b4*b4 + b6*b6; /* k */
831:   *eigr = (b0*b4 + b2*b6) / b7; /* eig_r */
832:   *eigi = (b2*b4 - b0*b6) / b7; /* eig_i */
833:   return(0);
834: }

838: PETSC_STATIC_INLINE PetscErrorCode dvd_compute_n_rr(PetscInt i_s,PetscInt n,PetscScalar *eigr,PetscScalar *eigi,Vec *u,Vec *Ax,Vec *Bx)
839: {
840:   PetscErrorCode  ierr;
841:   PetscInt        i;
842:   PetscScalar     b0,b1;

845:   for (i=0; i<n; i++) {
846: #if !defined(PETSC_USE_COMPLEX)
847:     if (eigi[i_s+i] != 0.0) {
848:       PetscScalar eigr0=0.0,eigi0=0.0;
849:       dvd_complex_rayleigh_quotient(u[i],u[i+1],Ax[i],Ax[i+1],Bx[i],Bx[i+1],&eigr0,&eigi0);
850:       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) {
851:         PetscInfo4(u[0], "The eigenvalue %g+%g is far from its Rayleigh quotient value %g+%g\n",(double)eigr[i_s+i],(double)eigi[i_s+i],(double)eigr0,(double)eigi0);
852:       }
853:       i++;
854:     } else
855: #endif
856:     {
857:       VecDotBegin(Ax[i],u[i],&b0);
858:       VecDotBegin(Bx[i],u[i],&b1);
859:       VecDotEnd(Ax[i],u[i],&b0);
860:       VecDotEnd(Bx[i],u[i],&b1);
861:       b0 = b0/b1;
862:       if (PetscAbsScalar(eigr[i_s+i]-b0)/PetscAbsScalar(eigr[i_s+i]) > 1e-10) {
863:         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));
864:       }
865:     }
866:   }
867:   return(0);
868: }

872: /*
873:   Compute: u <- X, v <- K*(theta[0]*A+theta[1]*B)*X,
874:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1], Y <- W*pY[i_s..i_e-1]
875:   where
876:   pX,pY, the right and left eigenvectors of the projected system
877:   ld, the leading dimension of pX and pY
878: */
879: 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)
880: {
881:   PetscErrorCode  ierr;
882:   PetscInt        n = i_e-i_s,i;
883:   PetscScalar     *b;
884:   Vec             *Ax,*Bx,*r;
885:   Mat             M;
886:   BV              X;

889:   BVDuplicateResize(d->eps->V,4,&X);
890:   MatCreateSeqDense(PETSC_COMM_SELF,4,4,NULL,&M);
891:   /* u <- X(i) */
892:   dvd_improvex_compute_X(d,i_s,i_e,u,pX,ld);

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

896:   /* Bx <- B*X(i) */
897:   Bx = kr;
898:   if (d->BX) {
899:     for (i=i_s; i<i_e; ++i) {
900:       BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
901:     }
902:   } else {
903:     for (i=0;i<n;i++) {
904:       if (d->B) {
905:         MatMult(d->B, u[i], Bx[i]);
906:       } else {
907:         VecCopy(u[i], Bx[i]);
908:       }
909:     }
910:   }

912:   /* Ax <- A*X(i) */
913:   SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
914:   Ax = r;
915:   for (i=i_s; i<i_e; ++i) {
916:     BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
917:   }

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

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

927:   for (i=0;i<n;i++) {
928: #if !defined(PETSC_USE_COMPLEX)
929:     if (d->eigi[i_s+i] != 0.0) {
930:       /* [r_i r_i+1 kr_i kr_i+1]*= [ theta_2i'    0            1        0
931:                                        0         theta_2i'     0        1
932:                                      theta_2i+1 -thetai_i   -eigr_i -eigi_i
933:                                      thetai_i    theta_2i+1  eigi_i -eigr_i ] */
934:       MatDenseGetArray(M,&b);
935:       b[0] = b[5] = PetscConj(theta[2*i]);
936:       b[2] = b[7] = -theta[2*i+1];
937:       b[6] = -(b[3] = thetai[i]);
938:       b[1] = b[4] = 0.0;
939:       b[8] = b[13] = 1.0/d->nX[i_s+i];
940:       b[10] = b[15] = -d->eigr[i_s+i]/d->nX[i_s+i];
941:       b[14] = -(b[11] = d->eigi[i_s+i]/d->nX[i_s+i]);
942:       b[9] = b[12] = 0.0;
943:       MatDenseRestoreArray(M,&b);
944:       BVInsertVec(X,0,Ax[i]);
945:       BVInsertVec(X,1,Ax[i+1]);
946:       BVInsertVec(X,2,Bx[i]);
947:       BVInsertVec(X,3,Bx[i+1]);
948:       BVSetActiveColumns(X,0,4);
949:       BVMultInPlace(X,M,0,4);
950:       BVCopyVec(X,0,Ax[i]);
951:       BVCopyVec(X,1,Ax[i+1]);
952:       BVCopyVec(X,2,Bx[i]);
953:       BVCopyVec(X,3,Bx[i+1]);
954:       i++;
955:     } else
956: #endif
957:     {
958:       /* [Ax_i Bx_i]*= [ theta_2i'    1/nX_i
959:                         theta_2i+1  -eig_i/nX_i ] */
960:       MatDenseGetArray(M,&b);
961:       b[0] = PetscConj(theta[i*2]);
962:       b[1] = theta[i*2+1];
963:       b[4] = 1.0/d->nX[i_s+i];
964:       b[5] = -d->eigr[i_s+i]/d->nX[i_s+i];
965:       MatDenseRestoreArray(M,&b);
966:       BVInsertVec(X,0,Ax[i]);
967:       BVInsertVec(X,1,Bx[i]);
968:       BVSetActiveColumns(X,0,2);
969:       BVMultInPlace(X,M,0,2);
970:       BVCopyVec(X,0,Ax[i]);
971:       BVCopyVec(X,1,Bx[i]);
972:     }
973:   }
974:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

976:   /* v <- K^{-1} r = K^{-1}(theta_2i'*Ax + theta_2i+1*Bx) */
977:   for (i=0;i<n;i++) {
978:     d->improvex_precond(d,i_s+i,r[i],v[i]);
979:   }
980:   SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);

982:   /* kr <- P*(Ax - eig_i*Bx) */
983:   d->calcpairs_proj_res(d,i_s,i_e,kr);
984:   BVDestroy(&X);
985:   MatDestroy(&M);
986:   return(0);
987: }

991: /*
992:   Compute: u <- K^{-1}*X, v <- X,
993:   kr <- K^{-1}*(A-eig*B)*X, being X <- V*pX[i_s..i_e-1]
994:   where
995:   pX,pY, the right and left eigenvectors of the projected system
996:   ld, the leading dimension of pX and pY
997: */
998: 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)
999: {
1000:   PetscErrorCode  ierr;
1001:   PetscInt        n = i_e - i_s,i;
1002:   PetscScalar     *b;
1003:   Vec             *Ax,*Bx,*r;
1004:   Mat             M;
1005:   BV              X;

1008:   BVDuplicateResize(d->eps->V,4,&X);
1009:   MatCreateSeqDense(PETSC_COMM_SELF,4,2,NULL,&M);
1010:   /* [v u] <- X(i) Y(i) */
1011:   dvd_improvex_compute_X(d,i_s,i_e,v,pX,ld);
1012:   for (i=i_s; i<i_e; ++i) {
1013:     BVMultVec(d->W?d->W:d->eps->V,1.0,0.0,u[i-i_s],&pY[ld*i]);
1014:   }

1016:   /* Bx <- B*X(i) */
1017:   SlepcVecPoolGetVecs(d->auxV,i_e-i_s,&r);
1018:   Bx = r;
1019:   if (d->BX) {
1020:     for (i=i_s; i<i_e; ++i) {
1021:       BVMultVec(d->BX,1.0,0.0,Bx[i-i_s],&pX[ld*i]);
1022:     }
1023:   } else {
1024:     if (d->B) {
1025:       for (i=0;i<n;i++) {
1026:         MatMult(d->B,v[i],Bx[i]);
1027:       }
1028:     } else Bx = v;
1029:   }

1031:   /* Ax <- A*X(i) */
1032:   Ax = kr;
1033:   for (i=i_s; i<i_e; ++i) {
1034:     BVMultVec(d->AX,1.0,0.0,Ax[i-i_s],&pX[ld*i]);
1035:   }

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

1040:   for (i=0;i<n;i++) {
1041:     if (d->eigi[i_s+i] == 0.0) {
1042:       /* kr <- Ax -eig*Bx */
1043:       VecAXPBY(kr[i],-d->eigr[i_s+i]/d->nX[i_s+i],1.0/d->nX[i_s+i],Bx[i]);
1044:     } else {
1045:       /* [kr_i kr_i+1 r_i r_i+1]*= [   1        0
1046:                                        0        1
1047:                                     -eigr_i -eigi_i
1048:                                      eigi_i -eigr_i] */
1049:       MatDenseGetArray(M,&b);
1050:       b[0] = b[5] = 1.0/d->nX[i_s+i];
1051:       b[2] = b[7] = -d->eigr[i_s+i]/d->nX[i_s+i];
1052:       b[6] = -(b[3] = d->eigi[i_s+i]/d->nX[i_s+i]);
1053:       b[1] = b[4] = 0.0;
1054:       MatDenseRestoreArray(M,&b);
1055:       BVInsertVec(X,0,kr[i]);
1056:       BVInsertVec(X,1,kr[i+1]);
1057:       BVInsertVec(X,2,r[i]);
1058:       BVInsertVec(X,3,r[i+1]);
1059:       BVSetActiveColumns(X,0,4);
1060:       BVMultInPlace(X,M,0,2);
1061:       BVCopyVec(X,0,kr[i]);
1062:       BVCopyVec(X,1,kr[i+1]);
1063:       i++;
1064:     }
1065:   }
1066:   for (i=0; i<n; i++) d->nX[i_s+i] = 1.0;

1068:   /* kr <- P*kr */
1069:   d->calcpairs_proj_res(d,i_s,i_e,r);
1070:   SlepcVecPoolRestoreVecs(d->auxV,i_e-i_s,&r);

1072:   /* u <- K^{-1} X(i) */
1073:   for (i=0;i<n;i++) {
1074:     d->improvex_precond(d,i_s+i,v[i],u[i]);
1075:   }
1076:   return(0);
1077: }

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

1086:   /* Setup the step */
1087:   if (b->state >= DVD_STATE_CONF) {
1088:     data->maxits = maxits;
1089:     data->tol = tol;
1090:     data->fix = fix;
1091:     d->improvex_jd_lit = dvd_improvex_jd_lit_const_0;
1092:   }
1093:   return(0);
1094: }

1098: static PetscErrorCode dvd_improvex_jd_lit_const_0(dvdDashboard *d,PetscInt i,PetscScalar* theta,PetscScalar* thetai,PetscInt *maxits,PetscReal *tol)
1099: {
1100:   dvdImprovex_jd  *data = (dvdImprovex_jd*)d->improveX_data;
1101:   PetscReal       a;

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

1106:   if (d->nR[i]/a < data->fix) {
1107:     theta[0] = d->eigr[i];
1108:     theta[1] = 1.0;
1109: #if !defined(PETSC_USE_COMPLEX)
1110:     *thetai = d->eigi[i];
1111: #endif
1112:   } else {
1113:     theta[0] = d->target[0];
1114:     theta[1] = d->target[1];
1115: #if !defined(PETSC_USE_COMPLEX)
1116:     *thetai = 0.0;
1117: #endif
1118: }

1120: #if defined(PETSC_USE_COMPLEX)
1121:   if (thetai) *thetai = 0.0;
1122: #endif
1123:   *maxits = data->maxits;
1124:   *tol = data->tol;
1125:   return(0);
1126: }

1128: /**** Patterns implementation *************************************************/

1130: typedef PetscInt (*funcV0_t)(dvdDashboard*,PetscInt,PetscInt,Vec*);
1131: typedef PetscInt (*funcV1_t)(dvdDashboard*,PetscInt,PetscInt,Vec*,PetscScalar*,Vec);

1135: /* Compute (I - KZ*iXKZ*X')*V where,
1136:    V, the vectors to apply the projector,
1137:    cV, the number of vectors in V,
1138: */
1139: static PetscErrorCode dvd_improvex_apply_proj(dvdDashboard *d,Vec *V,PetscInt cV)
1140: {
1141: #if defined(PETSC_MISSING_LAPACK_GETRS)
1143:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
1144: #else
1145:   PetscErrorCode    ierr;
1146:   dvdImprovex_jd    *data = (dvdImprovex_jd*)d->improveX_data;
1147:   PetscInt          i,ldh,k,l;
1148:   PetscScalar       *h;
1149:   PetscBLASInt      cV_,n,info,ld;
1150: #if defined(PETSC_USE_COMPLEX)
1151:   PetscInt          j;
1152: #endif

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

1157:   /* h <- X'*V */
1158:   PetscMalloc1(data->size_iXKZ*cV,&h);
1159:   ldh = data->size_iXKZ;
1160:   BVGetActiveColumns(data->U,&l,&k);
1161:   if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
1162:   BVSetActiveColumns(data->U,0,k);
1163:   for (i=0; i<cV; i++) {
1164:     BVDotVec(data->U,V[i],&h[ldh*i]);
1165: #if defined(PETSC_USE_COMPLEX)
1166:     for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
1167: #endif
1168:   }
1169:   BVSetActiveColumns(data->U,l,k);

1171:   /* h <- iXKZ\h */
1172:   PetscBLASIntCast(cV,&cV_);
1173:   PetscBLASIntCast(data->size_iXKZ,&n);
1174:   PetscBLASIntCast(data->ldiXKZ,&ld);
1176:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1177:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("N", &n, &cV_, data->iXKZ, &ld, data->iXKZPivots, h, &n, &info));
1178:   PetscFPTrapPop();
1179:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack XGETRS %d", info);

1181:   /* V <- V - KZ*h */
1182:   BVSetActiveColumns(data->KZ,0,k);
1183:   for (i=0; i<cV; i++) {
1184:     BVMultVec(data->KZ,-1.0,1.0,V[i],&h[ldh*i]);
1185:   }
1186:   BVSetActiveColumns(data->KZ,l,k);
1187:   PetscFree(h);
1188:   return(0);
1189: #endif
1190: }

1194: /* Compute (I - X*iXKZ*KZ')*V where,
1195:    V, the vectors to apply the projector,
1196:    cV, the number of vectors in V,
1197: */
1198: static PetscErrorCode dvd_improvex_applytrans_proj(dvdDashboard *d,Vec *V,PetscInt cV)
1199: {
1200: #if defined(PETSC_MISSING_LAPACK_GETRS)
1202:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GETRS - Lapack routines are unavailable");
1203:   return(0);
1204: #else
1205:   PetscErrorCode    ierr;
1206:   dvdImprovex_jd    *data = (dvdImprovex_jd*)d->improveX_data;
1207:   PetscInt          i,ldh,k,l;
1208:   PetscScalar       *h;
1209:   PetscBLASInt      cV_, n, info, ld;
1210: #if defined(PETSC_USE_COMPLEX)
1211:   PetscInt          j;
1212: #endif

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

1217:   /* h <- KZ'*V */
1218:   PetscMalloc1(data->size_iXKZ*cV,&h);
1219:   ldh = data->size_iXKZ;
1220:   BVGetActiveColumns(data->U,&l,&k);
1221:   if (ldh!=k) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
1222:   BVSetActiveColumns(data->KZ,0,k);
1223:   for (i=0; i<cV; i++) {
1224:     BVDotVec(data->KZ,V[i],&h[ldh*i]);
1225: #if defined(PETSC_USE_COMPLEX)
1226:     for (j=0; j<k; j++) h[ldh*i+j] = PetscConj(h[ldh*i+j]);
1227: #endif
1228:   }
1229:   BVSetActiveColumns(data->KZ,l,k);

1231:   /* h <- iXKZ\h */
1232:   PetscBLASIntCast(cV,&cV_);
1233:   PetscBLASIntCast(data->size_iXKZ,&n);
1234:   PetscBLASIntCast(data->ldiXKZ,&ld);
1236:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
1237:   PetscStackCallBLAS("LAPACKgetrs",LAPACKgetrs_("C", &n, &cV_, data->iXKZ, &ld, data->iXKZPivots, h, &n, &info));
1238:   PetscFPTrapPop();
1239:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "Error in Lapack XGETRS %d", info);

1241:   /* V <- V - U*h */
1242:   BVSetActiveColumns(data->U,0,k);
1243:   for (i=0; i<cV; i++) {
1244:     BVMultVec(data->U,-1.0,1.0,V[i],&h[ldh*i]);
1245:   }
1246:   BVSetActiveColumns(data->U,l,k);
1247:   PetscFree(h);
1248:   return(0);
1249: #endif
1250: }