Actual source code: dvd_utils.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:   SLEPc eigensolver: "davidson"

  4:   Some utils

  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

 28: typedef struct {
 29:   PC pc;
 30: } dvdPCWrapper;

 32: /*
 33:   Configure the harmonics.
 34:   switch (mode) {
 35:   DVD_HARM_RR:    harmonic RR
 36:   DVD_HARM_RRR:   relative harmonic RR
 37:   DVD_HARM_REIGS: rightmost eigenvalues
 38:   DVD_HARM_LEIGS: largest eigenvalues
 39:   }
 40:   fixedTarged, if true use the target instead of the best eigenvalue
 41:   target, the fixed target to be used
 42: */
 43: typedef struct {
 44:   PetscScalar
 45:     Wa, Wb,       /* span{W} = span{Wa*AV - Wb*BV} */
 46:     Pa, Pb;       /* H=W'*(Pa*AV - Pb*BV), G=W'*(Wa*AV - Wb*BV) */
 47:   PetscBool
 48:     withTarget;
 49:   HarmType_t
 50:     mode;
 51: } dvdHarmonic;

 53: static PetscErrorCode dvd_static_precond_PC_0(dvdDashboard*,PetscInt,Vec,Vec);
 54: static PetscErrorCode dvd_jacobi_precond_0(dvdDashboard*,PetscInt,Vec,Vec);
 55: static PetscErrorCode dvd_jacobi_precond_d(dvdDashboard*);
 56: static PetscErrorCode dvd_precond_none(dvdDashboard*,PetscInt,Vec,Vec);
 57: static PetscErrorCode dvd_improvex_precond_d(dvdDashboard*);
 58: static PetscErrorCode dvd_initV_prof(dvdDashboard*);
 59: static PetscErrorCode dvd_calcPairs_prof(dvdDashboard*);
 60: static PetscErrorCode dvd_improveX_prof(dvdDashboard*,PetscInt,PetscInt,PetscInt*);
 61: static PetscErrorCode dvd_updateV_prof(dvdDashboard*);
 62: static PetscErrorCode dvd_profiler_d(dvdDashboard*);
 63: static PetscErrorCode dvd_harm_d(dvdDashboard*);
 64: static PetscErrorCode dvd_harm_transf(dvdHarmonic*,PetscScalar);
 65: static PetscErrorCode dvd_harm_updateW(dvdDashboard*);
 66: static PetscErrorCode dvd_harm_proj(dvdDashboard*);
 67: static PetscErrorCode dvd_harm_eigs_trans(dvdDashboard*);
 68: static PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard*,PetscScalar,PetscScalar,PetscScalar*,PetscScalar*);
 69: static PetscErrorCode dvd_harm_backtrans(dvdHarmonic*,PetscScalar*,PetscScalar*);


 72: /*
 73:   Create a static preconditioner from a PC
 74: */
 77: PetscErrorCode dvd_static_precond_PC(dvdDashboard *d,dvdBlackboard *b,PC pc)
 78: {
 80:   dvdPCWrapper   *dvdpc;
 81:   Mat            P;
 82:   PetscBool      t0,t1,t2;

 85:   /* Setup the step */
 86:   if (b->state >= DVD_STATE_CONF) {
 87:     /* If the preconditioner is valid */
 88:     if (pc) {
 89:       PetscMalloc(sizeof(dvdPCWrapper),&dvdpc);
 90:       PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdPCWrapper));
 91:       dvdpc->pc = pc;
 92:       PetscObjectReference((PetscObject)pc);
 93:       d->improvex_precond_data = dvdpc;
 94:       d->improvex_precond = dvd_static_precond_PC_0;

 96:       /* PC saves the matrix associated with the linear system, and it has to
 97:          be initialize to a valid matrix */
 98:       PCGetOperatorsSet(pc,NULL,&t0);
 99:       PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t1);
100:       PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t2);
101:       if (t0 && !t1) {
102:         PCGetOperators(pc,NULL,&P);
103:         PetscObjectReference((PetscObject)P);
104:         PCSetOperators(pc,P,P);
105:         PCSetReusePreconditioner(pc,PETSC_TRUE);
106:         MatDestroy(&P);
107:       } else if (t2) {
108:         PCSetOperators(pc,d->A,d->A);
109:         PCSetReusePreconditioner(pc,PETSC_TRUE);
110:       } else {
111:         d->improvex_precond = dvd_precond_none;
112:       }

