Actual source code: ispai.c

  1: /* $Id: ispai.c,v 1.25 2001/04/10 19:36:18 bsmith Exp $*/

  3: /* 
  4:    3/99 Modified by Stephen Barnard to support SPAI version 3.0 
  5: */

  7: /*
  8:       Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner
  9:    Code written by Stephen Barnard.

 11:       Note: there is some BAD memory bleeding below!

 13:       This code needs work

 15:    1) get rid of all memory bleeding
 16:    2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix
 17:       rather than having the sp flag for PC_SPAI

 19: */

 21: #include "src/sles/pc/pcimpl.h"        /*I "petscpc.h" I*/

 23: /*
 24:     These are the SPAI include files
 25: */
 26: EXTERN_C_BEGIN
 27: #include "spai.h"
 28: #include "matrix.h"
 29: #include "read_mm_matrix.h"
 30: EXTERN_C_END

 32: EXTERN int ConvertMatToMatrix(Mat,Mat,matrix**);
 33: EXTERN int ConvertMatrixToMat(matrix *,Mat *);
 34: EXTERN int MM_to_PETSC(char *,char *,char *);

 36: typedef struct {

 38:   matrix *B;              /* matrix in SPAI format */
 39:   matrix *BT;             /* transpose of matrix in SPAI format */
 40:   matrix *M;              /* the approximate inverse in SPAI format */

 42:   Mat    PM;              /* the approximate inverse PETSc format */

 44:   double epsilon;         /* tolerance */
 45:   int    nbsteps;         /* max number of "improvement" steps per line */
 46:   int    max;             /* max dimensions of is_I, q, etc. */
 47:   int    maxnew;          /* max number of new entries per step */
 48:   int    block_size;      /* constant block size */
 49:   int    cache_size;      /* one of (1,2,3,4,5,6) indicting size of cache */
 50:   int    verbose;         /* SPAI prints timing and statistics */

 52:   int    sp;              /* symmetric nonzero pattern */

 54: } PC_SPAI;

 56: /**********************************************************************/

 58: static int PCSetUp_SPAI(PC pc)
 59: {
 60:   PC_SPAI *ispai = (PC_SPAI*)pc->data;
 61:   int      ierr;
 62:   Mat      AT;
 63:   Mat      PBT;
 64:   MPI_Comm comm;

 67:   PetscObjectGetComm((PetscObject)pc->pmat,&comm);

 69:   init_SPAI();

 71:   if (ispai->sp) {
 72:     ConvertMatToMatrix(pc->pmat,pc->pmat,&ispai->B);
 73:   } else {
 74:     /* Use the transpose to get the column nonzero structure. */
 75:     MatTranspose(pc->pmat,&AT);
 76:     ConvertMatToMatrix(pc->pmat,AT,&ispai->B);
 77:     MatDestroy(AT);
 78:   }

 80:   /* Destroy the transpose */
 81:   /* Don't know how to do it. PETSc developers? */
 82: 
 83:   /* construct SPAI preconditioner */
 84:   /* FILE *messages */     /* file for warning messages */
 85:   /* double epsilon */     /* tolerance */
 86:   /* int nbsteps */        /* max number of "improvement" steps per line */
 87:   /* int max */            /* max dimensions of is_I, q, etc. */
 88:   /* int maxnew */         /* max number of new entries per step */
 89:   /* int block_size */     /* block_size == 1 specifies scalar elments
 90:                               block_size == n specifies nxn constant-block elements
 91:                               block_size == 0 specifies variable-block elements */
 92:   /* int cache_size */     /* one of (1,2,3,4,5,6) indicting size of cache */
 93:                            /* cache_size == 0 indicates no caching */
 94:   /* int    verbose    */  /* verbose == 0 specifies that SPAI is silent
 95:                               verbose == 1 prints timing and matrix statistics */

 97:   ispai->M = bspai(ispai->B,
 98:                    stdout,
 99:                    ispai->epsilon,
100:                    ispai->nbsteps,
101:                    ispai->max,
102:                    ispai->maxnew,
103:                    ispai->block_size,
104:                    ispai->cache_size,
105:                    ispai->verbose);

107:   if (!ispai->M) SETERRQ(1,"Unable to create SPAI preconditioner");

109:   ConvertMatrixToMat(ispai->M,&ispai->PM);

111:   /* free the SPAI matrices */
112:   sp_free_matrix(ispai->B);
113:   sp_free_matrix(ispai->M);

115:   return(0);
116: }

