Actual source code: ispai.c

  1: /* $Id: ispai.c,v 1.28 2001/08/07 03:03:40 balay 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

 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: #undef __FUNCT__  
 60: static int PCSetUp_SPAI(PC pc)
 61: {
 62:   PC_SPAI *ispai = (PC_SPAI*)pc->data;
 63:   int      ierr;
 64:   Mat      AT;
 65:   MPI_Comm comm;

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

 70:   init_SPAI();

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

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

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

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

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

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

116:   return(0);
117: }

119: /**********************************************************************/

121: #undef __FUNCT__  
123: static int PCApply_SPAI(PC pc,Vec xx,Vec y)
124: {
125:   PC_SPAI *ispai = (PC_SPAI*)pc->data;
126:   int      ierr;

129:   /* Now using PETSc's multiply */
130:   MatMult(ispai->PM,xx,y);
131:   return(0);
132: }

134: /**********************************************************************/

136: #undef __FUNCT__  
138: static int PCDestroy_SPAI(PC pc)
139: {
140:   int     ierr;
141:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

144:   if (ispai->PM) {MatDestroy(ispai->PM);}
145:   PetscFree(ispai);
146:   return(0);
147: }

149: /**********************************************************************/

151: #undef __FUNCT__  
153: static int PCView_SPAI(PC pc,PetscViewer viewer)
154: {
155:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
156:   int        ierr;
157:   PetscTruth isascii;

160:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
161:   if (isascii) {
162:     PetscViewerASCIIPrintf(viewer,"    SPAI preconditionern");
163:     PetscViewerASCIIPrintf(viewer,"    epsilon %gn",   ispai->epsilon);
164:     PetscViewerASCIIPrintf(viewer,"    nbsteps %dn",   ispai->nbsteps);
165:     PetscViewerASCIIPrintf(viewer,"    max %dn",       ispai->max);
166:     PetscViewerASCIIPrintf(viewer,"    maxnew %dn",    ispai->maxnew);
167:     PetscViewerASCIIPrintf(viewer,"    block_size %dn",ispai->block_size);
168:     PetscViewerASCIIPrintf(viewer,"    cache_size %dn",ispai->cache_size);
169:     PetscViewerASCIIPrintf(viewer,"    verbose %dn",   ispai->verbose);
170:     PetscViewerASCIIPrintf(viewer,"    sp %dn",        ispai->sp);

172:   }
173:   return(0);
174: }

176: EXTERN_C_BEGIN
177: #undef __FUNCT__  
179: int PCSPAISetEpsilon_SPAI(PC pc,double epsilon1)
180: {
181:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
183:   ispai->epsilon = epsilon1;
184:   return(0);
185: }
186: EXTERN_C_END
187: 
188: /**********************************************************************/

190: EXTERN_C_BEGIN
191: #undef __FUNCT__  
193: int PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1)
194: {
195:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
197:   ispai->nbsteps = nbsteps1;
198:   return(0);
199: }
200: EXTERN_C_END

202: /**********************************************************************/

204: /* added 1/7/99 g.h. */
205: EXTERN_C_BEGIN
206: #undef __FUNCT__  
208: int PCSPAISetMax_SPAI(PC pc,int max1)
209: {
210:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
212:   ispai->max = max1;
213:   return(0);
214: }
215: EXTERN_C_END

217: /**********************************************************************/

219: EXTERN_C_BEGIN
220: #undef __FUNCT__  
222: int PCSPAISetMaxNew_SPAI(PC pc,int maxnew1)
223: {
224:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
226:   ispai->maxnew = maxnew1;
227:   return(0);
228: }
229: EXTERN_C_END

231: /**********************************************************************/

233: EXTERN_C_BEGIN
234: #undef __FUNCT__  
236: int PCSPAISetBlockSize_SPAI(PC pc,int block_size1)
237: {
238:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
240:   ispai->block_size = block_size1;
241:   return(0);
242: }
243: EXTERN_C_END

245: /**********************************************************************/

247: EXTERN_C_BEGIN
248: #undef __FUNCT__  
250: int PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
251: {
252:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
254:   ispai->cache_size = cache_size;
255:   return(0);
256: }
257: EXTERN_C_END

259: /**********************************************************************/

