Actual source code: ispai.c

petsc-3.9.4 2018-09-11
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,"    epsilon %g\n",   (double)ispai->epsilon);
152:     PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps);
153:     PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max);
154:     PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew);
155:     PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size);
156:     PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size);
157:     PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose);
158:     PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp);
159:   }
160:   return(0);
161: }

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

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

172: /**********************************************************************/

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

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

183: /**********************************************************************/

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

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

195: /**********************************************************************/

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

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

206: /**********************************************************************/

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

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

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

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

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

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

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

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

239: /**********************************************************************/

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

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

250: /* -------------------------------------------------------------------*/

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

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

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

267:   Level: intermediate

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

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

280: /**********************************************************************/

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

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

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

297:   Level: intermediate

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

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

310: /**********************************************************************/

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

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

321:   Level: intermediate

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

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

334: /**********************************************************************/

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

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

344:   Level: intermediate

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

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

357: /**********************************************************************/

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

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

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

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

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


384:   Level: intermediate

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

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

397: /**********************************************************************/

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

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

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

411:   Level: intermediate

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

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

424: /**********************************************************************/

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

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

433:   Notes: print parameters, timings and matrix statistics

435:   Level: intermediate

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

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

448: /**********************************************************************/

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

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

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


465:   Level: intermediate

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

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

478: /**********************************************************************/

480: /**********************************************************************/

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

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

529: /**********************************************************************/

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

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

545:    Notes: This only works with AIJ matrices.

547:    Level: beginner

549:    Concepts: approximate inverse

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

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

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

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

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

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

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

594: /**********************************************************************/

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

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

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

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

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

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

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

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

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

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

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

655:   rows = M->lines;

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

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

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

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

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


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

699:   if (AT) {

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

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

717:   PetscFree(mapping);

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

725: /**********************************************************************/

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

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

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

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

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

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

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

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

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

789: /**********************************************************************/

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

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

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


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

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

813:   PetscMalloc1(size,&start_indices);

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

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

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

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

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