Actual source code: dvd_utils.c

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

  4:   Some utils

  6:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  7:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  8:    Copyright (c) 2002-2012, 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

 28: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
 29:                                        Vec Px);
 30: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px);
 31: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px);
 32: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d);

 34: typedef struct {
 35:   PC pc;
 36: } dvdPCWrapper;
 37: /*
 38:   Create a static preconditioner from a PC
 39: */
 42: PetscErrorCode dvd_static_precond_PC(dvdDashboard *d, dvdBlackboard *b, PC pc)
 43: {
 44:   PetscErrorCode  ierr;
 45:   dvdPCWrapper    *dvdpc;
 46:   Mat             P;
 47:   PetscBool       t0,t1;


 51:   /* Setup the step */
 52:   if (b->state >= DVD_STATE_CONF) {
 53:     /* If the preconditioner is valid */
 54:     if (pc) {
 55:       PetscMalloc(sizeof(dvdPCWrapper), &dvdpc);
 56:       dvdpc->pc = pc;
 57:       PetscObjectReference((PetscObject)pc);
 58:       d->improvex_precond_data = dvdpc;
 59:       d->improvex_precond = dvd_static_precond_PC_0;

 61:       /* PC saves the matrix associated with the linear system, and it has to
 62:          be initialize to a valid matrix */
 63:       PCGetOperatorsSet(pc,NULL,&t0);
 64:       PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t1);
 65:       if (t0 && !t1) {
 66:         PCGetOperators(pc, NULL, &P, NULL);
 67:         PetscObjectReference((PetscObject)P);
 68:         PCSetOperators(pc, P, P, SAME_PRECONDITIONER);
 69:         MatDestroy(&P);
 70:       } else
 71:         d->improvex_precond = dvd_precond_none;

 73:       DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);

 75:     /* Else, use no preconditioner */
 76:     } else
 77:       d->improvex_precond = dvd_precond_none;
 78:   }

 80:   return(0);
 81: }

 85: PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
 86: {
 87:   PetscErrorCode  ierr;
 88:   dvdPCWrapper    *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;


 92:   /* Free local data */
 93:   if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
 94:   PetscFree(d->improvex_precond_data);

 96:   return(0);
 97: }


102: PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d, PetscInt i, Vec x,
103:                                        Vec Px)
104: {
105:   PetscErrorCode  ierr;
106:   dvdPCWrapper    *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;


110:   PCApply(dvdpc->pc, x, Px);
111: 
112:   return(0);
113: }

115: typedef struct {
116:   Vec diagA, diagB;
117: } dvdJacobiPrecond;

121: /*
122:   Create the Jacobi preconditioner for Generalized Eigenproblems
123: */
124: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d, dvdBlackboard *b)
125: {
126:   PetscErrorCode  ierr;
127:   dvdJacobiPrecond
128:                   *dvdjp;
129:   PetscBool       t;


133:   /* Check if the problem matrices support GetDiagonal */
134:   MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
135:   if (t && d->B) {
136:     MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
137:   }

139:   /* Setting configuration constrains */
140:   b->own_vecs+= t?( (d->B == 0)?1:2 ) : 0;

142:   /* Setup the step */
143:   if (b->state >= DVD_STATE_CONF) {
144:     if (t) {
145:       PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
146:       dvdjp->diagA = *b->free_vecs; b->free_vecs++;
147:       MatGetDiagonal(d->A,dvdjp->diagA);
148:       if (d->B) {
149:         dvdjp->diagB = *b->free_vecs; b->free_vecs++;
150:         MatGetDiagonal(d->B, dvdjp->diagB);
151:       } else
152:         dvdjp->diagB = 0;
153:       d->improvex_precond_data = dvdjp;
154:       d->improvex_precond = dvd_jacobi_precond_0;

156:       DVD_FL_ADD(d->destroyList, dvd_improvex_precond_d);

158:     /* Else, use no preconditioner */
159:     } else
160:       d->improvex_precond = dvd_precond_none;
161:   }

163:   return(0);
164: }