118: /**********************************************************************/

120: static int PCApply_SPAI(PC pc,Vec x,Vec y)
121: {
122:   PC_SPAI *ispai = (PC_SPAI*)pc->data;
123:   int      ierr;

126:   /* Now using PETSc's multiply */
127:   MatMult(ispai->PM,x,y);
128:   return(0);
129: }

131: /**********************************************************************/

133: static int PCDestroy_SPAI(PC pc)
134: {
135:   int     ierr;
136:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

139:   MatDestroy(ispai->PM);
140:   PetscFree(ispai);
141:   return(0);
142: }

144: /**********************************************************************/

146: static int PCView_SPAI(PC pc,PetscViewer viewer)
147: {
148:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
149:   int        ierr;
150:   PetscTruth isascii;

153:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
154:   if (isascii) {
155:     PetscViewerASCIIPrintf(viewer,"    SPAI preconditionern");
156:     PetscViewerASCIIPrintf(viewer,"    epsilon %gn",   ispai->epsilon);
157:     PetscViewerASCIIPrintf(viewer,"    nbsteps %dn",   ispai->nbsteps);
158:     PetscViewerASCIIPrintf(viewer,"    max %dn",       ispai->max);
159:     PetscViewerASCIIPrintf(viewer,"    maxnew %dn",    ispai->maxnew);
160:     PetscViewerASCIIPrintf(viewer,"    block_size %dn",ispai->block_size);
161:     PetscViewerASCIIPrintf(viewer,"    cache_size %dn",ispai->cache_size);
162:     PetscViewerASCIIPrintf(viewer,"    verbose %dn",   ispai->verbose);
163:     PetscViewerASCIIPrintf(viewer,"    sp %dn",        ispai->sp);

165:   }
166:   return(0);
167: }

169: EXTERN_C_BEGIN
170: int PCSPAISetEpsilon_SPAI(PC pc,double epsilon)
171: {
172:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
174:   ispai->epsilon = epsilon;
175:   return(0);
176: }
177: EXTERN_C_END
178: 
179: /**********************************************************************/

181: EXTERN_C_BEGIN
182: int PCSPAISetNBSteps_SPAI(PC pc,int nbsteps)
183: {
184:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
186:   ispai->nbsteps = nbsteps;
187:   return(0);
188: }
189: EXTERN_C_END

191: /**********************************************************************/

193: /* added 1/7/99 g.h. */
194: EXTERN_C_BEGIN
195: int PCSPAISetMax_SPAI(PC pc,int max)
196: {
197:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
199:   ispai->max = max;
200:   return(0);
201: }
202: EXTERN_C_END

204: /**********************************************************************/

206: EXTERN_C_BEGIN
207: int PCSPAISetMaxNew_SPAI(PC pc,int maxnew)
208: {
209:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
211:   ispai->maxnew = maxnew;
212:   return(0);
213: }
214: EXTERN_C_END

216: /**********************************************************************/

218: EXTERN_C_BEGIN
219: int PCSPAISetBlockSize_SPAI(PC pc,int block_size)
220: {
221:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
223:   ispai->block_size = block_size;
224:   return(0);
225: }
226: EXTERN_C_END

228: /**********************************************************************/

230: EXTERN_C_BEGIN
231: int PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
232: {
233:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
235:   ispai->cache_size = cache_size;
236:   return(0);
237: }
238: EXTERN_C_END

240: /**********************************************************************/

242: EXTERN_C_BEGIN
243: int PCSPAISetVerbose_SPAI(PC pc,int verbose)
244: {
245:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
247:   ispai->verbose = verbose;
248:   return(0);
249: }
250: EXTERN_C_END

252: /**********************************************************************/

254: EXTERN_C_BEGIN
255: int PCSPAISetSp_SPAI(PC pc,int sp)
256: {
257:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
259:   ispai->sp = sp;
260:   return(0);
261: }
262: EXTERN_C_END

