Actual source code: ispai.c

petsc-master 2017-01-18
Report Typos and Errors

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

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

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

 12:       This code needs work

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

 19: */
 20: #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */

 22:  #include <petsc/private/pcimpl.h>
 23:  #include <../src/ksp/pc/impls/spai/petscspai.h>

 25: /*
 26:     These are the SPAI include files
 27: */
 28: EXTERN_C_BEGIN
 29: #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */
 30: #include <spai.h>
 31: #include <matrix.h>
 32: EXTERN_C_END

 34: extern PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**);
 35: extern PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix*,Mat*);
 36: extern PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv);
 37: extern PetscErrorCode MM_to_PETSC(char*,char*,char*);

 39: typedef struct {

 41:   matrix *B;                /* matrix in SPAI format */
 42:   matrix *BT;               /* transpose of matrix in SPAI format */
 43:   matrix *M;                /* the approximate inverse in SPAI format */

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

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

 55:   int      sp;              /* symmetric nonzero pattern */
 56:   MPI_Comm comm_spai;     /* communicator to be used with spai */
 57: } PC_SPAI;

 59: /**********************************************************************/

 61: static PetscErrorCode PCSetUp_SPAI(PC pc)
 62: {
 63:   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
 65:   Mat            AT;

 68:   init_SPAI();

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

 79:   /* Destroy the transpose */
 80:   /* Don't know how to do it. PETSc developers? */

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

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

105:   ConvertMatrixToMat(PetscObjectComm((PetscObject)pc),ispai->M,&ispai->PM);

107:   /* free the SPAI matrices */
108:   sp_free_matrix(ispai->B);
109:   sp_free_matrix(ispai->M);
110:   return(0);
111: }

113: /**********************************************************************/

115: static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y)
116: {
117:   PC_SPAI        *ispai = (PC_SPAI*)pc->data;

121:   /* Now using PETSc's multiply */
122:   MatMult(ispai->PM,xx,y);
123:   return(0);
124: }

126: /**********************************************************************/

128: static PetscErrorCode PCDestroy_SPAI(PC pc)
129: {
131:   PC_SPAI        *ispai = (PC_SPAI*)pc->data;

134:   MatDestroy(&ispai->PM);
135:   MPI_Comm_free(&(ispai->comm_spai));
136:   PetscFree(pc->data);
137:   return(0);
138: }

140: /**********************************************************************/

142: static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer)
143: {
144:   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
146:   PetscBool      iascii;

149:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
150:   if (iascii) {
151:     PetscViewerASCIIPrintf(viewer,"    SPAI preconditioner\n");
152:     PetscViewerASCIIPrintf(viewer,"    epsilon %g\n",   (double)ispai->epsilon);
153:     PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps);
154:     PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max);
155:     PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew);
156:     PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size);
157:     PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size);
158:     PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose);
159:     PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp);
160:   }
161:   return(0);
162: }

164: static PetscErrorCode  PCSPAISetEpsilon_SPAI(PC pc,double epsilon1)
165: {
166:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

169:   ispai->epsilon = epsilon1;
170:   return(0);
171: }

173: /**********************************************************************/

175: static PetscErrorCode  PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1)
176: {
177:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

180:   ispai->nbsteps = nbsteps1;
181:   return(0);
182: }

184: /**********************************************************************/

186: /* added 1/7/99 g.h. */
187: static PetscErrorCode  PCSPAISetMax_SPAI(PC pc,int max1)
188: {
189:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

192:   ispai->max = max1;
193:   return(0);
194: }

196: /**********************************************************************/

198: static PetscErrorCode  PCSPAISetMaxNew_SPAI(PC pc,int maxnew1)
199: {
200:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

203:   ispai->maxnew = maxnew1;
204:   return(0);
205: }

207: /**********************************************************************/

209: static PetscErrorCode  PCSPAISetBlockSize_SPAI(PC pc,int block_size1)
210: {
211:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

214:   ispai->block_size = block_size1;
215:   return(0);
216: }

218: /**********************************************************************/

220: static PetscErrorCode  PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
221: {
222:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

225:   ispai->cache_size = cache_size;
226:   return(0);
227: }

229: /**********************************************************************/