168: PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
169: {
170:   PetscErrorCode  ierr;
171:   dvdJacobiPrecond
172:                   *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;


176:   /* Compute inv(D - eig)*x */
177:   if (dvdjp->diagB == 0) {
178:     /* Px <- diagA - l */
179:     VecCopy(dvdjp->diagA, Px);
180:     VecShift(Px, -d->eigr[i]);
181:   } else {
182:     /* Px <- diagA - l*diagB */
183:     VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
184: 
185:   }

187:   /* Px(i) <- x/Px(i) */
188:   VecPointwiseDivide(Px, x, Px);

190:   return(0);
191: }

195: /*
196:   Create a trivial preconditioner
197: */
198: PetscErrorCode dvd_precond_none(dvdDashboard *d, PetscInt i, Vec x, Vec Px)
199: {
200:   PetscErrorCode  ierr;


204:   VecCopy(x, Px);

206:   return(0);
207: }


210: /*
211:   Use of PETSc profiler functions
212: */

214: /* Define stages */
215: #define DVD_STAGE_INITV 0 
216: #define DVD_STAGE_NEWITER 1 
217: #define DVD_STAGE_CALCPAIRS 2 
218: #define DVD_STAGE_IMPROVEX 3
219: #define DVD_STAGE_UPDATEV 4
220: #define DVD_STAGE_ORTHV 5

222: PetscErrorCode dvd_profiler_d(dvdDashboard *d);

224: typedef struct {
225:   PetscErrorCode (*old_initV)(struct _dvdDashboard*);
226:   PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
227:   PetscErrorCode (*old_improveX)(struct _dvdDashboard*, Vec *D,
228:                                  PetscInt max_size_D, PetscInt r_s,
229:                                  PetscInt r_e, PetscInt *size_D);
230:   PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
231:   PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
232: } DvdProfiler;

234: static PetscLogStage stages[6] = {0,0,0,0,0,0};

236: /*** Other things ****/

240: PetscErrorCode dvd_prof_init()
241: {
242:   PetscErrorCode  ierr;

245:   if (!stages[0]) {
246:     PetscLogStageRegister("Dvd_step_initV", &stages[DVD_STAGE_INITV]);
247: 
248:     PetscLogStageRegister("Dvd_step_calcPairs",&stages[DVD_STAGE_CALCPAIRS]);
249:     PetscLogStageRegister("Dvd_step_improveX",&stages[DVD_STAGE_IMPROVEX]);
250:     PetscLogStageRegister("Dvd_step_updateV",&stages[DVD_STAGE_UPDATEV]);
251:     PetscLogStageRegister("Dvd_step_orthV",&stages[DVD_STAGE_ORTHV]);
252:   }
253:   return(0);
254: }

258: PetscErrorCode dvd_initV_prof(dvdDashboard* d)
259: {
260:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
261:   PetscErrorCode  ierr;


265:   PetscLogStagePush(stages[DVD_STAGE_INITV]);
266:   p->old_initV(d);
267:   PetscLogStagePop();

269:   return(0);
270: }

274: PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d)
275: {
276:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
277:   PetscErrorCode  ierr;


281:   PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
282:   p->old_calcPairs(d);
283:   PetscLogStagePop();

285:   return(0);
286: }

290: PetscErrorCode dvd_improveX_prof(dvdDashboard* d, Vec *D, PetscInt max_size_D,
291:                        PetscInt r_s, PetscInt r_e, PetscInt *size_D)
292: {
293:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
294:   PetscErrorCode  ierr;


298:   PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
299:   p->old_improveX(d, D, max_size_D, r_s, r_e, size_D);
300:   PetscLogStagePop();

302:   return(0);
303: }

307: PetscErrorCode dvd_updateV_prof(dvdDashboard *d)
308: {
309:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
310:   PetscErrorCode  ierr;


314:   PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
315:   p->old_updateV(d);
316:   PetscLogStagePop();

318:   return(0);
319: }

323: PetscErrorCode dvd_orthV_prof(dvdDashboard *d)
324: {
325:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
326:   PetscErrorCode  ierr;


330:   PetscLogStagePush(stages[DVD_STAGE_ORTHV]);
331:   p->old_orthV(d);
332:   PetscLogStagePop();

334:   return(0);
335: }