261: EXTERN_C_BEGIN
262: #undef __FUNCT__  
264: int PCSPAISetVerbose_SPAI(PC pc,int verbose)
265: {
266:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
268:   ispai->verbose = verbose;
269:   return(0);
270: }
271: EXTERN_C_END

273: /**********************************************************************/

275: EXTERN_C_BEGIN
276: #undef __FUNCT__  
278: int PCSPAISetSp_SPAI(PC pc,int sp)
279: {
280:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
282:   ispai->sp = sp;
283:   return(0);
284: }
285: EXTERN_C_END

287: /* -------------------------------------------------------------------*/

289: #undef __FUNCT__  
291: int PCSPAISetEpsilon(PC pc,double epsilon1)
292: {
293:   int ierr,(*f)(PC,double);
295:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetEpsilon_C",(void (**)(void))&f);
296:   if (f) {
297:     (*f)(pc,epsilon1);
298:   }
299:   return(0);
300: }
301: 
302: /**********************************************************************/

304: #undef __FUNCT__  
306: int PCSPAISetNBSteps(PC pc,int nbsteps1)
307: {
308:   int ierr,(*f)(PC,int);
310:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetNBSteps_C",(void (**)(void))&f);
311:   if (f) {
312:     (*f)(pc,nbsteps1);
313:   }
314:   return(0);
315: }

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

319: /* added 1/7/99 g.h. */
320: #undef __FUNCT__  
322: int PCSPAISetMax(PC pc,int max1)
323: {
324:   int ierr,(*f)(PC,int);
326:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMax_C",(void (**)(void))&f);
327:   if (f) {
328:     (*f)(pc,max1);
329:   }
330:   return(0);
331: }

333: /**********************************************************************/

335: #undef __FUNCT__  
337: int PCSPAISetMaxNew(PC pc,int maxnew1)
338: {
339:   int ierr,(*f)(PC,int);
341:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMaxNew_C",(void (**)(void))&f);
342:   if (f) {
343:     (*f)(pc,maxnew1);
344:   }
345:   return(0);
346: }

348: /**********************************************************************/

350: #undef __FUNCT__  
352: int PCSPAISetBlockSize(PC pc,int block_size1)
353: {
354:   int ierr,(*f)(PC,int);
356:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetBlockSize_C",(void (**)(void))&f);
357:   if (f) {
358:     (*f)(pc,block_size1);
359:   }
360:   return(0);
361: }

363: /**********************************************************************/

365: #undef __FUNCT__  
367: int PCSPAISetCacheSize(PC pc,int cache_size)
368: {
369:   int ierr,(*f)(PC,int);
371:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetCacheSize_C",(void (**)(void))&f);
372:   if (f) {
373:     (*f)(pc,cache_size);
374:   }
375:   return(0);
376: }

378: /**********************************************************************/

380: #undef __FUNCT__  
382: int PCSPAISetVerbose(PC pc,int verbose)
383: {
384:   int ierr,(*f)(PC,int);
386:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetVerbose_C",(void (**)(void))&f);
387:   if (f) {
388:     (*f)(pc,verbose);
389:   }
390:   return(0);
391: }

393: /**********************************************************************/

395: #undef __FUNCT__  
397: int PCSPAISetSp(PC pc,int sp)
398: {
399:   int ierr,(*f)(PC,int);
401:   PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetSp_C",(void (**)(void))&f);
402:   if (f) {
403:     (*f)(pc,sp);
404:   }
405:   return(0);
406: }

408: /**********************************************************************/

410: /**********************************************************************/

412: #undef __FUNCT__  
414: static int PCSetFromOptions_SPAI(PC pc)
415: {
416:   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
417:   int        ierr,nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
418:   double     epsilon1;
419:   PetscTruth flg;

422:   PetscOptionsHead("SPAI options");
423:     PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);
424:     if (flg) {
425:       PCSPAISetEpsilon(pc,epsilon1);
426:     }
427:     PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);
428:     if (flg) {
429:       PCSPAISetNBSteps(pc,nbsteps1);
430:     }
431:     /* added 1/7/99 g.h. */
432:     PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);
433:     if (flg) {
434:       PCSPAISetMax(pc,max1);
435:     }
436:     PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);
437:     if (flg) {
438:       PCSPAISetMaxNew(pc,maxnew1);
439:     }
440:     PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);
441:     if (flg) {
442:       PCSPAISetBlockSize(pc,block_size1);
443:     }
444:     PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);
445:     if (flg) {
446:       PCSPAISetCacheSize(pc,cache_size);
447:     }
448:     PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);
449:     if (flg) {
450:       PCSPAISetVerbose(pc,verbose);
451:     }
452:     PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);
453:     if (flg) {
454:       PCSPAISetSp(pc,sp);
455:     }
456:   PetscOptionsTail();
457:   return(0);
458: }