264: /* -------------------------------------------------------------------*/

266: int PCSPAISetEpsilon(PC pc,double epsilon)
267: {
268:   int ierr,(*f)(PC,double);
270:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetEpsilon_C",(void (**)())&f);
271:   if (f) {
272:     (*f)(pc,epsilon);
273:   }
274:   return(0);
275: }
276: 
277: /**********************************************************************/

279: int PCSPAISetNBSteps(PC pc,int nbsteps)
280: {
281:   int ierr,(*f)(PC,int);
283:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetNBSteps_C",(void (**)())&f);
284:   if (f) {
285:     (*f)(pc,nbsteps);
286:   }
287:   return(0);
288: }

290: /**********************************************************************/

292: /* added 1/7/99 g.h. */
293: int PCSPAISetMax(PC pc,int max)
294: {
295:   int ierr,(*f)(PC,int);
297:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMax_C",(void (**)())&f);
298:   if (f) {
299:     (*f)(pc,max);
300:   }
301:   return(0);
302: }

304: /**********************************************************************/

306: int PCSPAISetMaxNew(PC pc,int maxnew)
307: {
308:   int ierr,(*f)(PC,int);
310:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMaxNew_C",(void (**)())&f);
311:   if (f) {
312:     (*f)(pc,maxnew);
313:   }
314:   return(0);
315: }

317: /**********************************************************************/

319: int PCSPAISetBlockSize(PC pc,int block_size)
320: {
321:   int ierr,(*f)(PC,int);
323:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetBlockSize_C",(void (**)())&f);
324:   if (f) {
325:     (*f)(pc,block_size);
326:   }
327:   return(0);
328: }

330: /**********************************************************************/

332: int PCSPAISetCacheSize(PC pc,int cache_size)
333: {
334:   int ierr,(*f)(PC,int);
336:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetCacheSize_C",(void (**)())&f);
337:   if (f) {
338:     (*f)(pc,cache_size);
339:   }
340:   return(0);
341: }

343: /**********************************************************************/

345: int PCSPAISetVerbose(PC pc,int verbose)
346: {
347:   int ierr,(*f)(PC,int);
349:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetVerbose_C",(void (**)())&f);
350:   if (f) {
351:     (*f)(pc,verbose);
352:   }
353:   return(0);
354: }

356: /**********************************************************************/

358: int PCSPAISetSp(PC pc,int sp)
359: {
360:   int ierr,(*f)(PC,int);
362:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetSp_C",(void (**)())&f);
363:   if (f) {
364:     (*f)(pc,sp);
365:   }
366:   return(0);
367: }

369: /**********************************************************************/

371: /**********************************************************************/

373: static int PCSetFromOptions_SPAI(PC pc)
374: {
375:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
376:   int        ierr,nbsteps,max,maxnew,block_size,cache_size,verbose,sp;
377:   double     epsilon;
378:   PetscTruth flg;

381:   PetscOptionsHead("SPAI options");
382:     PetscOptionsDouble("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon,&flg);
383:     if (flg) {
384:       PCSPAISetEpsilon(pc,epsilon);
385:     }
386:     PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps,&flg);
387:     if (flg) {
388:       PCSPAISetNBSteps(pc,nbsteps);
389:     }
390:     /* added 1/7/99 g.h. */
391:     PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max,&flg);
392:     if (flg) {
393:       PCSPAISetMax(pc,max);
394:     }
395:     PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew,&flg);
396:     if (flg) {
397:       PCSPAISetMaxNew(pc,maxnew);
398:     }
399:     PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size,&flg);
400:     if (flg) {
401:       PCSPAISetBlockSize(pc,block_size);
402:     }
403:     PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);
404:     if (flg) {
405:       PCSPAISetCacheSize(pc,cache_size);
406:     }
407:     PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);
408:     if (flg) {
409:       PCSPAISetVerbose(pc,verbose);
410:     }
411:     PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);
412:     if (flg) {
413:       PCSPAISetSp(pc,sp);
414:     }
415:   PetscOptionsTail();
416:   return(0);
417: }

419: /**********************************************************************/