339: PetscErrorCode dvd_profiler(dvdDashboard *d, dvdBlackboard *b)
340: {
341:   PetscErrorCode  ierr;
342:   DvdProfiler     *p;


346:   /* Setup the step */
347:   if (b->state >= DVD_STATE_CONF) {
348:     PetscFree(d->prof_data);
349:     PetscMalloc(sizeof(DvdProfiler), &p);
350:     d->prof_data = p;
351:     p->old_initV = d->initV; d->initV = dvd_initV_prof;
352:     p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
353:     p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
354:     p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;

356:     DVD_FL_ADD(d->destroyList, dvd_profiler_d);
357:   }

359:   return(0);
360: }

364: PetscErrorCode dvd_profiler_d(dvdDashboard *d)
365: {
366:   PetscErrorCode  ierr;
367:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;


371:   /* Free local data */
372:   PetscFree(p);

374:   return(0);
375: }



379: /*
380:   Configure the harmonics.
381:   switch(mode) {
382:   DVD_HARM_RR:    harmonic RR
383:   DVD_HARM_RRR:   relative harmonic RR
384:   DVD_HARM_REIGS: rightmost eigenvalues
385:   DVD_HARM_LEIGS: largest eigenvalues
386:   }
387:   fixedTarged, if true use the target instead of the best eigenvalue
388:   target, the fixed target to be used
389: */
390: typedef struct {
391:   PetscScalar
392:     Wa, Wb,       /* span{W} = span{Wa*AV - Wb*BV} */
393:     Pa, Pb;       /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
394:   PetscBool
395:     withTarget;
396:   HarmType_t
397:     mode;
398: } dvdHarmonic;

400: PetscErrorCode dvd_harm_start(dvdDashboard *d);
401: PetscErrorCode dvd_harm_end(dvdDashboard *d);
402: PetscErrorCode dvd_harm_d(dvdDashboard *d);
403: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t);
404: PetscErrorCode dvd_harm_updateW(dvdDashboard *d);
405: PetscErrorCode dvd_harm_proj(dvdDashboard *d);
406: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d);
407: PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi);

411: PetscErrorCode dvd_harm_conf(dvdDashboard *d, dvdBlackboard *b,
412:                              HarmType_t mode, PetscBool fixedTarget,
413:                              PetscScalar t)
414: {
415:   PetscErrorCode  ierr;
416:   dvdHarmonic     *dvdh;


420:   /* Set the problem to GNHEP:
421:      d->G maybe is upper triangular due to biorthogonality of V and W */
422:   d->sEP = d->sA = d->sB = 0;

424:   /* Setup the step */
425:   if (b->state >= DVD_STATE_CONF) {
426:     PetscMalloc(sizeof(dvdHarmonic), &dvdh);
427:     dvdh->withTarget = fixedTarget;
428:     dvdh->mode = mode;
429:     if (fixedTarget) dvd_harm_transf(dvdh, t);
430:     d->calcpairs_W_data = dvdh;
431:     d->calcpairs_W = dvd_harm_updateW;
432:     d->calcpairs_proj_trans = dvd_harm_proj;
433:     d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
434:     d->calcpairs_eig_backtrans = dvd_harm_eig_backtrans;

436:     DVD_FL_ADD(d->destroyList, dvd_harm_d);
437:   }

439:   return(0);
440: }


445: PetscErrorCode dvd_harm_d(dvdDashboard *d)
446: {
447:   PetscErrorCode  ierr;

450:   /* Free local data */
451:   PetscFree(d->calcpairs_W_data);
452:   return(0);
453: }