231: static PetscErrorCode  PCSPAISetVerbose_SPAI(PC pc,int verbose)
232: {
233:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

236:   ispai->verbose = verbose;
237:   return(0);
238: }

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

242: static PetscErrorCode  PCSPAISetSp_SPAI(PC pc,int sp)
243: {
244:   PC_SPAI *ispai = (PC_SPAI*)pc->data;

247:   ispai->sp = sp;
248:   return(0);
249: }

251: /* -------------------------------------------------------------------*/

253: /*@
254:   PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner

256:   Input Parameters:
257: + pc - the preconditioner
258: - eps - epsilon (default .4)

260:   Notes:  Espilon must be between 0 and 1. It controls the
261:                  quality of the approximation of M to the inverse of
262:                  A. Higher values of epsilon lead to more work, more
263:                  fill, and usually better preconditioners. In many
264:                  cases the best choice of epsilon is the one that
265:                  divides the total solution time equally between the
266:                  preconditioner and the solver.

268:   Level: intermediate

270: .seealso: PCSPAI, PCSetType()
271:   @*/
272: PetscErrorCode  PCSPAISetEpsilon(PC pc,double epsilon1)
273: {

277:   PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));
278:   return(0);
279: }

281: /**********************************************************************/

283: /*@
284:   PCSPAISetNBSteps - set maximum number of improvement steps per row in
285:         the SPAI preconditioner

287:   Input Parameters:
288: + pc - the preconditioner
289: - n - number of steps (default 5)

291:   Notes:  SPAI constructs to approximation to every column of
292:                  the exact inverse of A in a series of improvement
293:                  steps. The quality of the approximation is determined
294:                  by epsilon. If an approximation achieving an accuracy
295:                  of epsilon is not obtained after ns steps, SPAI simply
296:                  uses the best approximation constructed so far.

298:   Level: intermediate

300: .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew()
301: @*/
302: PetscErrorCode  PCSPAISetNBSteps(PC pc,int nbsteps1)
303: {

307:   PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));
308:   return(0);
309: }

311: /**********************************************************************/

313: /* added 1/7/99 g.h. */
314: /*@
315:   PCSPAISetMax - set the size of various working buffers in
316:         the SPAI preconditioner

318:   Input Parameters:
319: + pc - the preconditioner
320: - n - size (default is 5000)

322:   Level: intermediate

324: .seealso: PCSPAI, PCSetType()
325: @*/
326: PetscErrorCode  PCSPAISetMax(PC pc,int max1)
327: {

331:   PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));
332:   return(0);
333: }

335: /**********************************************************************/

337: /*@
338:   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
339:    in SPAI preconditioner

341:   Input Parameters:
342: + pc - the preconditioner
343: - n - maximum number (default 5)

345:   Level: intermediate

347: .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps()
348: @*/
349: PetscErrorCode  PCSPAISetMaxNew(PC pc,int maxnew1)
350: {

354:   PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));
355:   return(0);
356: }

358: /**********************************************************************/

360: /*@
361:   PCSPAISetBlockSize - set the block size for the SPAI preconditioner

363:   Input Parameters:
364: + pc - the preconditioner
365: - n - block size (default 1)

367:   Notes: A block
368:                  size of 1 treats A as a matrix of scalar elements. A
369:                  block size of s > 1 treats A as a matrix of sxs
370:                  blocks. A block size of 0 treats A as a matrix with
371:                  variable sized blocks, which are determined by
372:                  searching for dense square diagonal blocks in A.
373:                  This can be very effective for finite-element
374:                  matrices.

376:                  SPAI will convert A to block form, use a block
377:                  version of the preconditioner algorithm, and then
378:                  convert the result back to scalar form.

380:                  In many cases the a block-size parameter other than 1
381:                  can lead to very significant improvement in
382:                  performance.


385:   Level: intermediate

387: .seealso: PCSPAI, PCSetType()
388: @*/
389: PetscErrorCode  PCSPAISetBlockSize(PC pc,int block_size1)
390: {

394:   PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));
395:   return(0);
396: }

398: /**********************************************************************/

