xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision feefa0e191a340680bb02e1467a36facdcb0b150)
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 */
2074c01ffaSSatish Balay #if !defined(PETSC_SKIP_COMPLEX)
21762360ffSBarry Smith   #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */
2274c01ffaSSatish Balay #endif
2301a81e61SBarry Smith 
24af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
251c7d2463SJed Brown #include <../src/ksp/pc/impls/spai/petscspai.h>
2601a81e61SBarry Smith 
2701a81e61SBarry Smith /*
2801a81e61SBarry Smith     These are the SPAI include files
2901a81e61SBarry Smith */
3001a81e61SBarry Smith EXTERN_C_BEGIN
31b1bda85cSSatish Balay #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */
32c6db04a5SJed Brown #include <spai.h>
33c6db04a5SJed Brown #include <matrix.h>
3401a81e61SBarry Smith EXTERN_C_END
3501a81e61SBarry Smith 
3609573ac7SBarry Smith extern PetscErrorCode ConvertMatToMatrix(MPI_Comm, Mat, Mat, matrix **);
3709573ac7SBarry Smith extern PetscErrorCode ConvertMatrixToMat(MPI_Comm, matrix *, Mat *);
38b03515a0SUmberto Zerbinati extern PetscErrorCode ConvertVectorToVec(MPI_Comm, vector *, Vec *);
3909573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char *, char *, char *);
4001a81e61SBarry Smith 
4101a81e61SBarry Smith typedef struct {
4201a81e61SBarry Smith   matrix *B;  /* matrix in SPAI format */
4301a81e61SBarry Smith   matrix *BT; /* transpose of matrix in SPAI format */
4401a81e61SBarry Smith   matrix *M;  /* the approximate inverse in SPAI format */
4501a81e61SBarry Smith 
4601a81e61SBarry Smith   Mat PM; /* the approximate inverse PETSc format */
4701a81e61SBarry Smith 
4801a81e61SBarry Smith   double epsilon;    /* tolerance */
4901a81e61SBarry Smith   int    nbsteps;    /* max number of "improvement" steps per line */
5001a81e61SBarry Smith   int    max;        /* max dimensions of is_I, q, etc. */
5101a81e61SBarry Smith   int    maxnew;     /* max number of new entries per step */
5201a81e61SBarry Smith   int    block_size; /* constant block size */
5301a81e61SBarry Smith   int    cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */
5401a81e61SBarry Smith   int    verbose;    /* SPAI prints timing and statistics */
5501a81e61SBarry Smith 
5601a81e61SBarry Smith   int      sp;        /* symmetric nonzero pattern */
5701a81e61SBarry Smith   MPI_Comm comm_spai; /* communicator to be used with spai */
5801a81e61SBarry Smith } PC_SPAI;
5901a81e61SBarry Smith 
60d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SPAI(PC pc)
61d71ae5a4SJacob Faibussowitsch {
6201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
6301a81e61SBarry Smith   Mat      AT;
6401a81e61SBarry Smith 
6501a81e61SBarry Smith   PetscFunctionBegin;
6601a81e61SBarry Smith   init_SPAI();
6701a81e61SBarry Smith 
6801a81e61SBarry Smith   if (ispai->sp) {
699566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, pc->pmat, &ispai->B));
7001a81e61SBarry Smith   } else {
7101a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
729566063dSJacob Faibussowitsch     PetscCall(MatTranspose(pc->pmat, MAT_INITIAL_MATRIX, &AT));
739566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, AT, &ispai->B));
749566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
7501a81e61SBarry Smith   }
7601a81e61SBarry Smith 
7701a81e61SBarry Smith   /* Destroy the transpose */
7801a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
7901a81e61SBarry Smith 
8001a81e61SBarry Smith   /* construct SPAI preconditioner */
8101a81e61SBarry Smith   /* FILE *messages */    /* file for warning messages */
8201a81e61SBarry Smith   /* double epsilon */    /* tolerance */
8301a81e61SBarry Smith   /* int nbsteps */       /* max number of "improvement" steps per line */
8401a81e61SBarry Smith   /* int max */           /* max dimensions of is_I, q, etc. */
8501a81e61SBarry Smith   /* int maxnew */        /* max number of new entries per step */
86d5b43468SJose E. Roman   /* int block_size */    /* block_size == 1 specifies scalar elements
8701a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
8801a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
892fa5cd67SKarl Rupp   /* int cache_size */    /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */
9001a81e61SBarry Smith   /* int    verbose    */ /* verbose == 0 specifies that SPAI is silent
9101a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9201a81e61SBarry Smith 
933ba16761SJacob Faibussowitsch   PetscCallExternal(bspai, ispai->B, &ispai->M, stdout, ispai->epsilon, ispai->nbsteps, ispai->max, ispai->maxnew, ispai->block_size, ispai->cache_size, ispai->verbose);
9401a81e61SBarry Smith 
959566063dSJacob Faibussowitsch   PetscCall(ConvertMatrixToMat(PetscObjectComm((PetscObject)pc), ispai->M, &ispai->PM));
9601a81e61SBarry Smith 
9701a81e61SBarry Smith   /* free the SPAI matrices */
9801a81e61SBarry Smith   sp_free_matrix(ispai->B);
9901a81e61SBarry Smith   sp_free_matrix(ispai->M);
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10101a81e61SBarry Smith }
10201a81e61SBarry Smith 
103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SPAI(PC pc, Vec xx, Vec y)
104d71ae5a4SJacob Faibussowitsch {
10501a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
10601a81e61SBarry Smith 
10701a81e61SBarry Smith   PetscFunctionBegin;
10801a81e61SBarry Smith   /* Now using PETSc's multiply */
1099566063dSJacob Faibussowitsch   PetscCall(MatMult(ispai->PM, xx, y));
1103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11101a81e61SBarry Smith }
11201a81e61SBarry Smith 
113d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SPAI(PC pc, Mat X, Mat Y)
114d71ae5a4SJacob Faibussowitsch {
11548e38a8aSPierre Jolivet   PC_SPAI *ispai = (PC_SPAI *)pc->data;
11648e38a8aSPierre Jolivet 
11748e38a8aSPierre Jolivet   PetscFunctionBegin;
11848e38a8aSPierre Jolivet   /* Now using PETSc's multiply */
1199566063dSJacob Faibussowitsch   PetscCall(MatMatMult(ispai->PM, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12148e38a8aSPierre Jolivet }
12248e38a8aSPierre Jolivet 
123d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SPAI(PC pc)
124d71ae5a4SJacob Faibussowitsch {
12501a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
12601a81e61SBarry Smith 
12701a81e61SBarry Smith   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ispai->PM));
1299566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&(ispai->comm_spai)));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", NULL));
1322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", NULL));
1332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", NULL));
1342e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", NULL));
1352e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", NULL));
1362e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", NULL));
1372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", NULL));
1382e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", NULL));
1393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14001a81e61SBarry Smith }
14101a81e61SBarry Smith 
142d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SPAI(PC pc, PetscViewer viewer)
143d71ae5a4SJacob Faibussowitsch {
14401a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
145ace3abfcSBarry Smith   PetscBool iascii;
14601a81e61SBarry Smith 
14701a81e61SBarry Smith   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
14901a81e61SBarry Smith   if (iascii) {
150b03515a0SUmberto Zerbinati     PetscCall(PetscViewerASCIIPrintf(viewer, "    epsilon %g\n", ispai->epsilon));
1519566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    nbsteps %d\n", ispai->nbsteps));
1529566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    max %d\n", ispai->max));
1539566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    maxnew %d\n", ispai->maxnew));
1549566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    block_size %d\n", ispai->block_size));
1559566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    cache_size %d\n", ispai->cache_size));
1569566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    verbose %d\n", ispai->verbose));
1579566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "    sp %d\n", ispai->sp));
15801a81e61SBarry Smith   }
1593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16001a81e61SBarry Smith }
16101a81e61SBarry Smith 
162d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1)
163d71ae5a4SJacob Faibussowitsch {
16401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1655fd66863SKarl Rupp 
16601a81e61SBarry Smith   PetscFunctionBegin;
167b03515a0SUmberto Zerbinati   ispai->epsilon = (double)epsilon1;
1683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16901a81e61SBarry Smith }
17001a81e61SBarry Smith 
171d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc, PetscInt nbsteps1)
172d71ae5a4SJacob Faibussowitsch {
17301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1745fd66863SKarl Rupp 
17501a81e61SBarry Smith   PetscFunctionBegin;
176b03515a0SUmberto Zerbinati   ispai->nbsteps = (int)nbsteps1;
1773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17801a81e61SBarry Smith }
17901a81e61SBarry Smith 
18001a81e61SBarry Smith /* added 1/7/99 g.h. */
181d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMax_SPAI(PC pc, PetscInt max1)
182d71ae5a4SJacob Faibussowitsch {
18301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1845fd66863SKarl Rupp 
18501a81e61SBarry Smith   PetscFunctionBegin;
186b03515a0SUmberto Zerbinati   ispai->max = (int)max1;
1873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18801a81e61SBarry Smith }
18901a81e61SBarry Smith 
190d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc, PetscInt maxnew1)
191d71ae5a4SJacob Faibussowitsch {
19201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
1935fd66863SKarl Rupp 
19401a81e61SBarry Smith   PetscFunctionBegin;
195b03515a0SUmberto Zerbinati   ispai->maxnew = (int)maxnew1;
1963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
19701a81e61SBarry Smith }
19801a81e61SBarry Smith 
199d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc, PetscInt block_size1)
200d71ae5a4SJacob Faibussowitsch {
20101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2025fd66863SKarl Rupp 
20301a81e61SBarry Smith   PetscFunctionBegin;
204b03515a0SUmberto Zerbinati   ispai->block_size = (int)block_size1;
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20601a81e61SBarry Smith }
20701a81e61SBarry Smith 
208d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc, PetscInt cache_size)
209d71ae5a4SJacob Faibussowitsch {
21001a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2115fd66863SKarl Rupp 
21201a81e61SBarry Smith   PetscFunctionBegin;
213b03515a0SUmberto Zerbinati   ispai->cache_size = (int)cache_size;
2143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21501a81e61SBarry Smith }
21601a81e61SBarry Smith 
217d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc, PetscInt verbose)
218d71ae5a4SJacob Faibussowitsch {
21901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2205fd66863SKarl Rupp 
22101a81e61SBarry Smith   PetscFunctionBegin;
222b03515a0SUmberto Zerbinati   ispai->verbose = (int)verbose;
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
22401a81e61SBarry Smith }
22501a81e61SBarry Smith 
226d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetSp_SPAI(PC pc, PetscInt sp)
227d71ae5a4SJacob Faibussowitsch {
22801a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI *)pc->data;
2295fd66863SKarl Rupp 
23001a81e61SBarry Smith   PetscFunctionBegin;
231b03515a0SUmberto Zerbinati   ispai->sp = (int)sp;
2323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23301a81e61SBarry Smith }
23401a81e61SBarry Smith 
23501a81e61SBarry Smith /*@
236f1580f4eSBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the `PCSPAI` preconditioner
23701a81e61SBarry Smith 
23801a81e61SBarry Smith   Input Parameters:
23901a81e61SBarry Smith + pc       - the preconditioner
240*feefa0e1SJacob Faibussowitsch - epsilon1 - epsilon (default .4)
24101a81e61SBarry Smith 
242f1580f4eSBarry Smith   Note:
24395452b02SPatrick Sanan   Espilon must be between 0 and 1. It controls the
24401a81e61SBarry Smith   quality of the approximation of M to the inverse of
24501a81e61SBarry Smith   A. Higher values of epsilon lead to more work, more
24601a81e61SBarry Smith   fill, and usually better preconditioners. In many
24701a81e61SBarry Smith   cases the best choice of epsilon is the one that
24801a81e61SBarry Smith   divides the total solution time equally between the
24901a81e61SBarry Smith   preconditioner and the solver.
25001a81e61SBarry Smith 
25101a81e61SBarry Smith   Level: intermediate
25201a81e61SBarry Smith 
253db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
25401a81e61SBarry Smith   @*/
255d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1)
256d71ae5a4SJacob Faibussowitsch {
25701a81e61SBarry Smith   PetscFunctionBegin;
258b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1));
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26001a81e61SBarry Smith }
26101a81e61SBarry Smith 
26201a81e61SBarry Smith /*@
26301a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
264f1580f4eSBarry Smith   the `PCSPAI` preconditioner
26501a81e61SBarry Smith 
26601a81e61SBarry Smith   Input Parameters:
26701a81e61SBarry Smith + pc       - the preconditioner
268*feefa0e1SJacob Faibussowitsch - nbsteps1 - number of steps (default 5)
26901a81e61SBarry Smith 
270f1580f4eSBarry Smith   Note:
271f1580f4eSBarry Smith   `PCSPAI` constructs to approximation to every column of
27201a81e61SBarry Smith   the exact inverse of A in a series of improvement
27301a81e61SBarry Smith   steps. The quality of the approximation is determined
27401a81e61SBarry Smith   by epsilon. If an approximation achieving an accuracy
27501a81e61SBarry Smith   of epsilon is not obtained after ns steps, SPAI simply
27601a81e61SBarry Smith   uses the best approximation constructed so far.
27701a81e61SBarry Smith 
27801a81e61SBarry Smith   Level: intermediate
27901a81e61SBarry Smith 
280db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
28101a81e61SBarry Smith @*/
282d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1)
283d71ae5a4SJacob Faibussowitsch {
28401a81e61SBarry Smith   PetscFunctionBegin;
285b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1));
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
28701a81e61SBarry Smith }
28801a81e61SBarry Smith 
28901a81e61SBarry Smith /* added 1/7/99 g.h. */
29001a81e61SBarry Smith /*@
29101a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
292f1580f4eSBarry Smith   the `PCSPAI` preconditioner
29301a81e61SBarry Smith 
29401a81e61SBarry Smith   Input Parameters:
29501a81e61SBarry Smith + pc   - the preconditioner
296*feefa0e1SJacob Faibussowitsch - max1 - size (default is 5000)
29701a81e61SBarry Smith 
29801a81e61SBarry Smith   Level: intermediate
29901a81e61SBarry Smith 
300db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
30101a81e61SBarry Smith @*/
302d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1)
303d71ae5a4SJacob Faibussowitsch {
30401a81e61SBarry Smith   PetscFunctionBegin;
305b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1));
3063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30701a81e61SBarry Smith }
30801a81e61SBarry Smith 
30901a81e61SBarry Smith /*@
31001a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
311f1580f4eSBarry Smith   in `PCSPAI` preconditioner
31201a81e61SBarry Smith 
31301a81e61SBarry Smith   Input Parameters:
31401a81e61SBarry Smith + pc      - the preconditioner
315*feefa0e1SJacob Faibussowitsch - maxnew1 - maximum number (default 5)
31601a81e61SBarry Smith 
31701a81e61SBarry Smith   Level: intermediate
31801a81e61SBarry Smith 
319db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
32001a81e61SBarry Smith @*/
321d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1)
322d71ae5a4SJacob Faibussowitsch {
32301a81e61SBarry Smith   PetscFunctionBegin;
324b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1));
3253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32601a81e61SBarry Smith }
32701a81e61SBarry Smith 
32801a81e61SBarry Smith /*@
329f1580f4eSBarry Smith   PCSPAISetBlockSize - set the block size for the `PCSPAI` preconditioner
33001a81e61SBarry Smith 
33101a81e61SBarry Smith   Input Parameters:
33201a81e61SBarry Smith + pc          - the preconditioner
333*feefa0e1SJacob Faibussowitsch - block_size1 - block size (default 1)
33401a81e61SBarry Smith 
33595452b02SPatrick Sanan   Notes:
33695452b02SPatrick Sanan   A block
33701a81e61SBarry Smith   size of 1 treats A as a matrix of scalar elements. A
33801a81e61SBarry Smith   block size of s > 1 treats A as a matrix of sxs
33901a81e61SBarry Smith   blocks. A block size of 0 treats A as a matrix with
34001a81e61SBarry Smith   variable sized blocks, which are determined by
34101a81e61SBarry Smith   searching for dense square diagonal blocks in A.
34201a81e61SBarry Smith   This can be very effective for finite-element
34301a81e61SBarry Smith   matrices.
34401a81e61SBarry Smith 
34501a81e61SBarry Smith   SPAI will convert A to block form, use a block
34601a81e61SBarry Smith   version of the preconditioner algorithm, and then
34701a81e61SBarry Smith   convert the result back to scalar form.
34801a81e61SBarry Smith 
34901a81e61SBarry Smith   In many cases the a block-size parameter other than 1
35001a81e61SBarry Smith   can lead to very significant improvement in
35101a81e61SBarry Smith   performance.
35201a81e61SBarry Smith 
35301a81e61SBarry Smith   Level: intermediate
35401a81e61SBarry Smith 
355db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
35601a81e61SBarry Smith @*/
357d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetBlockSize(PC pc, PetscInt block_size1)
358d71ae5a4SJacob Faibussowitsch {
35901a81e61SBarry Smith   PetscFunctionBegin;
360b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetBlockSize_C", (PC, PetscInt), (pc, block_size1));
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
36201a81e61SBarry Smith }
36301a81e61SBarry Smith 
36401a81e61SBarry Smith /*@
365f1580f4eSBarry Smith   PCSPAISetCacheSize - specify cache size in the `PCSPAI` preconditioner
36601a81e61SBarry Smith 
36701a81e61SBarry Smith   Input Parameters:
36801a81e61SBarry Smith + pc         - the preconditioner
369*feefa0e1SJacob Faibussowitsch - cache_size - cache size {0,1,2,3,4,5} (default 5)
37001a81e61SBarry Smith 
371f1580f4eSBarry Smith   Note:
372f1580f4eSBarry Smith   `PCSPAI` uses a hash table to cache messages and avoid
37301a81e61SBarry Smith   redundant communication. If suggest always using
37401a81e61SBarry Smith   5. This parameter is irrelevant in the serial
37501a81e61SBarry Smith   version.
37601a81e61SBarry Smith 
37701a81e61SBarry Smith   Level: intermediate
37801a81e61SBarry Smith 
379db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
38001a81e61SBarry Smith @*/
381d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetCacheSize(PC pc, PetscInt cache_size)
382d71ae5a4SJacob Faibussowitsch {
38301a81e61SBarry Smith   PetscFunctionBegin;
384b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetCacheSize_C", (PC, PetscInt), (pc, cache_size));
3853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38601a81e61SBarry Smith }
38701a81e61SBarry Smith 
38801a81e61SBarry Smith /*@
389f1580f4eSBarry Smith   PCSPAISetVerbose - verbosity level for the `PCSPAI` preconditioner
39001a81e61SBarry Smith 
39101a81e61SBarry Smith   Input Parameters:
39201a81e61SBarry Smith + pc      - the preconditioner
393*feefa0e1SJacob Faibussowitsch - verbose - level (default 1)
39401a81e61SBarry Smith 
395f1580f4eSBarry Smith   Note:
39695452b02SPatrick Sanan   print parameters, timings and matrix statistics
39701a81e61SBarry Smith 
39801a81e61SBarry Smith   Level: intermediate
39901a81e61SBarry Smith 
400db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
40101a81e61SBarry Smith @*/
402d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose)
403d71ae5a4SJacob Faibussowitsch {
40401a81e61SBarry Smith   PetscFunctionBegin;
405b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose));
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40701a81e61SBarry Smith }
40801a81e61SBarry Smith 
40901a81e61SBarry Smith /*@
410f1580f4eSBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the `PCSPAI` preconditioner
41101a81e61SBarry Smith 
41201a81e61SBarry Smith   Input Parameters:
41301a81e61SBarry Smith + pc - the preconditioner
414*feefa0e1SJacob Faibussowitsch - sp - 0 or 1
41501a81e61SBarry Smith 
416f1580f4eSBarry Smith   Note:
41795452b02SPatrick Sanan   If A has a symmetric nonzero pattern use -sp 1 to
41801a81e61SBarry Smith   improve performance by eliminating some communication
41901a81e61SBarry Smith   in the parallel version. Even if A does not have a
42001a81e61SBarry Smith   symmetric nonzero pattern -sp 1 may well lead to good
42101a81e61SBarry Smith   results, but the code will not follow the published
42201a81e61SBarry Smith   SPAI algorithm exactly.
42301a81e61SBarry Smith 
42401a81e61SBarry Smith   Level: intermediate
42501a81e61SBarry Smith 
426db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
42701a81e61SBarry Smith @*/
428d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp)
429d71ae5a4SJacob Faibussowitsch {
43001a81e61SBarry Smith   PetscFunctionBegin;
431b03515a0SUmberto Zerbinati   PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp));
4323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
43301a81e61SBarry Smith }
43401a81e61SBarry Smith 
435d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_SPAI(PC pc, PetscOptionItems *PetscOptionsObject)
436d71ae5a4SJacob Faibussowitsch {
43701a81e61SBarry Smith   PC_SPAI  *ispai = (PC_SPAI *)pc->data;
43801a81e61SBarry Smith   int       nbsteps1, max1, maxnew1, block_size1, cache_size, verbose, sp;
43901a81e61SBarry Smith   double    epsilon1;
440ace3abfcSBarry Smith   PetscBool flg;
44101a81e61SBarry Smith 
44201a81e61SBarry Smith   PetscFunctionBegin;
443d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "SPAI options");
4449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_spai_epsilon", "", "PCSPAISetEpsilon", ispai->epsilon, &epsilon1, &flg));
4451baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetEpsilon(pc, epsilon1));
4469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_nbsteps", "", "PCSPAISetNBSteps", ispai->nbsteps, &nbsteps1, &flg));
4471baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetNBSteps(pc, nbsteps1));
44801a81e61SBarry Smith   /* added 1/7/99 g.h. */
4499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_max", "", "PCSPAISetMax", ispai->max, &max1, &flg));
4501baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMax(pc, max1));
4519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_maxnew", "", "PCSPAISetMaxNew", ispai->maxnew, &maxnew1, &flg));
4521baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetMaxNew(pc, maxnew1));
4539566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_block_size", "", "PCSPAISetBlockSize", ispai->block_size, &block_size1, &flg));
4541baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetBlockSize(pc, block_size1));
4559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_cache_size", "", "PCSPAISetCacheSize", ispai->cache_size, &cache_size, &flg));
4561baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetCacheSize(pc, cache_size));
4579566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_verbose", "", "PCSPAISetVerbose", ispai->verbose, &verbose, &flg));
4581baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetVerbose(pc, verbose));
4599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_sp", "", "PCSPAISetSp", ispai->sp, &sp, &flg));
4601baa6e33SBarry Smith   if (flg) PetscCall(PCSPAISetSp(pc, sp));
461d0609cedSBarry Smith   PetscOptionsHeadEnd();
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
46301a81e61SBarry Smith }
46401a81e61SBarry Smith 
46501a81e61SBarry Smith /*MC
466f1580f4eSBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method
46701a81e61SBarry Smith 
46801a81e61SBarry Smith    Options Database Keys:
469c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
47001a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
47101a81e61SBarry Smith .  -pc_spai_max <m> - set max
47201a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
47301a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
47401a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
47501a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
47601a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
47701a81e61SBarry Smith 
47801a81e61SBarry Smith    Level: beginner
47901a81e61SBarry Smith 
480f1580f4eSBarry Smith    Note:
481f1580f4eSBarry Smith     This only works with `MATAIJ` matrices.
482f1580f4eSBarry Smith 
483f1580f4eSBarry Smith    References:
484f1580f4eSBarry Smith  . * -  Grote and Barnard (SIAM J. Sci. Comput.; vol 18, nr 3)
485f1580f4eSBarry Smith 
486db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
487db781477SPatrick Sanan           `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`,
488f1580f4eSBarry Smith           `PCSPAISetVerbose()`, `PCSPAISetSp()`, `PCSPAISetNBSteps()`, `PCSPAISetCacheSize()`
48901a81e61SBarry Smith M*/
49001a81e61SBarry Smith 
491d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc)
492d71ae5a4SJacob Faibussowitsch {
49301a81e61SBarry Smith   PC_SPAI *ispai;
49401a81e61SBarry Smith 
49501a81e61SBarry Smith   PetscFunctionBegin;
4964dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ispai));
4972a4967abSBarry Smith   pc->data = ispai;
49801a81e61SBarry Smith 
49901a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
50001a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
50148e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_SPAI;
50201a81e61SBarry Smith   pc->ops->applyrichardson = 0;
50301a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
50401a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
50501a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
50601a81e61SBarry Smith 
50701a81e61SBarry Smith   ispai->epsilon    = .4;
50801a81e61SBarry Smith   ispai->nbsteps    = 5;
50901a81e61SBarry Smith   ispai->max        = 5000;
51001a81e61SBarry Smith   ispai->maxnew     = 5;
51101a81e61SBarry Smith   ispai->block_size = 1;
51201a81e61SBarry Smith   ispai->cache_size = 5;
51301a81e61SBarry Smith   ispai->verbose    = 0;
51401a81e61SBarry Smith 
51501a81e61SBarry Smith   ispai->sp = 1;
5169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &(ispai->comm_spai)));
51701a81e61SBarry Smith 
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", PCSPAISetEpsilon_SPAI));
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", PCSPAISetNBSteps_SPAI));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", PCSPAISetMax_SPAI));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", PCSPAISetMaxNew_SPAI));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", PCSPAISetBlockSize_SPAI));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", PCSPAISetCacheSize_SPAI));
5249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", PCSPAISetVerbose_SPAI));
5259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", PCSPAISetSp_SPAI));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
52701a81e61SBarry Smith }
52801a81e61SBarry Smith 
52901a81e61SBarry Smith /*
53001a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
53101a81e61SBarry Smith */
532d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A, Mat AT, matrix **B)
533d71ae5a4SJacob Faibussowitsch {
53401a81e61SBarry Smith   matrix                  *M;
53501a81e61SBarry Smith   int                      i, j, col;
53601a81e61SBarry Smith   int                      row_indx;
53701a81e61SBarry Smith   int                      len, pe, local_indx, start_indx;
53801a81e61SBarry Smith   int                     *mapping;
53901a81e61SBarry Smith   const int               *cols;
54001a81e61SBarry Smith   const double            *vals;
5412122902bSSatish Balay   int                      n, mnl, nnl, nz, rstart, rend;
5422a4967abSBarry Smith   PetscMPIInt              size, rank;
54301a81e61SBarry Smith   struct compressed_lines *rows;
54401a81e61SBarry Smith 
54501a81e61SBarry Smith   PetscFunctionBegin;
5469566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
5479566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5489566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A, &n, &n));
5499566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &mnl, &nnl));
55001a81e61SBarry Smith 
55101a81e61SBarry Smith   /*
55201a81e61SBarry Smith     not sure why a barrier is required. commenting out
5539566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(comm));
55401a81e61SBarry Smith   */
55501a81e61SBarry Smith 
55629b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
55701a81e61SBarry Smith 
55801a81e61SBarry Smith   M->n              = n;
55901a81e61SBarry Smith   M->bs             = 1;
56001a81e61SBarry Smith   M->max_block_size = 1;
56101a81e61SBarry Smith 
56201a81e61SBarry Smith   M->mnls          = (int *)malloc(sizeof(int) * size);
56301a81e61SBarry Smith   M->start_indices = (int *)malloc(sizeof(int) * size);
56401a81e61SBarry Smith   M->pe            = (int *)malloc(sizeof(int) * n);
56501a81e61SBarry Smith   M->block_sizes   = (int *)malloc(sizeof(int) * n);
56601a81e61SBarry Smith   for (i = 0; i < n; i++) M->block_sizes[i] = 1;
56701a81e61SBarry Smith 
5689566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnl, 1, MPI_INT, M->mnls, 1, MPI_INT, comm));
56901a81e61SBarry Smith 
57001a81e61SBarry Smith   M->start_indices[0] = 0;
5712fa5cd67SKarl Rupp   for (i = 1; i < size; i++) M->start_indices[i] = M->start_indices[i - 1] + M->mnls[i - 1];
57201a81e61SBarry Smith 
57301a81e61SBarry Smith   M->mnl            = M->mnls[M->myid];
57401a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
57501a81e61SBarry Smith 
57601a81e61SBarry Smith   for (i = 0; i < size; i++) {
57701a81e61SBarry Smith     start_indx = M->start_indices[i];
5782fa5cd67SKarl Rupp     for (j = 0; j < M->mnls[i]; j++) M->pe[start_indx + j] = i;
57901a81e61SBarry Smith   }
58001a81e61SBarry Smith 
58101a81e61SBarry Smith   if (AT) {
5822f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 1);
58301a81e61SBarry Smith   } else {
5842f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank], 0);
58501a81e61SBarry Smith   }
58601a81e61SBarry Smith 
58701a81e61SBarry Smith   rows = M->lines;
58801a81e61SBarry Smith 
58901a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
5909566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M->n, &mapping));
59101a81e61SBarry Smith   pe         = 0;
59201a81e61SBarry Smith   local_indx = 0;
59301a81e61SBarry Smith   for (i = 0; i < M->n; i++) {
59401a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
59501a81e61SBarry Smith       pe++;
59601a81e61SBarry Smith       local_indx = 0;
59701a81e61SBarry Smith     }
59801a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
59901a81e61SBarry Smith     local_indx++;
60001a81e61SBarry Smith   }
60101a81e61SBarry Smith 
60201a81e61SBarry Smith   /************** Set up the row structure *****************/
60301a81e61SBarry Smith 
6049566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
60501a81e61SBarry Smith   for (i = rstart; i < rend; i++) {
60601a81e61SBarry Smith     row_indx = i - rstart;
6079566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A, i, &nz, &cols, &vals));
6082122902bSSatish Balay     /* allocate buffers */
6092122902bSSatish Balay     rows->ptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6102122902bSSatish Balay     rows->A[row_indx]    = (double *)malloc(nz * sizeof(double));
6112122902bSSatish Balay     /* copy the matrix */
61201a81e61SBarry Smith     for (j = 0; j < nz; j++) {
61301a81e61SBarry Smith       col = cols[j];
61401a81e61SBarry Smith       len = rows->len[row_indx]++;
6152fa5cd67SKarl Rupp 
61601a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
61701a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
61801a81e61SBarry Smith     }
61901a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
6202fa5cd67SKarl Rupp 
6219566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A, i, &nz, &cols, &vals));
62201a81e61SBarry Smith   }
62301a81e61SBarry Smith 
62401a81e61SBarry Smith   /************** Set up the column structure *****************/
62501a81e61SBarry Smith 
62601a81e61SBarry Smith   if (AT) {
62701a81e61SBarry Smith     for (i = rstart; i < rend; i++) {
62801a81e61SBarry Smith       row_indx = i - rstart;
6299566063dSJacob Faibussowitsch       PetscCall(MatGetRow(AT, i, &nz, &cols, &vals));
6302122902bSSatish Balay       /* allocate buffers */
6312122902bSSatish Balay       rows->rptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6322122902bSSatish Balay       /* copy the matrix (i.e., the structure) */
63301a81e61SBarry Smith       for (j = 0; j < nz; j++) {
63401a81e61SBarry Smith         col = cols[j];
63501a81e61SBarry Smith         len = rows->rlen[row_indx]++;
6362fa5cd67SKarl Rupp 
63701a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
63801a81e61SBarry Smith       }
6399566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(AT, i, &nz, &cols, &vals));
64001a81e61SBarry Smith     }
64101a81e61SBarry Smith   }
64201a81e61SBarry Smith 
6439566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping));
64401a81e61SBarry Smith 
64501a81e61SBarry Smith   order_pointers(M);
64601a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
64701a81e61SBarry Smith   *B       = M;
6483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64901a81e61SBarry Smith }
65001a81e61SBarry Smith 
65101a81e61SBarry Smith /*
65201a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
653df3898eeSBarry Smith    This assumes that the SPAI matrix B is stored in
65401a81e61SBarry Smith    COMPRESSED-ROW format.
65501a81e61SBarry Smith */
656d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertMatrixToMat(MPI_Comm comm, matrix *B, Mat *PB)
657d71ae5a4SJacob Faibussowitsch {
6584b6b5fe1SBarry Smith   PetscMPIInt size, rank;
65901a81e61SBarry Smith   int         m, n, M, N;
66001a81e61SBarry Smith   int         d_nz, o_nz;
66101a81e61SBarry Smith   int        *d_nnz, *o_nnz;
66201a81e61SBarry Smith   int         i, k, global_row, global_col, first_diag_col, last_diag_col;
66301a81e61SBarry Smith   PetscScalar val;
66401a81e61SBarry Smith 
66501a81e61SBarry Smith   PetscFunctionBegin;
6669566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6679566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
66801a81e61SBarry Smith 
66901a81e61SBarry Smith   m = n = B->mnls[rank];
67001a81e61SBarry Smith   d_nz = o_nz = 0;
67101a81e61SBarry Smith 
67206946f3aSJose E. Roman   /* Determine preallocation for MatCreateAIJ */
6739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &d_nnz));
6749566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m, &o_nnz));
67501a81e61SBarry Smith   for (i = 0; i < m; i++) d_nnz[i] = o_nnz[i] = 0;
67601a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
67701a81e61SBarry Smith   last_diag_col  = first_diag_col + B->mnls[rank];
67801a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
67901a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
68001a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
681db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
682db4deed7SKarl Rupp       else o_nnz[i]++;
68301a81e61SBarry Smith     }
68401a81e61SBarry Smith   }
68501a81e61SBarry Smith 
68601a81e61SBarry Smith   M = N = B->n;
68701a81e61SBarry Smith   /* Here we only know how to create AIJ format */
6889566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm, PB));
6899566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*PB, m, n, M, N));
6909566063dSJacob Faibussowitsch   PetscCall(MatSetType(*PB, MATAIJ));
6919566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*PB, d_nz, d_nnz));
6929566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*PB, d_nz, d_nnz, o_nz, o_nnz));
69301a81e61SBarry Smith 
69401a81e61SBarry Smith   for (i = 0; i < B->mnls[rank]; i++) {
69501a81e61SBarry Smith     global_row = B->start_indices[rank] + i;
69601a81e61SBarry Smith     for (k = 0; k < B->lines->len[i]; k++) {
69701a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
6982fa5cd67SKarl Rupp 
69901a81e61SBarry Smith       val = B->lines->A[i][k];
7009566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*PB, 1, &global_row, 1, &global_col, &val, ADD_VALUES));
70101a81e61SBarry Smith     }
70201a81e61SBarry Smith   }
70301a81e61SBarry Smith 
7049566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
7059566063dSJacob Faibussowitsch   PetscCall(PetscFree(o_nnz));
70601a81e61SBarry Smith 
7079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*PB, MAT_FINAL_ASSEMBLY));
7089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*PB, MAT_FINAL_ASSEMBLY));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
71001a81e61SBarry Smith }
71101a81e61SBarry Smith 
71201a81e61SBarry Smith /*
71301a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
71401a81e61SBarry Smith */
715d71ae5a4SJacob Faibussowitsch PetscErrorCode ConvertVectorToVec(MPI_Comm comm, vector *v, Vec *Pv)
716d71ae5a4SJacob Faibussowitsch {
7174b6b5fe1SBarry Smith   PetscMPIInt size, rank;
7184b6b5fe1SBarry Smith   int         m, M, i, *mnls, *start_indices, *global_indices;
71901a81e61SBarry Smith 
72001a81e61SBarry Smith   PetscFunctionBegin;
7219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
72301a81e61SBarry Smith 
72401a81e61SBarry Smith   m = v->mnl;
72501a81e61SBarry Smith   M = v->n;
72601a81e61SBarry Smith 
7279566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm, m, M, Pv));
72801a81e61SBarry Smith 
7299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &mnls));
7309566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&v->mnl, 1, MPI_INT, mnls, 1, MPI_INT, comm));
73101a81e61SBarry Smith 
7329566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size, &start_indices));
7332fa5cd67SKarl Rupp 
73401a81e61SBarry Smith   start_indices[0] = 0;
7352fa5cd67SKarl Rupp   for (i = 1; i < size; i++) start_indices[i] = start_indices[i - 1] + mnls[i - 1];
73601a81e61SBarry Smith 
7379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(v->mnl, &global_indices));
7382fa5cd67SKarl Rupp   for (i = 0; i < v->mnl; i++) global_indices[i] = start_indices[rank] + i;
73901a81e61SBarry Smith 
7409566063dSJacob Faibussowitsch   PetscCall(PetscFree(mnls));
7419566063dSJacob Faibussowitsch   PetscCall(PetscFree(start_indices));
74201a81e61SBarry Smith 
7439566063dSJacob Faibussowitsch   PetscCall(VecSetValues(*Pv, v->mnl, global_indices, v->v, INSERT_VALUES));
7449566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*Pv));
7459566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*Pv));
74601a81e61SBarry Smith 
7479566063dSJacob Faibussowitsch   PetscCall(PetscFree(global_indices));
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74901a81e61SBarry Smith }
750