114:       EPSDavidsonFLAdd(&d->destroyList,dvd_improvex_precond_d);

116:     /* Else, use no preconditioner */
117:     } else d->improvex_precond = dvd_precond_none;
118:   }
119:   return(0);
120: }

124: static PetscErrorCode dvd_improvex_precond_d(dvdDashboard *d)
125: {
127:   dvdPCWrapper   *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;

130:   /* Free local data */
131:   if (dvdpc->pc) { PCDestroy(&dvdpc->pc); }
132:   PetscFree(d->improvex_precond_data);
133:   return(0);
134: }

138: static PetscErrorCode dvd_static_precond_PC_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
139: {
141:   dvdPCWrapper   *dvdpc = (dvdPCWrapper*)d->improvex_precond_data;

144:   PCApply(dvdpc->pc, x, Px);
145:   return(0);
146: }

148: typedef struct {
149:   Vec diagA, diagB;
150: } dvdJacobiPrecond;

154: /*
155:   Create the Jacobi preconditioner for Generalized Eigenproblems
156: */
157: PetscErrorCode dvd_jacobi_precond(dvdDashboard *d,dvdBlackboard *b)
158: {
159:   PetscErrorCode   ierr;
160:   dvdJacobiPrecond *dvdjp;
161:   PetscBool        t;

164:   /* Check if the problem matrices support GetDiagonal */
165:   MatHasOperation(d->A, MATOP_GET_DIAGONAL, &t);
166:   if (t && d->B) {
167:     MatHasOperation(d->B, MATOP_GET_DIAGONAL, &t);
168:   }

170:   /* Setup the step */
171:   if (b->state >= DVD_STATE_CONF) {
172:     PetscMalloc(sizeof(dvdJacobiPrecond), &dvdjp);
173:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdJacobiPrecond));
174:     if (t) {
175:       MatGetVecs(d->A,&dvdjp->diagA,NULL);
176:       MatGetDiagonal(d->A,dvdjp->diagA);
177:       if (d->B) {
178:         MatGetVecs(d->B,&dvdjp->diagB,NULL);
179:         MatGetDiagonal(d->B,dvdjp->diagB);
180:       } else dvdjp->diagB = 0;
181:       d->improvex_precond_data = dvdjp;
182:       d->improvex_precond = dvd_jacobi_precond_0;

184:       EPSDavidsonFLAdd(&d->destroyList,dvd_jacobi_precond_d);

186:     /* Else, use no preconditioner */
187:     } else {
188:       dvdjp->diagA = dvdjp->diagB = 0;
189:       d->improvex_precond = dvd_precond_none;
190:     }
191:   }
192:   return(0);
193: }