460: /**********************************************************************/

462: EXTERN_C_BEGIN
463: /*
464:    PCCreate_SPAI - Creates the preconditioner context for the SPAI 
465:                    preconditioner written by Stephen Barnard.

467: */
468: #undef __FUNCT__  
470: int PCCreate_SPAI(PC pc)
471: {
472:   PC_SPAI *ispai;
473:   int     ierr;

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

479:   pc->ops->destroy         = PCDestroy_SPAI;
480:   pc->ops->apply           = PCApply_SPAI;
481:   pc->ops->applyrichardson = 0;
482:   pc->ops->setup           = PCSetUp_SPAI;
483:   pc->ops->view            = PCView_SPAI;
484:   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;

486:   pc->name          = 0;
487:   ispai->epsilon    = .4;
488:   ispai->nbsteps    = 5;
489:   ispai->max        = 5000;
490:   ispai->maxnew     = 5;
491:   ispai->block_size = 1;
492:   ispai->cache_size = 5;
493:   ispai->verbose    = 0;

495:   ispai->sp         = 1;

497:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C",
498:                     "PCSPAISetEpsilon_SPAI",
499:                      PCSPAISetEpsilon_SPAI);
500:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C",
501:                     "PCSPAISetNBSteps_SPAI",
502:                      PCSPAISetNBSteps_SPAI);
503:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C",
504:                     "PCSPAISetMax_SPAI",
505:                      PCSPAISetMax_SPAI);
506:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC",
507:                     "PCSPAISetMaxNew_SPAI",
508:                      PCSPAISetMaxNew_SPAI);
509:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C",
510:                     "PCSPAISetBlockSize_SPAI",
511:                      PCSPAISetBlockSize_SPAI);
512:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C",
513:                     "PCSPAISetCacheSize_SPAI",
514:                      PCSPAISetCacheSize_SPAI);
515:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C",
516:                     "PCSPAISetVerbose_SPAI",
517:                      PCSPAISetVerbose_SPAI);
518:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C",
519:                     "PCSPAISetSp_SPAI",
520:                      PCSPAISetSp_SPAI);

522:   return(0);
523: }
524: EXTERN_C_END

526: /**********************************************************************/

