xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 
58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SPAI(PC pc)
59d71ae5a4SJacob Faibussowitsch {
6001a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
6101a81e61SBarry Smith   Mat      AT;
6201a81e61SBarry Smith 
6301a81e61SBarry Smith   PetscFunctionBegin;
6401a81e61SBarry Smith   init_SPAI();
6501a81e61SBarry Smith 
6601a81e61SBarry Smith   if (ispai->sp) {
679566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, pc->pmat, &ispai->B));
6801a81e61SBarry Smith   } else {
6901a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
709566063dSJacob Faibussowitsch     PetscCall(MatTranspose(pc->pmat, MAT_INITIAL_MATRIX, &AT));
719566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, AT, &ispai->B));
729566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
7301a81e61SBarry Smith   }
7401a81e61SBarry Smith 
7501a81e61SBarry Smith   /* Destroy the transpose */
7601a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
7701a81e61SBarry Smith 
7801a81e61SBarry Smith   /* construct SPAI preconditioner */
7901a81e61SBarry Smith   /* FILE *messages */    /* file for warning messages */
8001a81e61SBarry Smith   /* double epsilon */    /* tolerance */
8101a81e61SBarry Smith   /* int nbsteps */       /* max number of "improvement" steps per line */
8201a81e61SBarry Smith   /* int max */           /* max dimensions of is_I, q, etc. */
8301a81e61SBarry Smith   /* int maxnew */        /* max number of new entries per step */
84*d5b43468SJose E. Roman   /* int block_size */    /* block_size == 1 specifies scalar elements
8501a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
8601a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
872fa5cd67SKarl Rupp   /* int cache_size */    /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */
8801a81e61SBarry Smith   /* int    verbose    */ /* verbose == 0 specifies that SPAI is silent
8901a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9001a81e61SBarry Smith 
919371c9d4SSatish Balay   PetscCall(bspai(ispai->B, &ispai->M, stdout, ispai->epsilon, ispai->nbsteps, ispai->max, ispai->maxnew, ispai->block_size, ispai->cache_size, ispai->verbose));
9201a81e61SBarry Smith 
939566063dSJacob Faibussowitsch   PetscCall(ConvertMatrixToMat(PetscObjectComm((PetscObject)pc), ispai->M, &ispai->PM));
9401a81e61SBarry Smith 
9501a81e61SBarry Smith   /* free the SPAI matrices */
9601a81e61SBarry Smith   sp_free_matrix(ispai->B);
9701a81e61SBarry Smith   sp_free_matrix(ispai->M);
9801a81e61SBarry Smith   PetscFunctionReturn(0);
9901a81e61SBarry Smith }
10001a81e61SBarry Smith 
101d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SPAI(PC pc, Vec xx, Vec y)
102d71ae5a4SJacob Faibussowitsch {
10301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
10401a81e61SBarry Smith 
10501a81e61SBarry Smith   PetscFunctionBegin;
10601a81e61SBarry Smith   /* Now using PETSc's multiply */
1079566063dSJacob Faibussowitsch   PetscCall(MatMult(ispai->PM, xx, y));
10801a81e61SBarry Smith   PetscFunctionReturn(0);
10901a81e61SBarry Smith }
11001a81e61SBarry Smith 
111d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SPAI(PC pc, Mat X, Mat Y)
112d71ae5a4SJacob Faibussowitsch {
11348e38a8aSPierre Jolivet   PC_SPAI *ispai = (PC_SPAI *)pc->data;
11448e38a8aSPierre Jolivet 
11548e38a8aSPierre Jolivet   PetscFunctionBegin;
11648e38a8aSPierre Jolivet   /* Now using PETSc's multiply */
1179566063dSJacob Faibussowitsch   PetscCall(MatMatMult(ispai->PM, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
11848e38a8aSPierre Jolivet   PetscFunctionReturn(0);
11948e38a8aSPierre Jolivet }
12048e38a8aSPierre Jolivet 
121d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SPAI(PC pc)
122d71ae5a4SJacob Faibussowitsch {
12301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
12401a81e61SBarry Smith 
12501a81e61SBarry Smith   PetscFunctionBegin;
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ispai->PM));
1279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&(ispai->comm_spai)));
1289566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", NULL));
1302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", NULL));
1312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", NULL));
1322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", NULL));
1332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", NULL));
1342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", NULL));
1352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", NULL));
1362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", NULL));
13701a81e61SBarry Smith   PetscFunctionReturn(0);
13801a81e61SBarry Smith }
13901a81e61SBarry Smith 
140d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SPAI(PC pc, PetscViewer viewer)
141d71ae5a4SJacob Faibussowitsch {
14201a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
143ace3abfcSBarry Smith   PetscBool iascii;
14401a81e61SBarry Smith 
14501a81e61SBarry Smith   PetscFunctionBegin;
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
14701a81e61SBarry Smith   if (iascii) {
148b03515a0SUmberto Zerbinati     PetscCall(PetscViewerASCIIPrintf(viewer, "    epsilon %g\n", ispai->epsilon));
1499566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    nbsteps %d\n", ispai->nbsteps));
1509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    max %d\n", ispai->max));
1519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    maxnew %d\n", ispai->maxnew));
1529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    block_size %d\n", ispai->block_size));
1539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    cache_size %d\n", ispai->cache_size));
1549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    verbose %d\n", ispai->verbose));
1559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    sp %d\n", ispai->sp));
15601a81e61SBarry Smith   }
15701a81e61SBarry Smith   PetscFunctionReturn(0);
15801a81e61SBarry Smith }
15901a81e61SBarry Smith 
160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1)
161d71ae5a4SJacob Faibussowitsch {
16201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1635fd66863SKarl Rupp 
16401a81e61SBarry Smith   PetscFunctionBegin;
165b03515a0SUmberto Zerbinati   ispai->epsilon = (double)epsilon1;
16601a81e61SBarry Smith   PetscFunctionReturn(0);
16701a81e61SBarry Smith }
16801a81e61SBarry Smith 
169d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc, PetscInt nbsteps1)
170d71ae5a4SJacob Faibussowitsch {
17101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1725fd66863SKarl Rupp 
17301a81e61SBarry Smith   PetscFunctionBegin;
174b03515a0SUmberto Zerbinati   ispai->nbsteps = (int)nbsteps1;
17501a81e61SBarry Smith   PetscFunctionReturn(0);
17601a81e61SBarry Smith }
17701a81e61SBarry Smith 
17801a81e61SBarry Smith /* added 1/7/99 g.h. */
179d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMax_SPAI(PC pc, PetscInt max1)
180d71ae5a4SJacob Faibussowitsch {
18101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1825fd66863SKarl Rupp 
18301a81e61SBarry Smith   PetscFunctionBegin;
184b03515a0SUmberto Zerbinati   ispai->max = (int)max1;
18501a81e61SBarry Smith   PetscFunctionReturn(0);
18601a81e61SBarry Smith }
18701a81e61SBarry Smith 
188d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc, PetscInt maxnew1)
189d71ae5a4SJacob Faibussowitsch {
19001a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1915fd66863SKarl Rupp 
19201a81e61SBarry Smith   PetscFunctionBegin;
193b03515a0SUmberto Zerbinati   ispai->maxnew = (int)maxnew1;
19401a81e61SBarry Smith   PetscFunctionReturn(0);
19501a81e61SBarry Smith }
19601a81e61SBarry Smith 
197d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc, PetscInt block_size1)
198d71ae5a4SJacob Faibussowitsch {
19901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2005fd66863SKarl Rupp 
20101a81e61SBarry Smith   PetscFunctionBegin;
202b03515a0SUmberto Zerbinati   ispai->block_size = (int)block_size1;
20301a81e61SBarry Smith   PetscFunctionReturn(0);
20401a81e61SBarry Smith }
20501a81e61SBarry Smith 
206d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc, PetscInt cache_size)
207d71ae5a4SJacob Faibussowitsch {
20801a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2095fd66863SKarl Rupp 
21001a81e61SBarry Smith   PetscFunctionBegin;
211b03515a0SUmberto Zerbinati   ispai->cache_size = (int)cache_size;
21201a81e61SBarry Smith   PetscFunctionReturn(0);
21301a81e61SBarry Smith }
21401a81e61SBarry Smith 
215d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc, PetscInt verbose)
216d71ae5a4SJacob Faibussowitsch {
21701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2185fd66863SKarl Rupp 
21901a81e61SBarry Smith   PetscFunctionBegin;
220b03515a0SUmberto Zerbinati   ispai->verbose = (int)verbose;
22101a81e61SBarry Smith   PetscFunctionReturn(0);
22201a81e61SBarry Smith }
22301a81e61SBarry Smith 
224d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetSp_SPAI(PC pc, PetscInt sp)
225d71ae5a4SJacob Faibussowitsch {
22601a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2275fd66863SKarl Rupp 
22801a81e61SBarry Smith   PetscFunctionBegin;
229b03515a0SUmberto Zerbinati   ispai->sp = (int)sp;
23001a81e61SBarry Smith   PetscFunctionReturn(0);
23101a81e61SBarry Smith }
23201a81e61SBarry Smith 
23301a81e61SBarry Smith /*@
234f1580f4eSBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the `PCSPAI` preconditioner
23501a81e61SBarry Smith 
23601a81e61SBarry Smith   Input Parameters:
23701a81e61SBarry Smith + pc - the preconditioner
23801a81e61SBarry Smith - eps - epsilon (default .4)
23901a81e61SBarry Smith 
240f1580f4eSBarry Smith   Note:
24195452b02SPatrick Sanan     Espilon must be between 0 and 1. It controls the
24201a81e61SBarry Smith                  quality of the approximation of M to the inverse of
24301a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
24401a81e61SBarry Smith                  fill, and usually better preconditioners. In many
24501a81e61SBarry Smith                  cases the best choice of epsilon is the one that
24601a81e61SBarry Smith                  divides the total solution time equally between the
24701a81e61SBarry Smith                  preconditioner and the solver.
24801a81e61SBarry Smith 
24901a81e61SBarry Smith   Level: intermediate
25001a81e61SBarry Smith 
251db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
25201a81e61SBarry Smith   @*/
253d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1)
254d71ae5a4SJacob Faibussowitsch {
25501a81e61SBarry Smith   PetscFunctionBegin;
256b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1));
25701a81e61SBarry Smith   PetscFunctionReturn(0);
25801a81e61SBarry Smith }
25901a81e61SBarry Smith 
26001a81e61SBarry Smith /*@
26101a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
262f1580f4eSBarry Smith         the `PCSPAI` preconditioner
26301a81e61SBarry Smith 
26401a81e61SBarry Smith   Input Parameters:
26501a81e61SBarry Smith + pc - the preconditioner
26601a81e61SBarry Smith - n - number of steps (default 5)
26701a81e61SBarry Smith 
268f1580f4eSBarry Smith   Note:
269f1580f4eSBarry Smith     `PCSPAI` constructs to approximation to every column of
27001a81e61SBarry Smith                  the exact inverse of A in a series of improvement
27101a81e61SBarry Smith                  steps. The quality of the approximation is determined
27201a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
27301a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
27401a81e61SBarry Smith                  uses the best approximation constructed so far.
27501a81e61SBarry Smith 
27601a81e61SBarry Smith   Level: intermediate
27701a81e61SBarry Smith 
278db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
27901a81e61SBarry Smith @*/
280d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1)
281d71ae5a4SJacob Faibussowitsch {
28201a81e61SBarry Smith   PetscFunctionBegin;
283b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1));
28401a81e61SBarry Smith   PetscFunctionReturn(0);
28501a81e61SBarry Smith }
28601a81e61SBarry Smith 
28701a81e61SBarry Smith /* added 1/7/99 g.h. */
28801a81e61SBarry Smith /*@
28901a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
290f1580f4eSBarry Smith         the `PCSPAI` preconditioner
29101a81e61SBarry Smith 
29201a81e61SBarry Smith   Input Parameters:
29301a81e61SBarry Smith + pc - the preconditioner
29401a81e61SBarry Smith - n - size (default is 5000)
29501a81e61SBarry Smith 
29601a81e61SBarry Smith   Level: intermediate
29701a81e61SBarry Smith 
298db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
29901a81e61SBarry Smith @*/
300d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1)
301d71ae5a4SJacob Faibussowitsch {
30201a81e61SBarry Smith   PetscFunctionBegin;
303b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1));
30401a81e61SBarry Smith   PetscFunctionReturn(0);
30501a81e61SBarry Smith }
30601a81e61SBarry Smith 
30701a81e61SBarry Smith /*@
30801a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
309f1580f4eSBarry Smith    in `PCSPAI` preconditioner
31001a81e61SBarry Smith 
31101a81e61SBarry Smith   Input Parameters:
31201a81e61SBarry Smith + pc - the preconditioner
31301a81e61SBarry Smith - n - maximum number (default 5)
31401a81e61SBarry Smith 
31501a81e61SBarry Smith   Level: intermediate
31601a81e61SBarry Smith 
317db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
31801a81e61SBarry Smith @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1)
320d71ae5a4SJacob Faibussowitsch {
32101a81e61SBarry Smith   PetscFunctionBegin;
322b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1));
32301a81e61SBarry Smith   PetscFunctionReturn(0);
32401a81e61SBarry Smith }
32501a81e61SBarry Smith 
32601a81e61SBarry Smith /*@
327f1580f4eSBarry Smith   PCSPAISetBlockSize - set the block size for the `PCSPAI` preconditioner
32801a81e61SBarry Smith 
32901a81e61SBarry Smith   Input Parameters:
33001a81e61SBarry Smith + pc - the preconditioner
33101a81e61SBarry Smith - n - block size (default 1)
33201a81e61SBarry Smith 
33395452b02SPatrick Sanan   Notes:
33495452b02SPatrick Sanan     A block
33501a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
33601a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
33701a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
33801a81e61SBarry Smith                  variable sized blocks, which are determined by
33901a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
34001a81e61SBarry Smith                  This can be very effective for finite-element
34101a81e61SBarry Smith                  matrices.
34201a81e61SBarry Smith 
34301a81e61SBarry Smith                  SPAI will convert A to block form, use a block
34401a81e61SBarry Smith                  version of the preconditioner algorithm, and then
34501a81e61SBarry Smith                  convert the result back to scalar form.
34601a81e61SBarry Smith 
34701a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
34801a81e61SBarry Smith                  can lead to very significant improvement in
34901a81e61SBarry Smith                  performance.
35001a81e61SBarry Smith 
35101a81e61SBarry Smith   Level: intermediate
35201a81e61SBarry Smith 
353db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
35401a81e61SBarry Smith @*/
355d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetBlockSize(PC pc, PetscInt block_size1)
356d71ae5a4SJacob Faibussowitsch {
35701a81e61SBarry Smith   PetscFunctionBegin;
358b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetBlockSize_C", (PC, PetscInt), (pc, block_size1));
35901a81e61SBarry Smith   PetscFunctionReturn(0);
36001a81e61SBarry Smith }
36101a81e61SBarry Smith 
36201a81e61SBarry Smith /*@
363f1580f4eSBarry Smith   PCSPAISetCacheSize - specify cache size in the `PCSPAI` preconditioner
36401a81e61SBarry Smith 
36501a81e61SBarry Smith   Input Parameters:
36601a81e61SBarry Smith + pc - the preconditioner
36701a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
36801a81e61SBarry Smith 
369f1580f4eSBarry Smith   Note:
370f1580f4eSBarry Smith     `PCSPAI` uses a hash table to cache messages and avoid
37101a81e61SBarry Smith                  redundant communication. If suggest always using
37201a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
37301a81e61SBarry Smith                  version.
37401a81e61SBarry Smith 
37501a81e61SBarry Smith   Level: intermediate
37601a81e61SBarry Smith 
377db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
37801a81e61SBarry Smith @*/
379d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetCacheSize(PC pc, PetscInt cache_size)
380d71ae5a4SJacob Faibussowitsch {
38101a81e61SBarry Smith   PetscFunctionBegin;
382b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetCacheSize_C", (PC, PetscInt), (pc, cache_size));
38301a81e61SBarry Smith   PetscFunctionReturn(0);
38401a81e61SBarry Smith }
38501a81e61SBarry Smith 
38601a81e61SBarry Smith /*@
387f1580f4eSBarry Smith   PCSPAISetVerbose - verbosity level for the `PCSPAI` preconditioner
38801a81e61SBarry Smith 
38901a81e61SBarry Smith   Input Parameters:
39001a81e61SBarry Smith + pc - the preconditioner
39101a81e61SBarry Smith - n - level (default 1)
39201a81e61SBarry Smith 
393f1580f4eSBarry Smith   Note:
39495452b02SPatrick Sanan     print parameters, timings and matrix statistics
39501a81e61SBarry Smith 
39601a81e61SBarry Smith   Level: intermediate
39701a81e61SBarry Smith 
398db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
39901a81e61SBarry Smith @*/
400d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose)
401d71ae5a4SJacob Faibussowitsch {
40201a81e61SBarry Smith   PetscFunctionBegin;
403b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose));
40401a81e61SBarry Smith   PetscFunctionReturn(0);
40501a81e61SBarry Smith }
40601a81e61SBarry Smith 
40701a81e61SBarry Smith /*@
408f1580f4eSBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the `PCSPAI` preconditioner
40901a81e61SBarry Smith 
41001a81e61SBarry Smith   Input Parameters:
41101a81e61SBarry Smith + pc - the preconditioner
41201a81e61SBarry Smith - n - 0 or 1
41301a81e61SBarry Smith 
414f1580f4eSBarry Smith   Note:
41595452b02SPatrick Sanan     If A has a symmetric nonzero pattern use -sp 1 to
41601a81e61SBarry Smith                  improve performance by eliminating some communication
41701a81e61SBarry Smith                  in the parallel version. Even if A does not have a
41801a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
41901a81e61SBarry Smith                  results, but the code will not follow the published
42001a81e61SBarry Smith                  SPAI algorithm exactly.
42101a81e61SBarry Smith 
42201a81e61SBarry Smith   Level: intermediate
42301a81e61SBarry Smith 
424db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
42501a81e61SBarry Smith @*/
426d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp)
427d71ae5a4SJacob Faibussowitsch {
42801a81e61SBarry Smith   PetscFunctionBegin;
429b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp));
43001a81e61SBarry Smith   PetscFunctionReturn(0);
43101a81e61SBarry Smith }
43201a81e61SBarry Smith 
433d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_SPAI(PC pc, PetscOptionItems *PetscOptionsObject)
434d71ae5a4SJacob Faibussowitsch {
43501a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
43601a81e61SBarry Smith   int       nbsteps1, max1, maxnew1, block_size1, cache_size, verbose, sp;
43701a81e61SBarry Smith   double    epsilon1;
438ace3abfcSBarry Smith   PetscBool flg;
43901a81e61SBarry Smith 
44001a81e61SBarry Smith   PetscFunctionBegin;
441d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SPAI options");
4429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_spai_epsilon", "", "PCSPAISetEpsilon", ispai->epsilon, &epsilon1, &flg));
4431baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetEpsilon(pc, epsilon1));
4449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_nbsteps", "", "PCSPAISetNBSteps", ispai->nbsteps, &nbsteps1, &flg));
4451baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetNBSteps(pc, nbsteps1));
44601a81e61SBarry Smith   /* added 1/7/99 g.h. */
4479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_max", "", "PCSPAISetMax", ispai->max, &max1, &flg));
4481baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMax(pc, max1));
4499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_maxnew", "", "PCSPAISetMaxNew", ispai->maxnew, &maxnew1, &flg));
4501baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMaxNew(pc, maxnew1));
4519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_block_size", "", "PCSPAISetBlockSize", ispai->block_size, &block_size1, &flg));
4521baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetBlockSize(pc, block_size1));
4539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_cache_size", "", "PCSPAISetCacheSize", ispai->cache_size, &cache_size, &flg));
4541baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetCacheSize(pc, cache_size));
4559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_verbose", "", "PCSPAISetVerbose", ispai->verbose, &verbose, &flg));
4561baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetVerbose(pc, verbose));
4579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_sp", "", "PCSPAISetSp", ispai->sp, &sp, &flg));
4581baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetSp(pc, sp));
459d0609cedSBarry Smith   PetscOptionsHeadEnd();
46001a81e61SBarry Smith   PetscFunctionReturn(0);
46101a81e61SBarry Smith }
46201a81e61SBarry Smith 
46301a81e61SBarry Smith /*MC
464f1580f4eSBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method
46501a81e61SBarry Smith 
46601a81e61SBarry Smith    Options Database Keys:
467c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
46801a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
46901a81e61SBarry Smith .  -pc_spai_max <m> - set max
47001a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
47101a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
47201a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
47301a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
47401a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
47501a81e61SBarry Smith 
47601a81e61SBarry Smith    Level: beginner
47701a81e61SBarry Smith 
478f1580f4eSBarry Smith    Note:
479f1580f4eSBarry Smith     This only works with `MATAIJ` matrices.
480f1580f4eSBarry Smith 
481f1580f4eSBarry Smith    References:
482f1580f4eSBarry Smith  . * -  Grote and Barnard (SIAM J. Sci. Comput.; vol 18, nr 3)
483f1580f4eSBarry Smith 
484db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
485db781477SPatrick Sanan           `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`,
486f1580f4eSBarry Smith           `PCSPAISetVerbose()`, `PCSPAISetSp()`, `PCSPAISetNBSteps()`, `PCSPAISetCacheSize()`
48701a81e61SBarry Smith M*/
48801a81e61SBarry Smith 
489d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc)
490d71ae5a4SJacob Faibussowitsch {
49101a81e61SBarry Smith   PC_SPAI *ispai;
49201a81e61SBarry Smith 
49301a81e61SBarry Smith   PetscFunctionBegin;
4944dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ispai));
4952a4967abSBarry Smith   pc->data = ispai;
49601a81e61SBarry Smith 
49701a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
49801a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
49948e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_SPAI;
50001a81e61SBarry Smith   pc->ops->applyrichardson = 0;
50101a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
50201a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
50301a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
50401a81e61SBarry Smith 
50501a81e61SBarry Smith   ispai->epsilon    = .4;
50601a81e61SBarry Smith   ispai->nbsteps    = 5;
50701a81e61SBarry Smith   ispai->max        = 5000;
50801a81e61SBarry Smith   ispai->maxnew     = 5;
50901a81e61SBarry Smith   ispai->block_size = 1;
51001a81e61SBarry Smith   ispai->cache_size = 5;
51101a81e61SBarry Smith   ispai->verbose    = 0;
51201a81e61SBarry Smith 
51301a81e61SBarry Smith   ispai->sp = 1;
5149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &(ispai->comm_spai)));
51501a81e61SBarry Smith 
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", PCSPAISetEpsilon_SPAI));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", PCSPAISetNBSteps_SPAI));
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", PCSPAISetMax_SPAI));
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", PCSPAISetMaxNew_SPAI));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", PCSPAISetBlockSize_SPAI));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", PCSPAISetCacheSize_SPAI));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", PCSPAISetVerbose_SPAI));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", PCSPAISetSp_SPAI));
52401a81e61SBarry Smith   PetscFunctionReturn(0);
52501a81e61SBarry Smith }
52601a81e61SBarry Smith 
52701a81e61SBarry Smith /*
52801a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
52901a81e61SBarry Smith */
530d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A, Mat AT, matrix **B)
531d71ae5a4SJacob Faibussowitsch {
53201a81e61SBarry Smith   matrix                  *M;
53301a81e61SBarry Smith   int                      i, j, col;
53401a81e61SBarry Smith   int                      row_indx;
53501a81e61SBarry Smith   int                      len, pe, local_indx, start_indx;
53601a81e61SBarry Smith   int                     *mapping;
53701a81e61SBarry Smith   const int               *cols;
53801a81e61SBarry Smith   const double            *vals;
5392122902bSSatish Balay   int                      n, mnl, nnl, nz, rstart, rend;
5402a4967abSBarry Smith   PetscMPIInt              size, rank;
54101a81e61SBarry Smith   struct compressed_lines *rows;
54201a81e61SBarry Smith 
54301a81e61SBarry Smith   PetscFunctionBegin;
5449566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5469566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &n, &n));
5479566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &mnl, &nnl));
54801a81e61SBarry Smith 
54901a81e61SBarry Smith   /*
55001a81e61SBarry Smith     not sure why a barrier is required. commenting out
5519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(comm));
55201a81e61SBarry Smith   */
55301a81e61SBarry Smith 
55429b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
55501a81e61SBarry Smith 
55601a81e61SBarry Smith   M->n              = n;
55701a81e61SBarry Smith   M->bs             = 1;
55801a81e61SBarry Smith   M->max_block_size = 1;
55901a81e61SBarry Smith 
56001a81e61SBarry Smith   M->mnls          = (int *)malloc(sizeof(int) * size);
56101a81e61SBarry Smith   M->start_indices = (int *)malloc(sizeof(int) * size);
56201a81e61SBarry Smith   M->pe            = (int *)malloc(sizeof(int) * n);
56301a81e61SBarry Smith   M->block_sizes   = (int *)malloc(sizeof(int) * n);
56401a81e61SBarry Smith   for (i = 0; i < n; i++) M->block_sizes[i] = 1;
56501a81e61SBarry Smith 
5669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnl, 1, MPI_INT, M->mnls, 1, MPI_INT, comm));
56701a81e61SBarry Smith 
56801a81e61SBarry Smith   M->start_indices[0] = 0;
5692fa5cd67SKarl Rupp   for (i = 1; i < size; i++) M->start_indices[i] = M->start_indices[i - 1] + M->mnls[i - 1];
57001a81e61SBarry Smith 
57101a81e61SBarry Smith   M->mnl            = M->mnls[M->myid];
57201a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
57301a81e61SBarry Smith 
57401a81e61SBarry Smith   for (i = 0; i < size; i++) {
57501a81e61SBarry Smith     start_indx = M->start_indices[i];
5762fa5cd67SKarl Rupp     for (j = 0; j < M->mnls[i]; j++) M->pe[start_indx + j] = i;
57701a81e61SBarry Smith   }
57801a81e61SBarry Smith 
57901a81e61SBarry Smith   if (AT) {
5802f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 1);
58101a81e61SBarry Smith   } else {
5822f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 0);
58301a81e61SBarry Smith   }
58401a81e61SBarry Smith 
58501a81e61SBarry Smith   rows = M->lines;
58601a81e61SBarry Smith 
58701a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
5889566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M->n, &mapping));
58901a81e61SBarry Smith   pe         = 0;
59001a81e61SBarry Smith   local_indx = 0;
59101a81e61SBarry Smith   for (i = 0; i < M->n; i++) {
59201a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
59301a81e61SBarry Smith       pe++;
59401a81e61SBarry Smith       local_indx = 0;
59501a81e61SBarry Smith     }
59601a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
59701a81e61SBarry Smith     local_indx++;
59801a81e61SBarry Smith   }
59901a81e61SBarry Smith 
60001a81e61SBarry Smith   /************** Set up the row structure *****************/
60101a81e61SBarry Smith 
6029566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
60301a81e61SBarry Smith   for (i = rstart; i < rend; i++) {
60401a81e61SBarry Smith     row_indx = i - rstart;
6059566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nz, &cols, &vals));
6062122902bSSatish Balay     /* allocate buffers */
6072122902bSSatish Balay     rows->ptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6082122902bSSatish Balay     rows->A[row_indx]    = (double *)malloc(nz * sizeof(double));
6092122902bSSatish Balay     /* copy the matrix */
61001a81e61SBarry Smith     for (j = 0; j < nz; j++) {
61101a81e61SBarry Smith       col = cols[j];
61201a81e61SBarry Smith       len = rows->len[row_indx]++;
6132fa5cd67SKarl Rupp 
61401a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
61501a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
61601a81e61SBarry Smith     }
61701a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
6182fa5cd67SKarl Rupp 
6199566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nz, &cols, &vals));
62001a81e61SBarry Smith   }
62101a81e61SBarry Smith 
62201a81e61SBarry Smith   /************** Set up the column structure *****************/
62301a81e61SBarry Smith 
62401a81e61SBarry Smith   if (AT) {
62501a81e61SBarry Smith     for (i = rstart; i < rend; i++) {
62601a81e61SBarry Smith       row_indx = i - rstart;
6279566063dSJacob Faibussowitsch       PetscCall(MatGetRow(AT, i, &nz, &cols, &vals));
6282122902bSSatish Balay       /* allocate buffers */
6292122902bSSatish Balay       rows->rptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6302122902bSSatish Balay       /* copy the matrix (i.e., the structure) */
63101a81e61SBarry Smith       for (j = 0; j < nz; j++) {
63201a81e61SBarry Smith         col = cols[j];
63301a81e61SBarry Smith         len = rows->rlen[row_indx]++;
6342fa5cd67SKarl Rupp 
63501a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
63601a81e61SBarry Smith       }
6379566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(AT, i, &nz, &cols, &vals));
63801a81e61SBarry Smith     }
63901a81e61SBarry Smith   }
64001a81e61SBarry Smith 
6419566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping));
64201a81e61SBarry Smith 
64301a81e61SBarry Smith   order_pointers(M);
64401a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
64501a81e61SBarry Smith   *B       = M;
64601a81e61SBarry Smith   PetscFunctionReturn(0);
64701a81e61SBarry Smith }
64801a81e61SBarry Smith 
64901a81e61SBarry Smith /*
65001a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
651df3898eeSBarry Smith    This assumes that the SPAI matrix B is stored in
65201a81e61SBarry Smith    COMPRESSED-ROW format.
65301a81e61SBarry Smith */
654d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertMatrixToMat(MPI_Comm comm, matrix *B, Mat *PB)
655d71ae5a4SJacob Faibussowitsch {
6564b6b5fe1SBarry Smith   PetscMPIInt size, rank;
65701a81e61SBarry Smith   int         m, n, M, N;
65801a81e61SBarry Smith   int         d_nz, o_nz;
65901a81e61SBarry Smith   int        *d_nnz, *o_nnz;
66001a81e61SBarry Smith   int         i, k, global_row, global_col, first_diag_col, last_diag_col;
66101a81e61SBarry Smith   PetscScalar val;
66201a81e61SBarry Smith 
66301a81e61SBarry Smith   PetscFunctionBegin;
6649566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
66601a81e61SBarry Smith 
66701a81e61SBarry Smith   m = n = B->mnls[rank];
66801a81e61SBarry Smith   d_nz = o_nz = 0;
66901a81e61SBarry Smith 
67006946f3aSJose E. Roman   /* Determine preallocation for MatCreateAIJ */
6719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &d_nnz));
6729566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &o_nnz));
67301a81e61SBarry Smith   for (i = 0; i < m; i++) d_nnz[i] = o_nnz[i] = 0;
67401a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
67501a81e61SBarry Smith   last_diag_col  = first_diag_col + B->mnls[rank];
67601a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
67701a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
67801a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
679db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
680db4deed7SKarl Rupp       else o_nnz[i]++;
68101a81e61SBarry Smith     }
68201a81e61SBarry Smith   }
68301a81e61SBarry Smith 
68401a81e61SBarry Smith   M = N = B->n;
68501a81e61SBarry Smith   /* Here we only know how to create AIJ format */
6869566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, PB));
6879566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*PB, m, n, M, N));
6889566063dSJacob Faibussowitsch   PetscCall(MatSetType(*PB, MATAIJ));
6899566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*PB, d_nz, d_nnz));
6909566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*PB, d_nz, d_nnz, o_nz, o_nnz));
69101a81e61SBarry Smith 
69201a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
69301a81e61SBarry Smith     global_row = B->start_indices[rank] + i;
69401a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
69501a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
6962fa5cd67SKarl Rupp 
69701a81e61SBarry Smith       val = B->lines->A[i][k];
6989566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*PB, 1, &global_row, 1, &global_col, &val, ADD_VALUES));
69901a81e61SBarry Smith     }
70001a81e61SBarry Smith   }
70101a81e61SBarry Smith 
7029566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
7039566063dSJacob Faibussowitsch   PetscCall(PetscFree(o_nnz));
70401a81e61SBarry Smith 
7059566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*PB, MAT_FINAL_ASSEMBLY));
7069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*PB, MAT_FINAL_ASSEMBLY));
70701a81e61SBarry Smith   PetscFunctionReturn(0);
70801a81e61SBarry Smith }
70901a81e61SBarry Smith 
71001a81e61SBarry Smith /*
71101a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
71201a81e61SBarry Smith */
713d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertVectorToVec(MPI_Comm comm, vector *v, Vec *Pv)
714d71ae5a4SJacob Faibussowitsch {
7154b6b5fe1SBarry Smith   PetscMPIInt size, rank;
7164b6b5fe1SBarry Smith   int         m, M, i, *mnls, *start_indices, *global_indices;
71701a81e61SBarry Smith 
71801a81e61SBarry Smith   PetscFunctionBegin;
7199566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
72101a81e61SBarry Smith 
72201a81e61SBarry Smith   m = v->mnl;
72301a81e61SBarry Smith   M = v->n;
72401a81e61SBarry Smith 
7259566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, m, M, Pv));
72601a81e61SBarry Smith 
7279566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &mnls));
7289566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&v->mnl, 1, MPI_INT, mnls, 1, MPI_INT, comm));
72901a81e61SBarry Smith 
7309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &start_indices));
7312fa5cd67SKarl Rupp 
73201a81e61SBarry Smith   start_indices[0] = 0;
7332fa5cd67SKarl Rupp   for (i = 1; i < size; i++) start_indices[i] = start_indices[i - 1] + mnls[i - 1];
73401a81e61SBarry Smith 
7359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(v->mnl, &global_indices));
7362fa5cd67SKarl Rupp   for (i = 0; i < v->mnl; i++) global_indices[i] = start_indices[rank] + i;
73701a81e61SBarry Smith 
7389566063dSJacob Faibussowitsch   PetscCall(PetscFree(mnls));
7399566063dSJacob Faibussowitsch   PetscCall(PetscFree(start_indices));
74001a81e61SBarry Smith 
7419566063dSJacob Faibussowitsch   PetscCall(VecSetValues(*Pv, v->mnl, global_indices, v->v, INSERT_VALUES));
7429566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*Pv));
7439566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*Pv));
74401a81e61SBarry Smith 
7459566063dSJacob Faibussowitsch   PetscCall(PetscFree(global_indices));
74601a81e61SBarry Smith   PetscFunctionReturn(0);
74701a81e61SBarry Smith }
748