421: EXTERN_C_BEGIN
422: /*
423:    PCCreate_SPAI - Creates the preconditioner context for the SPAI 
424:                    preconditioner written by Stephen Barnard.

426: */
427: int PCCreate_SPAI(PC pc)
428: {
429:   PC_SPAI *ispai;
430:   int     ierr;

433:   ierr               = PetscNew(PC_SPAI,&ispai);
434:   pc->data           = (void*)ispai;

436:   pc->ops->destroy         = PCDestroy_SPAI;
437:   pc->ops->apply           = PCApply_SPAI;
438:   pc->ops->applyrichardson = 0;
439:   pc->ops->setup           = PCSetUp_SPAI;
440:   pc->ops->view            = PCView_SPAI;
441:   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;

443:   pc->name          = 0;
444:   ispai->epsilon    = .4;
445:   ispai->nbsteps    = 5;
446:   ispai->max        = 5000;
447:   ispai->maxnew     = 5;
448:   ispai->block_size = 1;
449:   ispai->cache_size = 5;
450:   ispai->verbose    = 0;

452:   ispai->sp         = 1;

454:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C",
455:                     "PCSPAISetEpsilon_SPAI",
456:                      PCSPAISetEpsilon_SPAI);
457:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C",
458:                     "PCSPAISetNBSteps_SPAI",
459:                      PCSPAISetNBSteps_SPAI);
460:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C",
461:                     "PCSPAISetMax_SPAI",
462:                      PCSPAISetMax_SPAI);
463:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC",
464:                     "PCSPAISetMaxNew_SPAI",
465:                      PCSPAISetMaxNew_SPAI);
466:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C",
467:                     "PCSPAISetBlockSize_SPAI",
468:                      PCSPAISetBlockSize_SPAI);
469:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C",
470:                     "PCSPAISetCacheSize_SPAI",
471:                      PCSPAISetCacheSize_SPAI);
472:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C",
473:                     "PCSPAISetVerbose_SPAI",
474:                      PCSPAISetVerbose_SPAI);
475:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C",
476:                     "PCSPAISetSp_SPAI",
477:                      PCSPAISetSp_SPAI);

479:   return(0);
480: }
481: EXTERN_C_END

483: /**********************************************************************/

485: /*
486:    Converts from a PETSc matrix to an SPAI matrix 
487: */
488: int ConvertMatToMatrix(Mat A,Mat AT,matrix **B)
489: {
490:   matrix   *M;
491:   int      i,j,col;
492:   int      row_indx;
493:   int      len,pe,local_indx,start_indx;
494:   int      *mapping;
495:   int      ierr,*cols;
496:   double   *vals;
497:   int      *num_ptr,n,mnl,nnl,rank,size,nz,rstart,rend;
498:   MPI_Comm comm;
499:   struct compressed_lines *rows;

502:   PetscObjectGetComm((PetscObject)A,&comm);
503: 
504:   MPI_Comm_size(comm,&size);
505:   MPI_Comm_rank(comm,&rank);
506:   MatGetSize(A,&n,&n);
507:   MatGetLocalSize(A,&mnl,&nnl);

509:   MPI_Barrier(PETSC_COMM_WORLD);

511:   M = new_matrix((void *)comm);
512: 
513:   M->n = n;
514:   M->bs = 1;
515:   M->max_block_size = 1;

517:   M->mnls          = malloc(sizeof(int)*size);
518:   M->start_indices = malloc(sizeof(int)*size);
519:   M->pe            = malloc(sizeof(int)*n);
520:   M->block_sizes   = malloc(sizeof(int)*n);
521:   for (i=0; i<n; i++) M->block_sizes[i] = 1;

523:   MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,PETSC_COMM_WORLD);

525:   M->start_indices[0] = 0;
526:   for (i=1; i<size; i++) {
527:     M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
528:   }

530:   M->mnl = M->mnls[M->myid];
531:   M->my_start_index = M->start_indices[M->myid];

533:   for (i=0; i<size; i++) {
534:     start_indx = M->start_indices[i];
535:     for (j=0; j<M->mnls[i]; j++)
536:       M->pe[start_indx+j] = i;
537:   }