197: static PetscErrorCode dvd_jacobi_precond_0(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
198: {
199:   PetscErrorCode   ierr;
200:   dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;

203:   /* Compute inv(D - eig)*x */
204:   if (dvdjp->diagB == 0) {
205:     /* Px <- diagA - l */
206:     VecCopy(dvdjp->diagA, Px);
207:     VecShift(Px, -d->eigr[i]);
208:   } else {
209:     /* Px <- diagA - l*diagB */
210:     VecWAXPY(Px, -d->eigr[i], dvdjp->diagB, dvdjp->diagA);
211:   }

213:   /* Px(i) <- x/Px(i) */
214:   VecPointwiseDivide(Px, x, Px);
215:   return(0);
216: }

220: static PetscErrorCode dvd_jacobi_precond_d(dvdDashboard *d)
221: {
222:   PetscErrorCode   ierr;
223:   dvdJacobiPrecond *dvdjp = (dvdJacobiPrecond*)d->improvex_precond_data;

226:   if (dvdjp->diagA) {VecDestroy(&dvdjp->diagA);}
227:   if (dvdjp->diagB) {VecDestroy(&dvdjp->diagB);}
228:   PetscFree(d->improvex_precond_data);
229:   return(0);
230: }

234: /*
235:   Create a trivial preconditioner
236: */
237: static PetscErrorCode dvd_precond_none(dvdDashboard *d,PetscInt i,Vec x,Vec Px)
238: {

242:   VecCopy(x, Px);
243:   return(0);
244: }

246: /*
247:   Use of PETSc profiler functions
248: */

250: /* Define stages */
251: #define DVD_STAGE_INITV 0
252: #define DVD_STAGE_NEWITER 1
253: #define DVD_STAGE_CALCPAIRS 2
254: #define DVD_STAGE_IMPROVEX 3
255: #define DVD_STAGE_UPDATEV 4
256: #define DVD_STAGE_ORTHV 5

258: typedef struct {
259:   PetscErrorCode (*old_initV)(struct _dvdDashboard*);
260:   PetscErrorCode (*old_calcPairs)(struct _dvdDashboard*);
261:   PetscErrorCode (*old_improveX)(struct _dvdDashboard*,PetscInt r_s,PetscInt r_e,PetscInt *size_D);
262:   PetscErrorCode (*old_updateV)(struct _dvdDashboard*);
263:   PetscErrorCode (*old_orthV)(struct _dvdDashboard*);
264: } DvdProfiler;

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

268: /*** Other things ****/

272: PetscErrorCode dvd_prof_init()
273: {
274:   PetscErrorCode  ierr;

277:   if (!stages[0]) {
278:     PetscLogStageRegister("Dvd_step_initV",&stages[DVD_STAGE_INITV]);
279:     PetscLogStageRegister("Dvd_step_calcPairs",&stages[DVD_STAGE_CALCPAIRS]);
280:     PetscLogStageRegister("Dvd_step_improveX",&stages[DVD_STAGE_IMPROVEX]);
281:     PetscLogStageRegister("Dvd_step_updateV",&stages[DVD_STAGE_UPDATEV]);
282:     PetscLogStageRegister("Dvd_step_orthV",&stages[DVD_STAGE_ORTHV]);
283:   }
284:   return(0);
285: }

289: PetscErrorCode dvd_initV_prof(dvdDashboard* d)
290: {
291:   DvdProfiler     *p = (DvdProfiler*)d->prof_data;
292:   PetscErrorCode  ierr;

295:   PetscLogStagePush(stages[DVD_STAGE_INITV]);
296:   p->old_initV(d);
297:   PetscLogStagePop();
298:   return(0);
299: }

303: static PetscErrorCode dvd_calcPairs_prof(dvdDashboard* d)
304: {
305:   DvdProfiler    *p = (DvdProfiler*)d->prof_data;

309:   PetscLogStagePush(stages[DVD_STAGE_CALCPAIRS]);
310:   p->old_calcPairs(d);
311:   PetscLogStagePop();
312:   return(0);
313: }

317: static PetscErrorCode dvd_improveX_prof(dvdDashboard *d,PetscInt r_s,PetscInt r_e,PetscInt *size_D)
318: {
319:   DvdProfiler    *p = (DvdProfiler*)d->prof_data;

323:   PetscLogStagePush(stages[DVD_STAGE_IMPROVEX]);
324:   p->old_improveX(d, r_s, r_e, size_D);
325:   PetscLogStagePop();
326:   return(0);
327: }

331: static PetscErrorCode dvd_updateV_prof(dvdDashboard *d)
332: {
333:   DvdProfiler    *p = (DvdProfiler*)d->prof_data;

337:   PetscLogStagePush(stages[DVD_STAGE_UPDATEV]);
338:   p->old_updateV(d);
339:   PetscLogStagePop();
340:   return(0);
341: }

345: PetscErrorCode dvd_profiler(dvdDashboard *d,dvdBlackboard *b)
346: {
348:   DvdProfiler    *p;

351:   /* Setup the step */
352:   if (b->state >= DVD_STATE_CONF) {
353:     PetscFree(d->prof_data);
354:     PetscMalloc(sizeof(DvdProfiler),&p);
355:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(DvdProfiler));
356:     d->prof_data = p;
357:     p->old_initV = d->initV; d->initV = dvd_initV_prof;
358:     p->old_calcPairs = d->calcPairs; d->calcPairs = dvd_calcPairs_prof;
359:     p->old_improveX = d->improveX; d->improveX = dvd_improveX_prof;
360:     p->old_updateV = d->updateV; d->updateV = dvd_updateV_prof;

362:     EPSDavidsonFLAdd(&d->destroyList,dvd_profiler_d);
363:   }
364:   return(0);
365: }