458: PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh, PetscScalar t)
459: {

462:   switch(dvdh->mode) {
463:   case DVD_HARM_RR:    /* harmonic RR */
464:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
465:   case DVD_HARM_RRR:   /* relative harmonic RR */
466:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
467:   case DVD_HARM_REIGS: /* rightmost eigenvalues */
468:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
469:     break;
470:   case DVD_HARM_LEIGS: /* largest eigenvalues */
471:     dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
472:   case DVD_HARM_NONE:
473:   default:
474:     SETERRQ(PETSC_COMM_SELF,1, "Harmonic type not supported");
475:   }

477:   /* Check the transformation does not change the sign of the imaginary part */
478: #if !defined(PETSC_USE_COMPLEX)
479:   if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
480:     dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
481: #endif

483:   return(0);
484: }

488: PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
489: {
490:   dvdHarmonic     *data = (dvdHarmonic*)d->calcpairs_W_data;
491:   PetscErrorCode  ierr;
492:   PetscInt        i;

495:   /* Update the target if it is necessary */
496:   if (!data->withTarget) {
497:     dvd_harm_transf(data,d->eigr[0]);
498:   }

500:   for(i=d->V_new_s; i<d->V_new_e; i++) {
501:     /* W(i) <- Wa*AV(i) - Wb*BV(i) */
502:     VecAXPBYPCZ(d->W[i],data->Wa,-data->Wb,0.0,d->AV[i],(d->BV?d->BV:d->V)[i]);
503:   }
504:   return(0);
505: }

509: PetscErrorCode dvd_harm_proj(dvdDashboard *d)
510: {
511:   dvdHarmonic     *data = (dvdHarmonic*)d->calcpairs_W_data;
512:   PetscInt        i,j;

515:   if (d->sH != d->sG) SETERRQ(PETSC_COMM_SELF,1,"Projected matrices H and G must have the same structure");

517:   /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
518:   if (DVD_ISNOT(d->sH,DVD_MAT_LTRIANG))     /* Upper triangular part */
519:     for(i=d->V_new_s+d->cX_in_H; i<d->V_new_e+d->cX_in_H; i++)
520:       for(j=0; j<=i; j++) {
521:         PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
522:         d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
523:         d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
524:       }
525:   if (DVD_ISNOT(d->sH,DVD_MAT_UTRIANG))     /* Lower triangular part */
526:     for(i=0; i<d->V_new_e+d->cX_in_H; i++)
527:       for(j=PetscMax(d->V_new_s+d->cX_in_H,i+(DVD_ISNOT(d->sH,DVD_MAT_LTRIANG)?1:0));
528:           j<d->V_new_e+d->cX_in_H; j++) {
529:         PetscScalar h = d->H[d->ldH*i+j], g = d->G[d->ldH*i+j];
530:         d->H[d->ldH*i+j] = data->Pa*h - data->Pb*g;
531:         d->G[d->ldH*i+j] = data->Wa*h - data->Wb*g;
532:       }

534:   return(0);
535: }

539: PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data, PetscScalar *ar,
540:                                   PetscScalar *ai)
541: {
542:   PetscScalar xr;
543: #if !defined(PETSC_USE_COMPLEX)
544:   PetscScalar xi, k;
545: #endif

549:   xr = *ar;
550: #if !defined(PETSC_USE_COMPLEX)
552:   xi = *ai;

554:   if (xi != 0.0) {
555:     k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
556:         data->Wa*data->Wa*xi*xi;
557:     *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
558:            data->Wb*data->Wa*(xr*xr + xi*xi))/k;
559:     *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
560:   } else
561: #endif
562:     *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);

564:   return(0);
565: }

569: PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi)
570: {
571:   dvdHarmonic     *data = (dvdHarmonic*)d->calcpairs_W_data;
572:   PetscErrorCode  ierr;

575:   dvd_harm_backtrans(data,&ar,&ai);
576:   *br = ar;
577:   *bi = ai;
578:   return(0);
579: }

583: PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
584: {
585:   dvdHarmonic     *data = (dvdHarmonic*)d->calcpairs_W_data;
586:   PetscInt        i;
587:   PetscErrorCode  ierr;

590:   for(i=0; i<d->size_H; i++) {
591:     dvd_harm_backtrans(data, &d->eigr[i-d->cX_in_H], &d->eigi[i-d->cX_in_H]);
592:   }
593:   return(0);
594: }