539:   if (AT) {
540:     M->lines = new_compressed_lines(M->mnls[rank],1);
541:   } else {
542:     M->lines = new_compressed_lines(M->mnls[rank],0);
543:   }

545:   rows     = M->lines;

547:   /* Determine the mapping from global indices to pointers */
548:   ierr       = PetscMalloc(M->n*sizeof(int),&mapping);
549:   pe         = 0;
550:   local_indx = 0;
551:   for (i=0; i<M->n; i++) {
552:     if (local_indx >= M->mnls[pe]) {
553:       pe++;
554:       local_indx = 0;
555:     }
556:     mapping[i] = local_indx + M->start_indices[pe];
557:     local_indx++;
558:   }


561:   PetscMalloc(mnl*sizeof(int),&num_ptr);

563:   /*********************************************************/
564:   /************** Set up the row structure *****************/
565:   /*********************************************************/

567:   /* count number of nonzeros in every row */
568:   MatGetOwnershipRange(A,&rstart,&rend);
569:   for (i=rstart; i<rend; i++) {
570:     MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
571:     MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
572:   }

574:   /* allocate buffers */
575:   len = 0;
576:   for (i=0; i<n; i++) {
577:     if (len < num_ptr[i]) len = num_ptr[i];
578:   }

580:   for (i=rstart; i<rend; i++) {
581:     row_indx             = i-rstart;
582:     len                  = num_ptr[row_indx];
583:     rows->ptrs[row_indx] = malloc(len*sizeof(int));
584:     rows->A[row_indx]    = malloc(len*sizeof(double));
585:   }

587:   /* copy the matrix */
588:   for (i=rstart; i<rend; i++) {
589:     row_indx = i - rstart;
590:     ierr     = MatGetRow(A,i,&nz,&cols,&vals);
591:     for (j=0; j<nz; j++) {
592:       col = cols[j];
593:       len = rows->len[row_indx]++;
594:       rows->ptrs[row_indx][len] = mapping[col];
595:       rows->A[row_indx][len]    = vals[j];
596:     }
597:     rows->slen[row_indx] = rows->len[row_indx];
598:     MatRestoreRow(A,i,&nz,&cols,&vals);
599:   }


602:   /************************************************************/
603:   /************** Set up the column structure *****************/
604:   /*********************************************************/

606:   if (AT) {

608:     /* count number of nonzeros in every column */
609:     for (i=rstart; i<rend; i++) {
610:       MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
611:       MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
612:     }

614:     /* allocate buffers */
615:     len = 0;
616:     for (i=0; i<n; i++) {
617:       if (len < num_ptr[i]) len = num_ptr[i];
618:     }

620:     for (i=rstart; i<rend; i++) {
621:       row_indx = i-rstart;
622:       len      = num_ptr[row_indx];
623:       rows->rptrs[row_indx]     = malloc(len*sizeof(int));
624:     }

626:     /* copy the matrix (i.e., the structure) */
627:     for (i=rstart; i<rend; i++) {
628:       row_indx = i - rstart;
629:       ierr     = MatGetRow(AT,i,&nz,&cols,&vals);
630:       for (j=0; j<nz; j++) {
631:         col = cols[j];
632:         len = rows->rlen[row_indx]++;
633:         rows->rptrs[row_indx][len] = mapping[col];
634:       }
635:       ierr     = MatRestoreRow(AT,i,&nz,&cols,&vals);
636:     }
637:   }

639:   PetscFree(num_ptr);
640:   PetscFree(mapping);

642:   order_pointers(M);
643:   M->maxnz = calc_maxnz(M);

645:   *B = M;

647:   return(0);
648: }

650: /**********************************************************************/

