xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
101a81e61SBarry Smith 
201a81e61SBarry Smith /*
301a81e61SBarry Smith    3/99 Modified by Stephen Barnard to support SPAI version 3.0
401a81e61SBarry Smith */
501a81e61SBarry Smith 
601a81e61SBarry Smith /*
701a81e61SBarry Smith       Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner
801a81e61SBarry Smith    Code written by Stephen Barnard.
901a81e61SBarry Smith 
1001a81e61SBarry Smith       Note: there is some BAD memory bleeding below!
1101a81e61SBarry Smith 
1201a81e61SBarry Smith       This code needs work
1301a81e61SBarry Smith 
1401a81e61SBarry Smith    1) get rid of all memory bleeding
1501a81e61SBarry Smith    2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix
1601a81e61SBarry Smith       rather than having the sp flag for PC_SPAI
172a4967abSBarry Smith    3) fix to set the block size based on the matrix block size
1801a81e61SBarry Smith 
1901a81e61SBarry Smith */
20762360ffSBarry Smith #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */
2101a81e61SBarry Smith 
22af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
231c7d2463SJed Brown #include <../src/ksp/pc/impls/spai/petscspai.h>
2401a81e61SBarry Smith 
2501a81e61SBarry Smith /*
2601a81e61SBarry Smith     These are the SPAI include files
2701a81e61SBarry Smith */
2801a81e61SBarry Smith EXTERN_C_BEGIN
29b1bda85cSSatish Balay #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */
30c6db04a5SJed Brown #include <spai.h>
31c6db04a5SJed Brown #include <matrix.h>
3201a81e61SBarry Smith EXTERN_C_END
3301a81e61SBarry Smith 
3409573ac7SBarry Smith extern PetscErrorCode ConvertMatToMatrix(MPI_Comm, Mat, Mat, matrix **);
3509573ac7SBarry Smith extern PetscErrorCode ConvertMatrixToMat(MPI_Comm, matrix *, Mat *);
36b03515a0SUmberto Zerbinati extern PetscErrorCode ConvertVectorToVec(MPI_Comm, vector *, Vec *);
3709573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char *, char *, char *);
3801a81e61SBarry Smith 
3901a81e61SBarry Smith typedef struct {
4001a81e61SBarry Smith   matrix *B;  /* matrix in SPAI format */
4101a81e61SBarry Smith   matrix *BT; /* transpose of matrix in SPAI format */
4201a81e61SBarry Smith   matrix *M;  /* the approximate inverse in SPAI format */
4301a81e61SBarry Smith 
4401a81e61SBarry Smith   Mat PM; /* the approximate inverse PETSc format */
4501a81e61SBarry Smith 
4601a81e61SBarry Smith   double epsilon;    /* tolerance */
4701a81e61SBarry Smith   int    nbsteps;    /* max number of "improvement" steps per line */
4801a81e61SBarry Smith   int    max;        /* max dimensions of is_I, q, etc. */
4901a81e61SBarry Smith   int    maxnew;     /* max number of new entries per step */
5001a81e61SBarry Smith   int    block_size; /* constant block size */
5101a81e61SBarry Smith   int    cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */
5201a81e61SBarry Smith   int    verbose;    /* SPAI prints timing and statistics */
5301a81e61SBarry Smith 
5401a81e61SBarry Smith   int      sp;        /* symmetric nonzero pattern */
5501a81e61SBarry Smith   MPI_Comm comm_spai; /* communicator to be used with spai */
5601a81e61SBarry Smith } PC_SPAI;
5701a81e61SBarry Smith 
589371c9d4SSatish Balay static PetscErrorCode PCSetUp_SPAI(PC pc) {
5901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
6001a81e61SBarry Smith   Mat      AT;
6101a81e61SBarry Smith 
6201a81e61SBarry Smith   PetscFunctionBegin;
6301a81e61SBarry Smith   init_SPAI();
6401a81e61SBarry Smith 
6501a81e61SBarry Smith   if (ispai->sp) {
669566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, pc->pmat, &ispai->B));
6701a81e61SBarry Smith   } else {
6801a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
699566063dSJacob Faibussowitsch     PetscCall(MatTranspose(pc->pmat, MAT_INITIAL_MATRIX, &AT));
709566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, AT, &ispai->B));
719566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
7201a81e61SBarry Smith   }
7301a81e61SBarry Smith 
7401a81e61SBarry Smith   /* Destroy the transpose */
7501a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
7601a81e61SBarry Smith 
7701a81e61SBarry Smith   /* construct SPAI preconditioner */
7801a81e61SBarry Smith   /* FILE *messages */    /* file for warning messages */
7901a81e61SBarry Smith   /* double epsilon */    /* tolerance */
8001a81e61SBarry Smith   /* int nbsteps */       /* max number of "improvement" steps per line */
8101a81e61SBarry Smith   /* int max */           /* max dimensions of is_I, q, etc. */
8201a81e61SBarry Smith   /* int maxnew */        /* max number of new entries per step */
8301a81e61SBarry Smith   /* int block_size */    /* block_size == 1 specifies scalar elments
8401a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
8501a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
862fa5cd67SKarl Rupp   /* int cache_size */    /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */
8701a81e61SBarry Smith   /* int    verbose    */ /* verbose == 0 specifies that SPAI is silent
8801a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
8901a81e61SBarry Smith 
909371c9d4SSatish Balay   PetscCall(bspai(ispai->B, &ispai->M, stdout, ispai->epsilon, ispai->nbsteps, ispai->max, ispai->maxnew, ispai->block_size, ispai->cache_size, ispai->verbose));
9101a81e61SBarry Smith 
929566063dSJacob Faibussowitsch   PetscCall(ConvertMatrixToMat(PetscObjectComm((PetscObject)pc), ispai->M, &ispai->PM));
9301a81e61SBarry Smith 
9401a81e61SBarry Smith   /* free the SPAI matrices */
9501a81e61SBarry Smith   sp_free_matrix(ispai->B);
9601a81e61SBarry Smith   sp_free_matrix(ispai->M);
9701a81e61SBarry Smith   PetscFunctionReturn(0);
9801a81e61SBarry Smith }
9901a81e61SBarry Smith 
1009371c9d4SSatish Balay static PetscErrorCode PCApply_SPAI(PC pc, Vec xx, Vec y) {
10101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
10201a81e61SBarry Smith 
10301a81e61SBarry Smith   PetscFunctionBegin;
10401a81e61SBarry Smith   /* Now using PETSc's multiply */
1059566063dSJacob Faibussowitsch   PetscCall(MatMult(ispai->PM, xx, y));
10601a81e61SBarry Smith   PetscFunctionReturn(0);
10701a81e61SBarry Smith }
10801a81e61SBarry Smith 
1099371c9d4SSatish Balay static PetscErrorCode PCMatApply_SPAI(PC pc, Mat X, Mat Y) {
11048e38a8aSPierre Jolivet   PC_SPAI *ispai = (PC_SPAI *)pc->data;
11148e38a8aSPierre Jolivet 
11248e38a8aSPierre Jolivet   PetscFunctionBegin;
11348e38a8aSPierre Jolivet   /* Now using PETSc's multiply */
1149566063dSJacob Faibussowitsch   PetscCall(MatMatMult(ispai->PM, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
11548e38a8aSPierre Jolivet   PetscFunctionReturn(0);
11648e38a8aSPierre Jolivet }
11748e38a8aSPierre Jolivet 
1189371c9d4SSatish Balay static PetscErrorCode PCDestroy_SPAI(PC pc) {
11901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
12001a81e61SBarry Smith 
12101a81e61SBarry Smith   PetscFunctionBegin;
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ispai->PM));
1239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&(ispai->comm_spai)));
1249566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1252e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", NULL));
1262e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", NULL));
1272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", NULL));
1282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", NULL));
1292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", NULL));
1302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", NULL));
1312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", NULL));
1322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", NULL));
13301a81e61SBarry Smith   PetscFunctionReturn(0);
13401a81e61SBarry Smith }
13501a81e61SBarry Smith 
1369371c9d4SSatish Balay static PetscErrorCode PCView_SPAI(PC pc, PetscViewer viewer) {
13701a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
138ace3abfcSBarry Smith   PetscBool iascii;
13901a81e61SBarry Smith 
14001a81e61SBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
14201a81e61SBarry Smith   if (iascii) {
143b03515a0SUmberto Zerbinati     PetscCall(PetscViewerASCIIPrintf(viewer, "    epsilon %g\n", ispai->epsilon));
1449566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    nbsteps %d\n", ispai->nbsteps));
1459566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    max %d\n", ispai->max));
1469566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    maxnew %d\n", ispai->maxnew));
1479566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    block_size %d\n", ispai->block_size));
1489566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    cache_size %d\n", ispai->cache_size));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    verbose %d\n", ispai->verbose));
1509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    sp %d\n", ispai->sp));
15101a81e61SBarry Smith   }
15201a81e61SBarry Smith   PetscFunctionReturn(0);
15301a81e61SBarry Smith }
15401a81e61SBarry Smith 
1559371c9d4SSatish Balay static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1) {
15601a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1575fd66863SKarl Rupp 
15801a81e61SBarry Smith   PetscFunctionBegin;
159b03515a0SUmberto Zerbinati   ispai->epsilon = (double)epsilon1;
16001a81e61SBarry Smith   PetscFunctionReturn(0);
16101a81e61SBarry Smith }
16201a81e61SBarry Smith 
1639371c9d4SSatish Balay static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc, PetscInt nbsteps1) {
16401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1655fd66863SKarl Rupp 
16601a81e61SBarry Smith   PetscFunctionBegin;
167b03515a0SUmberto Zerbinati   ispai->nbsteps = (int)nbsteps1;
16801a81e61SBarry Smith   PetscFunctionReturn(0);
16901a81e61SBarry Smith }
17001a81e61SBarry Smith 
17101a81e61SBarry Smith /* added 1/7/99 g.h. */
1729371c9d4SSatish Balay static PetscErrorCode PCSPAISetMax_SPAI(PC pc, PetscInt max1) {
17301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1745fd66863SKarl Rupp 
17501a81e61SBarry Smith   PetscFunctionBegin;
176b03515a0SUmberto Zerbinati   ispai->max = (int)max1;
17701a81e61SBarry Smith   PetscFunctionReturn(0);
17801a81e61SBarry Smith }
17901a81e61SBarry Smith 
1809371c9d4SSatish Balay static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc, PetscInt maxnew1) {
18101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1825fd66863SKarl Rupp 
18301a81e61SBarry Smith   PetscFunctionBegin;
184b03515a0SUmberto Zerbinati   ispai->maxnew = (int)maxnew1;
18501a81e61SBarry Smith   PetscFunctionReturn(0);
18601a81e61SBarry Smith }
18701a81e61SBarry Smith 
1889371c9d4SSatish Balay static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc, PetscInt block_size1) {
18901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1905fd66863SKarl Rupp 
19101a81e61SBarry Smith   PetscFunctionBegin;
192b03515a0SUmberto Zerbinati   ispai->block_size = (int)block_size1;
19301a81e61SBarry Smith   PetscFunctionReturn(0);
19401a81e61SBarry Smith }
19501a81e61SBarry Smith 
1969371c9d4SSatish Balay static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc, PetscInt cache_size) {
19701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1985fd66863SKarl Rupp 
19901a81e61SBarry Smith   PetscFunctionBegin;
200b03515a0SUmberto Zerbinati   ispai->cache_size = (int)cache_size;
20101a81e61SBarry Smith   PetscFunctionReturn(0);
20201a81e61SBarry Smith }
20301a81e61SBarry Smith 
2049371c9d4SSatish Balay static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc, PetscInt verbose) {
20501a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2065fd66863SKarl Rupp 
20701a81e61SBarry Smith   PetscFunctionBegin;
208b03515a0SUmberto Zerbinati   ispai->verbose = (int)verbose;
20901a81e61SBarry Smith   PetscFunctionReturn(0);
21001a81e61SBarry Smith }
21101a81e61SBarry Smith 
2129371c9d4SSatish Balay static PetscErrorCode PCSPAISetSp_SPAI(PC pc, PetscInt sp) {
21301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2145fd66863SKarl Rupp 
21501a81e61SBarry Smith   PetscFunctionBegin;
216b03515a0SUmberto Zerbinati   ispai->sp = (int)sp;
21701a81e61SBarry Smith   PetscFunctionReturn(0);
21801a81e61SBarry Smith }
21901a81e61SBarry Smith 
22001a81e61SBarry Smith /*@
221*f1580f4eSBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the `PCSPAI` preconditioner
22201a81e61SBarry Smith 
22301a81e61SBarry Smith   Input Parameters:
22401a81e61SBarry Smith + pc - the preconditioner
22501a81e61SBarry Smith - eps - epsilon (default .4)
22601a81e61SBarry Smith 
227*f1580f4eSBarry Smith   Note:
22895452b02SPatrick Sanan     Espilon must be between 0 and 1. It controls the
22901a81e61SBarry Smith                  quality of the approximation of M to the inverse of
23001a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
23101a81e61SBarry Smith                  fill, and usually better preconditioners. In many
23201a81e61SBarry Smith                  cases the best choice of epsilon is the one that
23301a81e61SBarry Smith                  divides the total solution time equally between the
23401a81e61SBarry Smith                  preconditioner and the solver.
23501a81e61SBarry Smith 
23601a81e61SBarry Smith   Level: intermediate
23701a81e61SBarry Smith 
238db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
23901a81e61SBarry Smith   @*/
2409371c9d4SSatish Balay PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1) {
24101a81e61SBarry Smith   PetscFunctionBegin;
242b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1));
24301a81e61SBarry Smith   PetscFunctionReturn(0);
24401a81e61SBarry Smith }
24501a81e61SBarry Smith 
24601a81e61SBarry Smith /*@
24701a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
248*f1580f4eSBarry Smith         the `PCSPAI` preconditioner
24901a81e61SBarry Smith 
25001a81e61SBarry Smith   Input Parameters:
25101a81e61SBarry Smith + pc - the preconditioner
25201a81e61SBarry Smith - n - number of steps (default 5)
25301a81e61SBarry Smith 
254*f1580f4eSBarry Smith   Note:
255*f1580f4eSBarry Smith     `PCSPAI` constructs to approximation to every column of
25601a81e61SBarry Smith                  the exact inverse of A in a series of improvement
25701a81e61SBarry Smith                  steps. The quality of the approximation is determined
25801a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
25901a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
26001a81e61SBarry Smith                  uses the best approximation constructed so far.
26101a81e61SBarry Smith 
26201a81e61SBarry Smith   Level: intermediate
26301a81e61SBarry Smith 
264db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
26501a81e61SBarry Smith @*/
2669371c9d4SSatish Balay PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1) {
26701a81e61SBarry Smith   PetscFunctionBegin;
268b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1));
26901a81e61SBarry Smith   PetscFunctionReturn(0);
27001a81e61SBarry Smith }
27101a81e61SBarry Smith 
27201a81e61SBarry Smith /* added 1/7/99 g.h. */
27301a81e61SBarry Smith /*@
27401a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
275*f1580f4eSBarry Smith         the `PCSPAI` preconditioner
27601a81e61SBarry Smith 
27701a81e61SBarry Smith   Input Parameters:
27801a81e61SBarry Smith + pc - the preconditioner
27901a81e61SBarry Smith - n - size (default is 5000)
28001a81e61SBarry Smith 
28101a81e61SBarry Smith   Level: intermediate
28201a81e61SBarry Smith 
283db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
28401a81e61SBarry Smith @*/
2859371c9d4SSatish Balay PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1) {
28601a81e61SBarry Smith   PetscFunctionBegin;
287b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1));
28801a81e61SBarry Smith   PetscFunctionReturn(0);
28901a81e61SBarry Smith }
29001a81e61SBarry Smith 
29101a81e61SBarry Smith /*@
29201a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
293*f1580f4eSBarry Smith    in `PCSPAI` preconditioner
29401a81e61SBarry Smith 
29501a81e61SBarry Smith   Input Parameters:
29601a81e61SBarry Smith + pc - the preconditioner
29701a81e61SBarry Smith - n - maximum number (default 5)
29801a81e61SBarry Smith 
29901a81e61SBarry Smith   Level: intermediate
30001a81e61SBarry Smith 
301db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
30201a81e61SBarry Smith @*/
3039371c9d4SSatish Balay PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1) {
30401a81e61SBarry Smith   PetscFunctionBegin;
305b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1));
30601a81e61SBarry Smith   PetscFunctionReturn(0);
30701a81e61SBarry Smith }
30801a81e61SBarry Smith 
30901a81e61SBarry Smith /*@
310*f1580f4eSBarry Smith   PCSPAISetBlockSize - set the block size for the `PCSPAI` preconditioner
31101a81e61SBarry Smith 
31201a81e61SBarry Smith   Input Parameters:
31301a81e61SBarry Smith + pc - the preconditioner
31401a81e61SBarry Smith - n - block size (default 1)
31501a81e61SBarry Smith 
31695452b02SPatrick Sanan   Notes:
31795452b02SPatrick Sanan     A block
31801a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
31901a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
32001a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
32101a81e61SBarry Smith                  variable sized blocks, which are determined by
32201a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
32301a81e61SBarry Smith                  This can be very effective for finite-element
32401a81e61SBarry Smith                  matrices.
32501a81e61SBarry Smith 
32601a81e61SBarry Smith                  SPAI will convert A to block form, use a block
32701a81e61SBarry Smith                  version of the preconditioner algorithm, and then
32801a81e61SBarry Smith                  convert the result back to scalar form.
32901a81e61SBarry Smith 
33001a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
33101a81e61SBarry Smith                  can lead to very significant improvement in
33201a81e61SBarry Smith                  performance.
33301a81e61SBarry Smith 
33401a81e61SBarry Smith   Level: intermediate
33501a81e61SBarry Smith 
336db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
33701a81e61SBarry Smith @*/
3389371c9d4SSatish Balay PetscErrorCode PCSPAISetBlockSize(PC pc, PetscInt block_size1) {
33901a81e61SBarry Smith   PetscFunctionBegin;
340b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetBlockSize_C", (PC, PetscInt), (pc, block_size1));
34101a81e61SBarry Smith   PetscFunctionReturn(0);
34201a81e61SBarry Smith }
34301a81e61SBarry Smith 
34401a81e61SBarry Smith /*@
345*f1580f4eSBarry Smith   PCSPAISetCacheSize - specify cache size in the `PCSPAI` preconditioner
34601a81e61SBarry Smith 
34701a81e61SBarry Smith   Input Parameters:
34801a81e61SBarry Smith + pc - the preconditioner
34901a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
35001a81e61SBarry Smith 
351*f1580f4eSBarry Smith   Note:
352*f1580f4eSBarry Smith     `PCSPAI` uses a hash table to cache messages and avoid
35301a81e61SBarry Smith                  redundant communication. If suggest always using
35401a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
35501a81e61SBarry Smith                  version.
35601a81e61SBarry Smith 
35701a81e61SBarry Smith   Level: intermediate
35801a81e61SBarry Smith 
359db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
36001a81e61SBarry Smith @*/
3619371c9d4SSatish Balay PetscErrorCode PCSPAISetCacheSize(PC pc, PetscInt cache_size) {
36201a81e61SBarry Smith   PetscFunctionBegin;
363b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetCacheSize_C", (PC, PetscInt), (pc, cache_size));
36401a81e61SBarry Smith   PetscFunctionReturn(0);
36501a81e61SBarry Smith }
36601a81e61SBarry Smith 
36701a81e61SBarry Smith /*@
368*f1580f4eSBarry Smith   PCSPAISetVerbose - verbosity level for the `PCSPAI` preconditioner
36901a81e61SBarry Smith 
37001a81e61SBarry Smith   Input Parameters:
37101a81e61SBarry Smith + pc - the preconditioner
37201a81e61SBarry Smith - n - level (default 1)
37301a81e61SBarry Smith 
374*f1580f4eSBarry Smith   Note:
37595452b02SPatrick Sanan     print parameters, timings and matrix statistics
37601a81e61SBarry Smith 
37701a81e61SBarry Smith   Level: intermediate
37801a81e61SBarry Smith 
379db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
38001a81e61SBarry Smith @*/
3819371c9d4SSatish Balay PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose) {
38201a81e61SBarry Smith   PetscFunctionBegin;
383b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose));
38401a81e61SBarry Smith   PetscFunctionReturn(0);
38501a81e61SBarry Smith }
38601a81e61SBarry Smith 
38701a81e61SBarry Smith /*@
388*f1580f4eSBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the `PCSPAI` preconditioner
38901a81e61SBarry Smith 
39001a81e61SBarry Smith   Input Parameters:
39101a81e61SBarry Smith + pc - the preconditioner
39201a81e61SBarry Smith - n - 0 or 1
39301a81e61SBarry Smith 
394*f1580f4eSBarry Smith   Note:
39595452b02SPatrick Sanan     If A has a symmetric nonzero pattern use -sp 1 to
39601a81e61SBarry Smith                  improve performance by eliminating some communication
39701a81e61SBarry Smith                  in the parallel version. Even if A does not have a
39801a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
39901a81e61SBarry Smith                  results, but the code will not follow the published
40001a81e61SBarry Smith                  SPAI algorithm exactly.
40101a81e61SBarry Smith 
40201a81e61SBarry Smith   Level: intermediate
40301a81e61SBarry Smith 
404db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
40501a81e61SBarry Smith @*/
4069371c9d4SSatish Balay PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp) {
40701a81e61SBarry Smith   PetscFunctionBegin;
408b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp));
40901a81e61SBarry Smith   PetscFunctionReturn(0);
41001a81e61SBarry Smith }
41101a81e61SBarry Smith 
4129371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_SPAI(PC pc, PetscOptionItems *PetscOptionsObject) {
41301a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
41401a81e61SBarry Smith   int       nbsteps1, max1, maxnew1, block_size1, cache_size, verbose, sp;
41501a81e61SBarry Smith   double    epsilon1;
416ace3abfcSBarry Smith   PetscBool flg;
41701a81e61SBarry Smith 
41801a81e61SBarry Smith   PetscFunctionBegin;
419d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SPAI options");
4209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_spai_epsilon", "", "PCSPAISetEpsilon", ispai->epsilon, &epsilon1, &flg));
4211baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetEpsilon(pc, epsilon1));
4229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_nbsteps", "", "PCSPAISetNBSteps", ispai->nbsteps, &nbsteps1, &flg));
4231baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetNBSteps(pc, nbsteps1));
42401a81e61SBarry Smith   /* added 1/7/99 g.h. */
4259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_max", "", "PCSPAISetMax", ispai->max, &max1, &flg));
4261baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMax(pc, max1));
4279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_maxnew", "", "PCSPAISetMaxNew", ispai->maxnew, &maxnew1, &flg));
4281baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMaxNew(pc, maxnew1));
4299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_block_size", "", "PCSPAISetBlockSize", ispai->block_size, &block_size1, &flg));
4301baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetBlockSize(pc, block_size1));
4319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_cache_size", "", "PCSPAISetCacheSize", ispai->cache_size, &cache_size, &flg));
4321baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetCacheSize(pc, cache_size));
4339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_verbose", "", "PCSPAISetVerbose", ispai->verbose, &verbose, &flg));
4341baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetVerbose(pc, verbose));
4359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_sp", "", "PCSPAISetSp", ispai->sp, &sp, &flg));
4361baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetSp(pc, sp));
437d0609cedSBarry Smith   PetscOptionsHeadEnd();
43801a81e61SBarry Smith   PetscFunctionReturn(0);
43901a81e61SBarry Smith }
44001a81e61SBarry Smith 
44101a81e61SBarry Smith /*MC
442*f1580f4eSBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method
44301a81e61SBarry Smith 
44401a81e61SBarry Smith    Options Database Keys:
445c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
44601a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
44701a81e61SBarry Smith .  -pc_spai_max <m> - set max
44801a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
44901a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
45001a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
45101a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
45201a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
45301a81e61SBarry Smith 
45401a81e61SBarry Smith    Level: beginner
45501a81e61SBarry Smith 
456*f1580f4eSBarry Smith    Note:
457*f1580f4eSBarry Smith     This only works with `MATAIJ` matrices.
458*f1580f4eSBarry Smith 
459*f1580f4eSBarry Smith    References:
460*f1580f4eSBarry Smith  . * -  Grote and Barnard (SIAM J. Sci. Comput.; vol 18, nr 3)
461*f1580f4eSBarry Smith 
462db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
463db781477SPatrick Sanan           `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`,
464*f1580f4eSBarry Smith           `PCSPAISetVerbose()`, `PCSPAISetSp()`, `PCSPAISetNBSteps()`, `PCSPAISetCacheSize()`
46501a81e61SBarry Smith M*/
46601a81e61SBarry Smith 
4679371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc) {
46801a81e61SBarry Smith   PC_SPAI *ispai;
46901a81e61SBarry Smith 
47001a81e61SBarry Smith   PetscFunctionBegin;
4719566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &ispai));
4722a4967abSBarry Smith   pc->data = ispai;
47301a81e61SBarry Smith 
47401a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
47501a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
47648e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_SPAI;
47701a81e61SBarry Smith   pc->ops->applyrichardson = 0;
47801a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
47901a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
48001a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
48101a81e61SBarry Smith 
48201a81e61SBarry Smith   ispai->epsilon    = .4;
48301a81e61SBarry Smith   ispai->nbsteps    = 5;
48401a81e61SBarry Smith   ispai->max        = 5000;
48501a81e61SBarry Smith   ispai->maxnew     = 5;
48601a81e61SBarry Smith   ispai->block_size = 1;
48701a81e61SBarry Smith   ispai->cache_size = 5;
48801a81e61SBarry Smith   ispai->verbose    = 0;
48901a81e61SBarry Smith 
49001a81e61SBarry Smith   ispai->sp = 1;
4919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &(ispai->comm_spai)));
49201a81e61SBarry Smith 
4939566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", PCSPAISetEpsilon_SPAI));
4949566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", PCSPAISetNBSteps_SPAI));
4959566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", PCSPAISetMax_SPAI));
4969566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", PCSPAISetMaxNew_SPAI));
4979566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", PCSPAISetBlockSize_SPAI));
4989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", PCSPAISetCacheSize_SPAI));
4999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", PCSPAISetVerbose_SPAI));
5009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", PCSPAISetSp_SPAI));
50101a81e61SBarry Smith   PetscFunctionReturn(0);
50201a81e61SBarry Smith }
50301a81e61SBarry Smith 
50401a81e61SBarry Smith /*
50501a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
50601a81e61SBarry Smith */
5079371c9d4SSatish Balay PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A, Mat AT, matrix **B) {
50801a81e61SBarry Smith   matrix                  *M;
50901a81e61SBarry Smith   int                      i, j, col;
51001a81e61SBarry Smith   int                      row_indx;
51101a81e61SBarry Smith   int                      len, pe, local_indx, start_indx;
51201a81e61SBarry Smith   int                     *mapping;
51301a81e61SBarry Smith   const int               *cols;
51401a81e61SBarry Smith   const double            *vals;
5152122902bSSatish Balay   int                      n, mnl, nnl, nz, rstart, rend;
5162a4967abSBarry Smith   PetscMPIInt              size, rank;
51701a81e61SBarry Smith   struct compressed_lines *rows;
51801a81e61SBarry Smith 
51901a81e61SBarry Smith   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5229566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &n, &n));
5239566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &mnl, &nnl));
52401a81e61SBarry Smith 
52501a81e61SBarry Smith   /*
52601a81e61SBarry Smith     not sure why a barrier is required. commenting out
5279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(comm));
52801a81e61SBarry Smith   */
52901a81e61SBarry Smith 
53029b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
53101a81e61SBarry Smith 
53201a81e61SBarry Smith   M->n              = n;
53301a81e61SBarry Smith   M->bs             = 1;
53401a81e61SBarry Smith   M->max_block_size = 1;
53501a81e61SBarry Smith 
53601a81e61SBarry Smith   M->mnls          = (int *)malloc(sizeof(int) * size);
53701a81e61SBarry Smith   M->start_indices = (int *)malloc(sizeof(int) * size);
53801a81e61SBarry Smith   M->pe            = (int *)malloc(sizeof(int) * n);
53901a81e61SBarry Smith   M->block_sizes   = (int *)malloc(sizeof(int) * n);
54001a81e61SBarry Smith   for (i = 0; i < n; i++) M->block_sizes[i] = 1;
54101a81e61SBarry Smith 
5429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnl, 1, MPI_INT, M->mnls, 1, MPI_INT, comm));
54301a81e61SBarry Smith 
54401a81e61SBarry Smith   M->start_indices[0] = 0;
5452fa5cd67SKarl Rupp   for (i = 1; i < size; i++) M->start_indices[i] = M->start_indices[i - 1] + M->mnls[i - 1];
54601a81e61SBarry Smith 
54701a81e61SBarry Smith   M->mnl            = M->mnls[M->myid];
54801a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
54901a81e61SBarry Smith 
55001a81e61SBarry Smith   for (i = 0; i < size; i++) {
55101a81e61SBarry Smith     start_indx = M->start_indices[i];
5522fa5cd67SKarl Rupp     for (j = 0; j < M->mnls[i]; j++) M->pe[start_indx + j] = i;
55301a81e61SBarry Smith   }
55401a81e61SBarry Smith 
55501a81e61SBarry Smith   if (AT) {
5562f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 1);
55701a81e61SBarry Smith   } else {
5582f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 0);
55901a81e61SBarry Smith   }
56001a81e61SBarry Smith 
56101a81e61SBarry Smith   rows = M->lines;
56201a81e61SBarry Smith 
56301a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
5649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M->n, &mapping));
56501a81e61SBarry Smith   pe         = 0;
56601a81e61SBarry Smith   local_indx = 0;
56701a81e61SBarry Smith   for (i = 0; i < M->n; i++) {
56801a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
56901a81e61SBarry Smith       pe++;
57001a81e61SBarry Smith       local_indx = 0;
57101a81e61SBarry Smith     }
57201a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
57301a81e61SBarry Smith     local_indx++;
57401a81e61SBarry Smith   }
57501a81e61SBarry Smith 
57601a81e61SBarry Smith   /************** Set up the row structure *****************/
57701a81e61SBarry Smith 
5789566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
57901a81e61SBarry Smith   for (i = rstart; i < rend; i++) {
58001a81e61SBarry Smith     row_indx = i - rstart;
5819566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nz, &cols, &vals));
5822122902bSSatish Balay     /* allocate buffers */
5832122902bSSatish Balay     rows->ptrs[row_indx] = (int *)malloc(nz * sizeof(int));
5842122902bSSatish Balay     rows->A[row_indx]    = (double *)malloc(nz * sizeof(double));
5852122902bSSatish Balay     /* copy the matrix */
58601a81e61SBarry Smith     for (j = 0; j < nz; j++) {
58701a81e61SBarry Smith       col = cols[j];
58801a81e61SBarry Smith       len = rows->len[row_indx]++;
5892fa5cd67SKarl Rupp 
59001a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
59101a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
59201a81e61SBarry Smith     }
59301a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
5942fa5cd67SKarl Rupp 
5959566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nz, &cols, &vals));
59601a81e61SBarry Smith   }
59701a81e61SBarry Smith 
59801a81e61SBarry Smith   /************** Set up the column structure *****************/
59901a81e61SBarry Smith 
60001a81e61SBarry Smith   if (AT) {
60101a81e61SBarry Smith     for (i = rstart; i < rend; i++) {
60201a81e61SBarry Smith       row_indx = i - rstart;
6039566063dSJacob Faibussowitsch       PetscCall(MatGetRow(AT, i, &nz, &cols, &vals));
6042122902bSSatish Balay       /* allocate buffers */
6052122902bSSatish Balay       rows->rptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6062122902bSSatish Balay       /* copy the matrix (i.e., the structure) */
60701a81e61SBarry Smith       for (j = 0; j < nz; j++) {
60801a81e61SBarry Smith         col = cols[j];
60901a81e61SBarry Smith         len = rows->rlen[row_indx]++;
6102fa5cd67SKarl Rupp 
61101a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
61201a81e61SBarry Smith       }
6139566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(AT, i, &nz, &cols, &vals));
61401a81e61SBarry Smith     }
61501a81e61SBarry Smith   }
61601a81e61SBarry Smith 
6179566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping));
61801a81e61SBarry Smith 
61901a81e61SBarry Smith   order_pointers(M);
62001a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
62101a81e61SBarry Smith   *B       = M;
62201a81e61SBarry Smith   PetscFunctionReturn(0);
62301a81e61SBarry Smith }
62401a81e61SBarry Smith 
62501a81e61SBarry Smith /*
62601a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
627df3898eeSBarry Smith    This assumes that the SPAI matrix B is stored in
62801a81e61SBarry Smith    COMPRESSED-ROW format.
62901a81e61SBarry Smith */
6309371c9d4SSatish Balay PetscErrorCode ConvertMatrixToMat(MPI_Comm comm, matrix *B, Mat *PB) {
6314b6b5fe1SBarry Smith   PetscMPIInt size, rank;
63201a81e61SBarry Smith   int         m, n, M, N;
63301a81e61SBarry Smith   int         d_nz, o_nz;
63401a81e61SBarry Smith   int        *d_nnz, *o_nnz;
63501a81e61SBarry Smith   int         i, k, global_row, global_col, first_diag_col, last_diag_col;
63601a81e61SBarry Smith   PetscScalar val;
63701a81e61SBarry Smith 
63801a81e61SBarry Smith   PetscFunctionBegin;
6399566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
64101a81e61SBarry Smith 
64201a81e61SBarry Smith   m = n = B->mnls[rank];
64301a81e61SBarry Smith   d_nz = o_nz = 0;
64401a81e61SBarry Smith 
64506946f3aSJose E. Roman   /* Determine preallocation for MatCreateAIJ */
6469566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &d_nnz));
6479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &o_nnz));
64801a81e61SBarry Smith   for (i = 0; i < m; i++) d_nnz[i] = o_nnz[i] = 0;
64901a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
65001a81e61SBarry Smith   last_diag_col  = first_diag_col + B->mnls[rank];
65101a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
65201a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
65301a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
654db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
655db4deed7SKarl Rupp       else o_nnz[i]++;
65601a81e61SBarry Smith     }
65701a81e61SBarry Smith   }
65801a81e61SBarry Smith 
65901a81e61SBarry Smith   M = N = B->n;
66001a81e61SBarry Smith   /* Here we only know how to create AIJ format */
6619566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, PB));
6629566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*PB, m, n, M, N));
6639566063dSJacob Faibussowitsch   PetscCall(MatSetType(*PB, MATAIJ));
6649566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*PB, d_nz, d_nnz));
6659566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*PB, d_nz, d_nnz, o_nz, o_nnz));
66601a81e61SBarry Smith 
66701a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
66801a81e61SBarry Smith     global_row = B->start_indices[rank] + i;
66901a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
67001a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
6712fa5cd67SKarl Rupp 
67201a81e61SBarry Smith       val = B->lines->A[i][k];
6739566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*PB, 1, &global_row, 1, &global_col, &val, ADD_VALUES));
67401a81e61SBarry Smith     }
67501a81e61SBarry Smith   }
67601a81e61SBarry Smith 
6779566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
6789566063dSJacob Faibussowitsch   PetscCall(PetscFree(o_nnz));
67901a81e61SBarry Smith 
6809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*PB, MAT_FINAL_ASSEMBLY));
6819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*PB, MAT_FINAL_ASSEMBLY));
68201a81e61SBarry Smith   PetscFunctionReturn(0);
68301a81e61SBarry Smith }
68401a81e61SBarry Smith 
68501a81e61SBarry Smith /*
68601a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
68701a81e61SBarry Smith */
6889371c9d4SSatish Balay PetscErrorCode ConvertVectorToVec(MPI_Comm comm, vector *v, Vec *Pv) {
6894b6b5fe1SBarry Smith   PetscMPIInt size, rank;
6904b6b5fe1SBarry Smith   int         m, M, i, *mnls, *start_indices, *global_indices;
69101a81e61SBarry Smith 
69201a81e61SBarry Smith   PetscFunctionBegin;
6939566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6949566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
69501a81e61SBarry Smith 
69601a81e61SBarry Smith   m = v->mnl;
69701a81e61SBarry Smith   M = v->n;
69801a81e61SBarry Smith 
6999566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, m, M, Pv));
70001a81e61SBarry Smith 
7019566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &mnls));
7029566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&v->mnl, 1, MPI_INT, mnls, 1, MPI_INT, comm));
70301a81e61SBarry Smith 
7049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &start_indices));
7052fa5cd67SKarl Rupp 
70601a81e61SBarry Smith   start_indices[0] = 0;
7072fa5cd67SKarl Rupp   for (i = 1; i < size; i++) start_indices[i] = start_indices[i - 1] + mnls[i - 1];
70801a81e61SBarry Smith 
7099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(v->mnl, &global_indices));
7102fa5cd67SKarl Rupp   for (i = 0; i < v->mnl; i++) global_indices[i] = start_indices[rank] + i;
71101a81e61SBarry Smith 
7129566063dSJacob Faibussowitsch   PetscCall(PetscFree(mnls));
7139566063dSJacob Faibussowitsch   PetscCall(PetscFree(start_indices));
71401a81e61SBarry Smith 
7159566063dSJacob Faibussowitsch   PetscCall(VecSetValues(*Pv, v->mnl, global_indices, v->v, INSERT_VALUES));
7169566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*Pv));
7179566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*Pv));
71801a81e61SBarry Smith 
7199566063dSJacob Faibussowitsch   PetscCall(PetscFree(global_indices));
72001a81e61SBarry Smith   PetscFunctionReturn(0);
72101a81e61SBarry Smith }
722