400: /*@
401:   PCSPAISetCacheSize - specify cache size in the SPAI preconditioner

403:   Input Parameters:
404: + pc - the preconditioner
405: - n -  cache size {0,1,2,3,4,5} (default 5)

407:   Notes:    SPAI uses a hash table to cache messages and avoid
408:                  redundant communication. If suggest always using
409:                  5. This parameter is irrelevant in the serial
410:                  version.

412:   Level: intermediate

414: .seealso: PCSPAI, PCSetType()
415: @*/
416: PetscErrorCode  PCSPAISetCacheSize(PC pc,int cache_size)
417: {

421:   PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));
422:   return(0);
423: }

425: /**********************************************************************/

427: /*@
428:   PCSPAISetVerbose - verbosity level for the SPAI preconditioner

430:   Input Parameters:
431: + pc - the preconditioner
432: - n - level (default 1)

434:   Notes: print parameters, timings and matrix statistics

436:   Level: intermediate

438: .seealso: PCSPAI, PCSetType()
439: @*/
440: PetscErrorCode  PCSPAISetVerbose(PC pc,int verbose)
441: {

445:   PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));
446:   return(0);
447: }

449: /**********************************************************************/

451: /*@
452:   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner

454:   Input Parameters:
455: + pc - the preconditioner
456: - n - 0 or 1

458:   Notes: If A has a symmetric nonzero pattern use -sp 1 to
459:                  improve performance by eliminating some communication
460:                  in the parallel version. Even if A does not have a
461:                  symmetric nonzero pattern -sp 1 may well lead to good
462:                  results, but the code will not follow the published
463:                  SPAI algorithm exactly.


466:   Level: intermediate

468: .seealso: PCSPAI, PCSetType()
469: @*/
470: PetscErrorCode  PCSPAISetSp(PC pc,int sp)
471: {

475:   PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));
476:   return(0);
477: }

479: /**********************************************************************/

481: /**********************************************************************/

483: static PetscErrorCode PCSetFromOptions_SPAI(PetscOptionItems *PetscOptionsObject,PC pc)
484: {
485:   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
487:   int            nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
488:   double         epsilon1;
489:   PetscBool      flg;

492:   PetscOptionsHead(PetscOptionsObject,"SPAI options");
493:   PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);
494:   if (flg) {
495:     PCSPAISetEpsilon(pc,epsilon1);
496:   }
497:   PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);
498:   if (flg) {
499:     PCSPAISetNBSteps(pc,nbsteps1);
500:   }
501:   /* added 1/7/99 g.h. */
502:   PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);
503:   if (flg) {
504:     PCSPAISetMax(pc,max1);
505:   }
506:   PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);
507:   if (flg) {
508:     PCSPAISetMaxNew(pc,maxnew1);
509:   }
510:   PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);
511:   if (flg) {
512:     PCSPAISetBlockSize(pc,block_size1);
513:   }
514:   PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);
515:   if (flg) {
516:     PCSPAISetCacheSize(pc,cache_size);
517:   }
518:   PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);
519:   if (flg) {
520:     PCSPAISetVerbose(pc,verbose);
521:   }
522:   PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);
523:   if (flg) {
524:     PCSPAISetSp(pc,sp);
525:   }
526:   PetscOptionsTail();
527:   return(0);
528: }

530: /**********************************************************************/

532: /*MC
533:    PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard
534:      as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)

536:    Options Database Keys:
537: +  -pc_spai_epsilon <eps> - set tolerance
538: .  -pc_spai_nbstep <n> - set nbsteps
539: .  -pc_spai_max <m> - set max
540: .  -pc_spai_max_new <m> - set maxnew
541: .  -pc_spai_block_size <n> - set block size
542: .  -pc_spai_cache_size <n> - set cache size
543: .  -pc_spai_sp <m> - set sp
544: -  -pc_spai_set_verbose <true,false> - verbose output

546:    Notes: This only works with AIJ matrices.

548:    Level: beginner

550:    Concepts: approximate inverse

552: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
553:     PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(),
554:     PCSPAISetVerbose(), PCSPAISetSp()
555: M*/