369: static PetscErrorCode dvd_profiler_d(dvdDashboard *d)
370: {
372:   DvdProfiler    *p = (DvdProfiler*)d->prof_data;

375:   /* Free local data */
376:   PetscFree(p);
377:   return(0);
378: }

382: PetscErrorCode dvd_harm_conf(dvdDashboard *d,dvdBlackboard *b,HarmType_t mode,PetscBool fixedTarget,PetscScalar t)
383: {
385:   dvdHarmonic    *dvdh;

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

392:   /* Setup the step */
393:   if (b->state >= DVD_STATE_CONF) {
394:     PetscMalloc(sizeof(dvdHarmonic),&dvdh);
395:     PetscLogObjectMemory((PetscObject)d->eps,sizeof(dvdHarmonic));
396:     dvdh->withTarget = fixedTarget;
397:     dvdh->mode = mode;
398:     if (fixedTarget) dvd_harm_transf(dvdh, t);
399:     d->calcpairs_W_data = dvdh;
400:     d->calcpairs_W = dvd_harm_updateW;
401:     d->calcpairs_proj_trans = dvd_harm_proj;
402:     d->calcpairs_eigs_trans = dvd_harm_eigs_trans;
403:     d->calcpairs_eig_backtrans = dvd_harm_eig_backtrans;

405:     EPSDavidsonFLAdd(&d->destroyList,dvd_harm_d);
406:   }
407:   return(0);
408: }

412: static PetscErrorCode dvd_harm_d(dvdDashboard *d)
413: {

417:   /* Free local data */
418:   PetscFree(d->calcpairs_W_data);
419:   return(0);
420: }

424: static PetscErrorCode dvd_harm_transf(dvdHarmonic *dvdh,PetscScalar t)
425: {
427:   switch (dvdh->mode) {
428:   case DVD_HARM_RR:    /* harmonic RR */
429:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 0.0; dvdh->Pb = -1.0; break;
430:   case DVD_HARM_RRR:   /* relative harmonic RR */
431:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
432:   case DVD_HARM_REIGS: /* rightmost eigenvalues */
433:     dvdh->Wa = 1.0; dvdh->Wb = t;   dvdh->Pa = 1.0; dvdh->Pb = -PetscConj(t);
434:     break;
435:   case DVD_HARM_LEIGS: /* largest eigenvalues */
436:     dvdh->Wa = 0.0; dvdh->Wb = 1.0; dvdh->Pa = 1.0; dvdh->Pb = 0.0; break;
437:   case DVD_HARM_NONE:
438:   default:
439:     SETERRQ(PETSC_COMM_SELF,1, "Harmonic type not supported");
440:   }

442:   /* Check the transformation does not change the sign of the imaginary part */
443: #if !defined(PETSC_USE_COMPLEX)
444:   if (dvdh->Pb*dvdh->Wa - dvdh->Wb*dvdh->Pa < 0.0)
445:     dvdh->Pa*= -1.0, dvdh->Pb*= -1.0;
446: #endif
447:   return(0);
448: }

452: static PetscErrorCode dvd_harm_updateW(dvdDashboard *d)
453: {
454:   dvdHarmonic    *data = (dvdHarmonic*)d->calcpairs_W_data;
456:   PetscInt       l,k,i;
457:   BV             BX = d->BX?d->BX:d->eps->V;

460:   /* Update the target if it is necessary */
461:   if (!data->withTarget) {
462:     dvd_harm_transf(data,d->eigr[0]);
463:   }

465:   /* W(i) <- Wa*AV(i) - Wb*BV(i) */
466:   BVGetActiveColumns(d->eps->V,&l,&k);
467:   if (k != l+d->V_new_s) SETERRQ(PETSC_COMM_SELF,1, "Consistency broken");
468:   BVSetActiveColumns(d->W,l+d->V_new_s,l+d->V_new_e);
469:   BVSetActiveColumns(d->AX,l+d->V_new_s,l+d->V_new_e);
470:   BVSetActiveColumns(BX,l+d->V_new_s,l+d->V_new_e);
471:   BVCopy(d->AX,d->W);
472:   /* Work around bug in BVScale
473:   BVScale(d->W,data->Wa); */
474:   for (i=l+d->V_new_s;i<l+d->V_new_e; ++i) {
475:     BVScaleColumn(d->W,i,data->Wa);
476:   }
477:   BVAXPY(d->W,-data->Wb,BX);
478:   BVSetActiveColumns(d->W,l,k);
479:   BVSetActiveColumns(d->AX,l,k);
480:   BVSetActiveColumns(BX,l,k);
481:   return(0);
482: }