528: /*
529:    Converts from a PETSc matrix to an SPAI matrix 
530: */
531: #undef __FUNCT__  
533: int ConvertMatToMatrix(Mat A,Mat AT,matrix **B)
534: {
535:   matrix   *M;
536:   int      i,j,col;
537:   int      row_indx;
538:   int      len,pe,local_indx,start_indx;
539:   int      *mapping;
540:   int      ierr,*cols;
541:   double   *vals;
542:   int      *num_ptr,n,mnl,nnl,rank,size,nz,rstart,rend;
543:   MPI_Comm comm;
544:   struct compressed_lines *rows;

547:   PetscObjectGetComm((PetscObject)A,&comm);
548: 
549:   MPI_Comm_size(comm,&size);
550:   MPI_Comm_rank(comm,&rank);
551:   MatGetSize(A,&n,&n);
552:   MatGetLocalSize(A,&mnl,&nnl);

554:   MPI_Barrier(PETSC_COMM_WORLD);

556:   M = new_matrix((void *)comm);
557: 
558:   M->n = n;
559:   M->bs = 1;
560:   M->max_block_size = 1;

562:   M->mnls          = (int*)malloc(sizeof(int)*size);
563:   M->start_indices = (int*)malloc(sizeof(int)*size);
564:   M->pe            = (int*)malloc(sizeof(int)*n);
565:   M->block_sizes   = (int*)malloc(sizeof(int)*n);
566:   for (i=0; i<n; i++) M->block_sizes[i] = 1;

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

570:   M->start_indices[0] = 0;
571:   for (i=1; i<size; i++) {
572:     M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
573:   }

575:   M->mnl = M->mnls[M->myid];
576:   M->my_start_index = M->start_indices[M->myid];

578:   for (i=0; i<size; i++) {
579:     start_indx = M->start_indices[i];
580:     for (j=0; j<M->mnls[i]; j++)
581:       M->pe[start_indx+j] = i;
582:   }

584:   if (AT) {
585:     M->lines = new_compressed_lines(M->mnls[rank],1);
586:   } else {
587:     M->lines = new_compressed_lines(M->mnls[rank],0);
588:   }

590:   rows     = M->lines;

592:   /* Determine the mapping from global indices to pointers */
593:   ierr       = PetscMalloc(M->n*sizeof(int),&mapping);
594:   pe         = 0;
595:   local_indx = 0;
596:   for (i=0; i<M->n; i++) {
597:     if (local_indx >= M->mnls[pe]) {
598:       pe++;
599:       local_indx = 0;
600:     }
601:     mapping[i] = local_indx + M->start_indices[pe];
602:     local_indx++;
603:   }


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

608:   /*********************************************************/
609:   /************** Set up the row structure *****************/
610:   /*********************************************************/

612:   /* count number of nonzeros in every row */
613:   MatGetOwnershipRange(A,&rstart,&rend);
614:   for (i=rstart; i<rend; i++) {
615:     MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
616:     MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
617:   }

619:   /* allocate buffers */
620:   len = 0;
621:   for (i=0; i<mnl; i++) {
622:     if (len < num_ptr[i]) len = num_ptr[i];
623:   }

625:   for (i=rstart; i<rend; i++) {
626:     row_indx             = i-rstart;
627:     len                  = num_ptr[row_indx];
628:     rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int));
629:     rows->A[row_indx]    = (double*)malloc(len*sizeof(double));
630:   }

632:   /* copy the matrix */
633:   for (i=rstart; i<rend; i++) {
634:     row_indx = i - rstart;
635:     ierr     = MatGetRow(A,i,&nz,&cols,&vals);
636:     for (j=0; j<nz; j++) {
637:       col = cols[j];
638:       len = rows->len[row_indx]++;
639:       rows->ptrs[row_indx][len] = mapping[col];
640:       rows->A[row_indx][len]    = vals[j];
641:     }
642:     rows->slen[row_indx] = rows->len[row_indx];
643:     MatRestoreRow(A,i,&nz,&cols,&vals);
644:   }


647:   /************************************************************/
648:   /************** Set up the column structure *****************/
649:   /*********************************************************/

651:   if (AT) {

653:     /* count number of nonzeros in every column */
654:     for (i=rstart; i<rend; i++) {
655:       MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
656:       MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);
657:     }

659:     /* allocate buffers */
660:     len = 0;
661:     for (i=0; i<mnl; i++) {
662:       if (len < num_ptr[i]) len = num_ptr[i];
663:     }

665:     for (i=rstart; i<rend; i++) {
666:       row_indx = i-rstart;
667:       len      = num_ptr[row_indx];
668:       rows->rptrs[row_indx]     = (int*)malloc(len*sizeof(int));
669:     }

671:     /* copy the matrix (i.e., the structure) */
672:     for (i=rstart; i<rend; i++) {
673:       row_indx = i - rstart;
674:       ierr     = MatGetRow(AT,i,&nz,&cols,&vals);
675:       for (j=0; j<nz; j++) {
676:         col = cols[j];
677:         len = rows->rlen[row_indx]++;
678:         rows->rptrs[row_indx][len] = mapping[col];
679:       }
680:       ierr     = MatRestoreRow(AT,i,&nz,&cols,&vals);
681:     }
682:   }

684:   PetscFree(num_ptr);
685:   PetscFree(mapping);

687:   order_pointers(M);
688:   M->maxnz = calc_maxnz(M);

690:   *B = M;

692:   return(0);
693: }

695: /**********************************************************************/