652: /*
653:    Converts from an SPAI matrix B  to a PETSc matrix PB.
654:    This assumes that the the SPAI matrix B is stored in
655:    COMPRESSED-ROW format.
656: */
657: int ConvertMatrixToMat(matrix *B,Mat *PB)
658: {
659:   int    size,rank;
660:   int    ierr;
661:   int    m,n,M,N;
662:   int    d_nz,o_nz;
663:   int    *d_nnz,*o_nnz;
664:   int    i,k,global_row,global_col,first_diag_col,last_diag_col;
665:   Scalar val;

668:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
669:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
670: 
671:   m = n = B->mnls[rank];
672:   d_nz = o_nz = 0;

674:   /* Determine preallocation for MatCreateMPIAIJ */
675:   PetscMalloc(m*sizeof(int),&d_nnz);
676:   PetscMalloc(m*sizeof(int),&o_nnz);
677:   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
678:   first_diag_col = B->start_indices[rank];
679:   last_diag_col = first_diag_col + B->mnls[rank];
680:   for (i=0; i<B->mnls[rank]; i++) {
681:     for (k=0; k<B->lines->len[i]; k++) {
682:       global_col = B->lines->ptrs[i][k];
683:       if ((global_col >= first_diag_col) && (global_col <= last_diag_col))
684:         d_nnz[i]++;
685:       else
686:         o_nnz[i]++;
687:     }
688:   }

690:   M = N = B->n;
691:   MatCreateMPIAIJ(PETSC_COMM_WORLD,m,n,M,N,d_nz,d_nnz,o_nz,o_nnz,PB);

693:   for (i=0; i<B->mnls[rank]; i++) {
694:     global_row = B->start_indices[rank]+i;
695:     for (k=0; k<B->lines->len[i]; k++) {
696:       global_col = B->lines->ptrs[i][k];
697:       val = B->lines->A[i][k];
698:       MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);
699:     }
700:   }

702:   PetscFree(d_nnz);
703:   PetscFree(o_nnz);

705:   MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);
706:   MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);

708:   return(0);
709: }

711: /**********************************************************************/

713: /*
714:    Converts from an SPAI vector v  to a PETSc vec Pv.
715: */
716: int ConvertVectorToVec(vector *v,Vec *Pv)
717: {
718:   int size,rank,ierr,m,M,i,*mnls,*start_indices,*global_indices;
719: 
721:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
722:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
723: 
724:   m = v->mnl;
725:   M = v->n;
726: 
727: 
728:   VecCreateMPI(PETSC_COMM_WORLD,m,M,Pv);

730:   PetscMalloc(size*sizeof(int),&mnls);
731:   MPI_Allgather((void*)&v->mnl,1,MPI_INT,(void*)mnls,1,MPI_INT,PETSC_COMM_WORLD);
732: 
733:   PetscMalloc(size*sizeof(int),&start_indices);
734:   start_indices[0] = 0;
735:   for (i=1; i<size; i++)
736:     start_indices[i] = start_indices[i-1] +mnls[i-1];
737: 
738:   PetscMalloc(v->mnl*sizeof(int),&global_indices);
739:   for (i=0; i<v->mnl; i++)
740:     global_indices[i] = start_indices[rank] + i;

742:   PetscFree(mnls);
743:   PetscFree(start_indices);
744: 
745:   VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);

747:   PetscFree(global_indices);

749:   VecAssemblyBegin(*Pv);
750:   VecAssemblyEnd(*Pv);
751: 
752:   return(0);
753: }

755: /**********************************************************************/

757: /*

759:   Reads a matrix and a RHS vector in Matrix Market format and writes them
760:   in PETSc binary format.

762:   !!!! Warning!!!!! PETSc supports only serial execution of this routine.
763:   
764:   f0 <input_file> : matrix in Matrix Market format
765:   f1 <input_file> : vector in Matrix Market array format
766:   f2 <output_file> : matrix and vector in PETSc format

768: */
769: int MM_to_PETSC(char *f0,char *f1,char *f2)
770: {
771:   Mat         A_PETSC;          /* matrix */
772:   Vec         b_PETSC;          /* RHS */
773:   PetscViewer fd;               /* viewer */
774:   int         ierr;
775:   matrix      *A_spai;
776:   vector      *b_spai;

779:   A_spai = read_mm_matrix(f0,1,1,1,0,0,0,NULL);
780:   b_spai = read_rhs_for_matrix(f1,A_spai);

782:   ConvertMatrixToMat(A_spai,&A_PETSC);
783:   ConvertVectorToVec(b_spai,&b_PETSC);

785:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,f1,PETSC_BINARY_CREATE,&fd);
786:   MatView(A_PETSC,fd);
787:   VecView(b_PETSC,fd);

789:   return(0);
790: }