486: static PetscErrorCode dvd_harm_proj(dvdDashboard *d)
487: {
489:   dvdHarmonic    *data = (dvdHarmonic*)d->calcpairs_W_data;
490:   PetscInt       i,j,l0,l,k,ld;
491:   PetscScalar    h,g,*H,*G;

494:   BVGetActiveColumns(d->eps->V,&l0,&k);
495:   l = l0 + d->V_new_s;
496:   k = l0 + d->V_new_e;
497:   MatGetSize(d->H,&ld,NULL);
498:   MatDenseGetArray(d->H,&H);
499:   MatDenseGetArray(d->G,&G);
500:   /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
501:   /* Right part */
502:   for (i=l; i<k; i++) {
503:     for (j=l0; j<k; j++) {
504:       h = H[ld*i+j]; g = G[ld*i+j];
505:       H[ld*i+j] = data->Pa*h - data->Pb*g;
506:       G[ld*i+j] = data->Wa*h - data->Wb*g;
507:     }
508:   }
509:   /* Left part */
510:   for (i=l0; i<l; i++) {
511:     for (j=l; j<k; j++) {
512:       h = H[ld*i+j]; g = G[ld*i+j];
513:       H[ld*i+j] = data->Pa*h - data->Pb*g;
514:       G[ld*i+j] = data->Wa*h - data->Wb*g;
515:     }
516:   }
517:   MatDenseRestoreArray(d->H,&H);
518:   MatDenseRestoreArray(d->G,&G);
519:   return(0);
520: }

524: PetscErrorCode dvd_harm_updateproj(dvdDashboard *d)
525: {
527:   dvdHarmonic    *data = (dvdHarmonic*)d->calcpairs_W_data;
528:   PetscInt       i,j,l,k,ld;
529:   PetscScalar    h,g,*H,*G;

532:   BVGetActiveColumns(d->eps->V,&l,&k);
533:   k = l + d->V_tra_s;
534:   MatGetSize(d->H,&ld,NULL);
535:   MatDenseGetArray(d->H,&H);
536:   MatDenseGetArray(d->G,&G);
537:   /* [H G] <- [Pa*H - Pb*G, Wa*H - Wb*G] */
538:   /* Right part */
539:   for (i=l; i<k; i++) {
540:     for (j=0; j<l; j++) {
541:       h = H[ld*i+j]; g = G[ld*i+j];
542:       H[ld*i+j] = data->Pa*h - data->Pb*g;
543:       G[ld*i+j] = data->Wa*h - data->Wb*g;
544:     }
545:   }
546:   /* Lower triangular part */
547:   for (i=0; i<l; i++) {
548:     for (j=l; j<k; j++) {
549:       h = H[ld*i+j]; g = G[ld*i+j];
550:       H[ld*i+j] = data->Pa*h - data->Pb*g;
551:       G[ld*i+j] = data->Wa*h - data->Wb*g;
552:     }
553:   }
554:   MatDenseRestoreArray(d->H,&H);
555:   MatDenseRestoreArray(d->G,&G);
556:   return(0);
557: }

561: static PetscErrorCode dvd_harm_backtrans(dvdHarmonic *data,PetscScalar *ar,PetscScalar *ai)
562: {
563:   PetscScalar xr;
564: #if !defined(PETSC_USE_COMPLEX)
565:   PetscScalar xi, k;
566: #endif

570:   xr = *ar;
571: #if !defined(PETSC_USE_COMPLEX)
573:   xi = *ai;

575:   if (xi != 0.0) {
576:     k = (data->Pa - data->Wa*xr)*(data->Pa - data->Wa*xr) +
577:         data->Wa*data->Wa*xi*xi;
578:     *ar = (data->Pb*data->Pa - (data->Pb*data->Wa + data->Wb*data->Pa)*xr +
579:            data->Wb*data->Wa*(xr*xr + xi*xi))/k;
580:     *ai = (data->Pb*data->Wa - data->Wb*data->Pa)*xi/k;
581:   } else
582: #endif
583:     *ar = (data->Pb - data->Wb*xr) / (data->Pa - data->Wa*xr);
584:   return(0);
585: }