557: PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc)
558: {
559:   PC_SPAI        *ispai;

563:   PetscNewLog(pc,&ispai);
564:   pc->data = ispai;

566:   pc->ops->destroy         = PCDestroy_SPAI;
567:   pc->ops->apply           = PCApply_SPAI;
568:   pc->ops->applyrichardson = 0;
569:   pc->ops->setup           = PCSetUp_SPAI;
570:   pc->ops->view            = PCView_SPAI;
571:   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;

573:   ispai->epsilon    = .4;
574:   ispai->nbsteps    = 5;
575:   ispai->max        = 5000;
576:   ispai->maxnew     = 5;
577:   ispai->block_size = 1;
578:   ispai->cache_size = 5;
579:   ispai->verbose    = 0;

581:   ispai->sp = 1;
582:   MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai));

584:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",PCSPAISetEpsilon_SPAI);
585:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",PCSPAISetNBSteps_SPAI);
586:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",PCSPAISetMax_SPAI);
587:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",PCSPAISetMaxNew_SPAI);
588:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",PCSPAISetBlockSize_SPAI);
589:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",PCSPAISetCacheSize_SPAI);
590:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",PCSPAISetVerbose_SPAI);
591:   PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",PCSPAISetSp_SPAI);
592:   return(0);
593: }

595: /**********************************************************************/

597: /*
598:    Converts from a PETSc matrix to an SPAI matrix
599: */
600: PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B)
601: {
602:   matrix                  *M;
603:   int                     i,j,col;
604:   int                     row_indx;
605:   int                     len,pe,local_indx,start_indx;
606:   int                     *mapping;
607:   PetscErrorCode          ierr;
608:   const int               *cols;
609:   const double            *vals;
610:   int                     n,mnl,nnl,nz,rstart,rend;
611:   PetscMPIInt             size,rank;
612:   struct compressed_lines *rows;

615:   MPI_Comm_size(comm,&size);
616:   MPI_Comm_rank(comm,&rank);
617:   MatGetSize(A,&n,&n);
618:   MatGetLocalSize(A,&mnl,&nnl);

620:   /*
621:     not sure why a barrier is required. commenting out
622:   MPI_Barrier(comm);
623:   */

625:   M = new_matrix((SPAI_Comm)comm);

627:   M->n              = n;
628:   M->bs             = 1;
629:   M->max_block_size = 1;

631:   M->mnls          = (int*)malloc(sizeof(int)*size);
632:   M->start_indices = (int*)malloc(sizeof(int)*size);
633:   M->pe            = (int*)malloc(sizeof(int)*n);
634:   M->block_sizes   = (int*)malloc(sizeof(int)*n);
635:   for (i=0; i<n; i++) M->block_sizes[i] = 1;

637:   MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);

639:   M->start_indices[0] = 0;
640:   for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];

642:   M->mnl            = M->mnls[M->myid];
643:   M->my_start_index = M->start_indices[M->myid];

645:   for (i=0; i<size; i++) {
646:     start_indx = M->start_indices[i];
647:     for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i;
648:   }

650:   if (AT) {
651:     M->lines = new_compressed_lines(M->mnls[rank],1);
652:   } else {
653:     M->lines = new_compressed_lines(M->mnls[rank],0);
654:   }

656:   rows = M->lines;

658:   /* Determine the mapping from global indices to pointers */
659:   PetscMalloc1(M->n,&mapping);
660:   pe         = 0;
661:   local_indx = 0;
662:   for (i=0; i<M->n; i++) {
663:     if (local_indx >= M->mnls[pe]) {
664:       pe++;
665:       local_indx = 0;
666:     }
667:     mapping[i] = local_indx + M->start_indices[pe];
668:     local_indx++;
669:   }

671:   /*********************************************************/
672:   /************** Set up the row structure *****************/
673:   /*********************************************************/

675:   MatGetOwnershipRange(A,&rstart,&rend);
676:   for (i=rstart; i<rend; i++) {
677:     row_indx = i - rstart;
678:     MatGetRow(A,i,&nz,&cols,&vals);
679:     /* allocate buffers */
680:     rows->ptrs[row_indx] = (int*)malloc(nz*sizeof(int));
681:     rows->A[row_indx]    = (double*)malloc(nz*sizeof(double));
682:     /* copy the matrix */
683:     for (j=0; j<nz; j++) {
684:       col = cols[j];
685:       len = rows->len[row_indx]++;

687:       rows->ptrs[row_indx][len] = mapping[col];
688:       rows->A[row_indx][len]    = vals[j];
689:     }
690:     rows->slen[row_indx] = rows->len[row_indx];

692:     MatRestoreRow(A,i,&nz,&cols,&vals);
693:   }