697: /*
698:    Converts from an SPAI matrix B  to a PETSc matrix PB.
699:    This assumes that the the SPAI matrix B is stored in
700:    COMPRESSED-ROW format.
701: */
702: #undef __FUNCT__  
704: int ConvertMatrixToMat(matrix *B,Mat *PB)
705: {
706:   int         size,rank;
707:   int         ierr;
708:   int         m,n,M,N;
709:   int         d_nz,o_nz;
710:   int         *d_nnz,*o_nnz;
711:   int         i,k,global_row,global_col,first_diag_col,last_diag_col;
712:   PetscScalar val;

715:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
716:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
717: 
718:   m = n = B->mnls[rank];
719:   d_nz = o_nz = 0;

721:   /* Determine preallocation for MatCreateMPIAIJ */
722:   PetscMalloc(m*sizeof(int),&d_nnz);
723:   PetscMalloc(m*sizeof(int),&o_nnz);
724:   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
725:   first_diag_col = B->start_indices[rank];
726:   last_diag_col = first_diag_col + B->mnls[rank];
727:   for (i=0; i<B->mnls[rank]; i++) {
728:     for (k=0; k<B->lines->len[i]; k++) {
729:       global_col = B->lines->ptrs[i][k];
730:       if ((global_col >= first_diag_col) && (global_col <= last_diag_col))
731:         d_nnz[i]++;
732:       else
733:         o_nnz[i]++;
734:     }
735:   }

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

740:   for (i=0; i<B->mnls[rank]; i++) {
741:     global_row = B->start_indices[rank]+i;
742:     for (k=0; k<B->lines->len[i]; k++) {
743:       global_col = B->lines->ptrs[i][k];
744:       val = B->lines->A[i][k];
745:       MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);
746:     }
747:   }

749:   PetscFree(d_nnz);
750:   PetscFree(o_nnz);

752:   MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);
753:   MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);

755:   return(0);
756: }

758: /**********************************************************************/

760: /*
761:    Converts from an SPAI vector v  to a PETSc vec Pv.
762: */
763: #undef __FUNCT__  
765: int ConvertVectorToVec(vector *v,Vec *Pv)
766: {
767:   int size,rank,ierr,m,M,i,*mnls,*start_indices,*global_indices;
768: 
770:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
771:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
772: 
773:   m = v->mnl;
774:   M = v->n;
775: 
776: 
777:   VecCreateMPI(PETSC_COMM_WORLD,m,M,Pv);

779:   PetscMalloc(size*sizeof(int),&mnls);
780:   MPI_Allgather((void*)&v->mnl,1,MPI_INT,(void*)mnls,1,MPI_INT,PETSC_COMM_WORLD);
781: 
782:   PetscMalloc(size*sizeof(int),&start_indices);
783:   start_indices[0] = 0;
784:   for (i=1; i<size; i++)
785:     start_indices[i] = start_indices[i-1] +mnls[i-1];
786: 
787:   PetscMalloc(v->mnl*sizeof(int),&global_indices);
788:   for (i=0; i<v->mnl; i++)
789:     global_indices[i] = start_indices[rank] + i;

791:   PetscFree(mnls);
792:   PetscFree(start_indices);
793: 
794:   VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);

796:   PetscFree(global_indices);

798:   VecAssemblyBegin(*Pv);
799:   VecAssemblyEnd(*Pv);
800: 
801:   return(0);
802: }

804: /**********************************************************************/

806: /*

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

811:   !!!! Warning!!!!! PETSc supports only serial execution of this routine.
812:   
813:   f0 <input_file> : matrix in Matrix Market format
814:   f1 <input_file> : vector in Matrix Market array format
815:   f2 <output_file> : matrix and vector in PETSc format

817: */
818: #undef __FUNCT__  
820: int MM_to_PETSC(char *f0,char *f1,char *f2)
821: {
822:   Mat         A_PETSC;          /* matrix */
823:   Vec         b_PETSC;          /* RHS */
824:   PetscViewer fd;               /* viewer */
825:   int         ierr;
826:   matrix      *A_spai;
827:   vector      *b_spai;

830:   A_spai = read_mm_matrix(f0,1,1,1,0,0,0,NULL);
831:   b_spai = read_rhs_for_matrix(f1,A_spai);

833:   ConvertMatrixToMat(A_spai,&A_PETSC);
834:   ConvertVectorToVec(b_spai,&b_PETSC);

836:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,f1,PETSC_BINARY_CREATE,&fd);
837:   MatView(A_PETSC,fd);
838:   VecView(b_PETSC,fd);

840:   return(0);
841: }