589: static PetscErrorCode dvd_harm_eig_backtrans(dvdDashboard *d,PetscScalar ar,PetscScalar ai,PetscScalar *br,PetscScalar *bi)
590: {
591:   dvdHarmonic    *data = (dvdHarmonic*)d->calcpairs_W_data;

595:   dvd_harm_backtrans(data,&ar,&ai);
596:   *br = ar;
597:   *bi = ai;
598:   return(0);
599: }

603: static PetscErrorCode dvd_harm_eigs_trans(dvdDashboard *d)
604: {
605:   dvdHarmonic    *data = (dvdHarmonic*)d->calcpairs_W_data;
606:   PetscInt       i,l,k;

610:   BVGetActiveColumns(d->eps->V,&l,&k);
611:   for (i=0;i<k-l;i++) {
612:     dvd_harm_backtrans(data,&d->eigr[i],&d->eigi[i]);
613:   }
614:   return(0);
615: }

619: /* H = [H              Y(old)'*X(new);
620:         Y(new)'*X(old) Y(new)'*X(new) ],
621:      being old=0:l-1, new=l:k-1 */
622: PetscErrorCode BVMultS(BV X,BV Y,PetscScalar *H,PetscInt ldh)
623: {
625:   PetscInt       j,lx,ly,kx,ky;
626:   PetscScalar    *array;
627:   Mat            M;

630:   BVGetActiveColumns(X,&lx,&kx);
631:   BVGetActiveColumns(Y,&ly,&ky);
632:   MatCreateSeqDense(PETSC_COMM_SELF,ky,kx,NULL,&M);
633:   BVMatProject(X,NULL,Y,M);
634:   MatDenseGetArray(M,&array);
635:   /* upper part */
636:   for (j=lx;j<kx;j++) {
637:     PetscMemcpy(&H[ldh*j],&array[j*ky],ly*sizeof(PetscScalar));
638:   }
639:   /* lower part */
640:   for (j=0;j<kx;j++) {
641:     PetscMemcpy(&H[ldh*j+ly],&array[j*ky+ly],(ky-ly)*sizeof(PetscScalar));
642:   }
643:   MatDenseRestoreArray(M,&array);
644:   MatDestroy(&M);
645:   return(0);
646: }

650: /*
651:    SlepcMatDenseCopy - Copy a submatrix from A to B.

653:    Not Collective

655:    Input Parameters:
656: +  A    - source seq dense matrix
657: .  Ar0  - first row to copy from A
658: .  Ac0  - first column to copy from A
659: .  Br0  - first row to copy on B
660: .  Bc0  - first column to copy on B
661: .  rows - number of rows to copy
662: -  cols - number of columns to copy

664:    Level: advanced
665: */
666: PetscErrorCode SlepcMatDenseCopy(Mat A,PetscInt Ar0,PetscInt Ac0,Mat B,PetscInt Br0,PetscInt Bc0,PetscInt rows,PetscInt cols)
667: {
669:   PetscInt       i,n,m,ldA,ldB;
670:   PetscScalar    *pA,*pB;

673:   if (!rows || !cols) return(0);
674:   MatGetSize(A,&m,&n); ldA=m;
675:   if (Ar0<0 || Ar0>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial row in A");
676:   if (Ac0<0 || Ac0>=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in A");
677:   if (Ar0+rows>m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid number of rows");
678:   if (Ac0+cols>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid number of columns");
679:   MatGetSize(B,&m,&n); ldB=m;
680:   if (Br0<0 || Br0>=m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial row in B");
681:   if (Bc0<0 || Bc0>=n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid initial column in B");
682:   if (Br0+rows>m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid number of rows");
683:   if (Bc0+cols>n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid number of columns");
684:   MatDenseGetArray(A,&pA);
685:   MatDenseGetArray(B,&pB);
686:   for (i=0;i<cols;i++) {
687:     PetscMemcpy(&pB[ldB*(Bc0+i)+Br0],&pA[ldA*(Ac0+i)+Ar0],sizeof(PetscScalar)*rows);
688:   }
689:   MatDenseRestoreArray(A,&pA);
690:   MatDenseRestoreArray(B,&pB);
691:   return(0);
692: }