696:   /************************************************************/
697:   /************** Set up the column structure *****************/
698:   /*********************************************************/

700:   if (AT) {

702:     for (i=rstart; i<rend; i++) {
703:       row_indx = i - rstart;
704:       MatGetRow(AT,i,&nz,&cols,&vals);
705:       /* allocate buffers */
706:       rows->rptrs[row_indx] = (int*)malloc(nz*sizeof(int));
707:       /* copy the matrix (i.e., the structure) */
708:       for (j=0; j<nz; j++) {
709:         col = cols[j];
710:         len = rows->rlen[row_indx]++;

712:         rows->rptrs[row_indx][len] = mapping[col];
713:       }
714:       MatRestoreRow(AT,i,&nz,&cols,&vals);
715:     }
716:   }

718:   PetscFree(mapping);

720:   order_pointers(M);
721:   M->maxnz = calc_maxnz(M);
722:   *B       = M;
723:   return(0);
724: }

726: /**********************************************************************/

728: /*
729:    Converts from an SPAI matrix B  to a PETSc matrix PB.
730:    This assumes that the SPAI matrix B is stored in
731:    COMPRESSED-ROW format.
732: */
733: PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB)
734: {
735:   PetscMPIInt    size,rank;
737:   int            m,n,M,N;
738:   int            d_nz,o_nz;
739:   int            *d_nnz,*o_nnz;
740:   int            i,k,global_row,global_col,first_diag_col,last_diag_col;
741:   PetscScalar    val;

744:   MPI_Comm_size(comm,&size);
745:   MPI_Comm_rank(comm,&rank);

747:   m    = n = B->mnls[rank];
748:   d_nz = o_nz = 0;

750:   /* Determine preallocation for MatCreateMPIAIJ */
751:   PetscMalloc1(m,&d_nnz);
752:   PetscMalloc1(m,&o_nnz);
753:   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
754:   first_diag_col = B->start_indices[rank];
755:   last_diag_col  = first_diag_col + B->mnls[rank];
756:   for (i=0; i<B->mnls[rank]; i++) {
757:     for (k=0; k<B->lines->len[i]; k++) {
758:       global_col = B->lines->ptrs[i][k];
759:       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
760:       else o_nnz[i]++;
761:     }
762:   }

764:   M = N = B->n;
765:   /* Here we only know how to create AIJ format */
766:   MatCreate(comm,PB);
767:   MatSetSizes(*PB,m,n,M,N);
768:   MatSetType(*PB,MATAIJ);
769:   MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);
770:   MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);

772:   for (i=0; i<B->mnls[rank]; i++) {
773:     global_row = B->start_indices[rank]+i;
774:     for (k=0; k<B->lines->len[i]; k++) {
775:       global_col = B->lines->ptrs[i][k];

777:       val  = B->lines->A[i][k];
778:       MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);
779:     }
780:   }

782:   PetscFree(d_nnz);
783:   PetscFree(o_nnz);

785:   MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);
786:   MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);
787:   return(0);
788: }

790: /**********************************************************************/

792: /*
793:    Converts from an SPAI vector v  to a PETSc vec Pv.
794: */
795: PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv)
796: {
798:   PetscMPIInt    size,rank;
799:   int            m,M,i,*mnls,*start_indices,*global_indices;

802:   MPI_Comm_size(comm,&size);
803:   MPI_Comm_rank(comm,&rank);

805:   m = v->mnl;
806:   M = v->n;


809:   VecCreateMPI(comm,m,M,Pv);

811:   PetscMalloc1(size,&mnls);
812:   MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);

814:   PetscMalloc1(size,&start_indices);

816:   start_indices[0] = 0;
817:   for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1];

819:   PetscMalloc1(v->mnl,&global_indices);
820:   for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i;

822:   PetscFree(mnls);
823:   PetscFree(start_indices);

825:   VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);
826:   VecAssemblyBegin(*Pv);
827:   VecAssemblyEnd(*Pv);

829:   PetscFree(global_indices);
830:   return(0);
831: }