101a81e61SBarry Smith /*
201a81e61SBarry Smith 3/99 Modified by Stephen Barnard to support SPAI version 3.0
301a81e61SBarry Smith */
401a81e61SBarry Smith
501a81e61SBarry Smith /*
601a81e61SBarry Smith Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner
701a81e61SBarry Smith Code written by Stephen Barnard.
801a81e61SBarry Smith
901a81e61SBarry Smith Note: there is some BAD memory bleeding below!
1001a81e61SBarry Smith
1101a81e61SBarry Smith This code needs work
1201a81e61SBarry Smith
1301a81e61SBarry Smith 1) get rid of all memory bleeding
1401a81e61SBarry Smith 2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix
1501a81e61SBarry Smith rather than having the sp flag for PC_SPAI
162a4967abSBarry Smith 3) fix to set the block size based on the matrix block size
1701a81e61SBarry Smith
1801a81e61SBarry Smith */
1974c01ffaSSatish Balay #if !defined(PETSC_SKIP_COMPLEX)
20762360ffSBarry Smith #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */
2174c01ffaSSatish Balay #endif
2201a81e61SBarry Smith
23af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
241c7d2463SJed Brown #include <../src/ksp/pc/impls/spai/petscspai.h>
2501a81e61SBarry Smith
2601a81e61SBarry Smith /*
2701a81e61SBarry Smith These are the SPAI include files
2801a81e61SBarry Smith */
2901a81e61SBarry Smith EXTERN_C_BEGIN
30b1bda85cSSatish Balay #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */
31c6db04a5SJed Brown #include <spai.h>
32c6db04a5SJed Brown #include <matrix.h>
3301a81e61SBarry Smith EXTERN_C_END
3401a81e61SBarry Smith
35ba38deedSJacob Faibussowitsch static PetscErrorCode ConvertMatToMatrix(MPI_Comm, Mat, Mat, matrix **);
36ba38deedSJacob Faibussowitsch static PetscErrorCode ConvertMatrixToMat(MPI_Comm, matrix *, Mat *);
3701a81e61SBarry Smith
3801a81e61SBarry Smith typedef struct {
3901a81e61SBarry Smith matrix *B; /* matrix in SPAI format */
4001a81e61SBarry Smith matrix *BT; /* transpose of matrix in SPAI format */
4101a81e61SBarry Smith matrix *M; /* the approximate inverse in SPAI format */
4201a81e61SBarry Smith
4301a81e61SBarry Smith Mat PM; /* the approximate inverse PETSc format */
4401a81e61SBarry Smith
4501a81e61SBarry Smith double epsilon; /* tolerance */
4601a81e61SBarry Smith int nbsteps; /* max number of "improvement" steps per line */
4701a81e61SBarry Smith int max; /* max dimensions of is_I, q, etc. */
4801a81e61SBarry Smith int maxnew; /* max number of new entries per step */
4901a81e61SBarry Smith int block_size; /* constant block size */
5001a81e61SBarry Smith int cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */
5101a81e61SBarry Smith int verbose; /* SPAI prints timing and statistics */
5201a81e61SBarry Smith
5301a81e61SBarry Smith int sp; /* symmetric nonzero pattern */
5401a81e61SBarry Smith MPI_Comm comm_spai; /* communicator to be used with spai */
5501a81e61SBarry Smith } PC_SPAI;
5601a81e61SBarry Smith
PCSetUp_SPAI(PC pc)57d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_SPAI(PC pc)
58d71ae5a4SJacob Faibussowitsch {
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 */
83d5b43468SJose E. Roman /* int block_size */ /* block_size == 1 specifies scalar elements
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
903ba16761SJacob Faibussowitsch PetscCallExternal(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);
973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
9801a81e61SBarry Smith }
9901a81e61SBarry Smith
PCApply_SPAI(PC pc,Vec xx,Vec y)100d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SPAI(PC pc, Vec xx, Vec y)
101d71ae5a4SJacob Faibussowitsch {
10201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
10301a81e61SBarry Smith
10401a81e61SBarry Smith PetscFunctionBegin;
10501a81e61SBarry Smith /* Now using PETSc's multiply */
1069566063dSJacob Faibussowitsch PetscCall(MatMult(ispai->PM, xx, y));
1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
10801a81e61SBarry Smith }
10901a81e61SBarry Smith
PCMatApply_SPAI(PC pc,Mat X,Mat Y)110d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_SPAI(PC pc, Mat X, Mat Y)
111d71ae5a4SJacob Faibussowitsch {
11248e38a8aSPierre Jolivet PC_SPAI *ispai = (PC_SPAI *)pc->data;
11348e38a8aSPierre Jolivet
11448e38a8aSPierre Jolivet PetscFunctionBegin;
11548e38a8aSPierre Jolivet /* Now using PETSc's multiply */
116fb842aefSJose E. Roman PetscCall(MatMatMult(ispai->PM, X, MAT_REUSE_MATRIX, PETSC_CURRENT, &Y));
1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
11848e38a8aSPierre Jolivet }
11948e38a8aSPierre Jolivet
PCDestroy_SPAI(PC pc)120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SPAI(PC pc)
121d71ae5a4SJacob Faibussowitsch {
12201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
12301a81e61SBarry Smith
12401a81e61SBarry Smith PetscFunctionBegin;
1259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ispai->PM));
126f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_free(&ispai->comm_spai));
1279566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data));
1282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", NULL));
1292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", NULL));
1302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", NULL));
1312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", NULL));
1322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", NULL));
1332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", NULL));
1342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", NULL));
1352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", NULL));
1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
13701a81e61SBarry Smith }
13801a81e61SBarry Smith
PCView_SPAI(PC pc,PetscViewer viewer)139d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_SPAI(PC pc, PetscViewer viewer)
140d71ae5a4SJacob Faibussowitsch {
14101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
142*9f196a02SMartin Diehl PetscBool isascii;
14301a81e61SBarry Smith
14401a81e61SBarry Smith PetscFunctionBegin;
145*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
146*9f196a02SMartin Diehl if (isascii) {
147b03515a0SUmberto Zerbinati PetscCall(PetscViewerASCIIPrintf(viewer, " epsilon %g\n", ispai->epsilon));
1489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " nbsteps %d\n", ispai->nbsteps));
1499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " max %d\n", ispai->max));
1509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxnew %d\n", ispai->maxnew));
1519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " block_size %d\n", ispai->block_size));
1529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " cache_size %d\n", ispai->cache_size));
1539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " verbose %d\n", ispai->verbose));
1549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " sp %d\n", ispai->sp));
15501a81e61SBarry Smith }
1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
15701a81e61SBarry Smith }
15801a81e61SBarry Smith
PCSPAISetEpsilon_SPAI(PC pc,PetscReal epsilon1)159d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1)
160d71ae5a4SJacob Faibussowitsch {
16101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
1625fd66863SKarl Rupp
16301a81e61SBarry Smith PetscFunctionBegin;
164b03515a0SUmberto Zerbinati ispai->epsilon = (double)epsilon1;
1653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
16601a81e61SBarry Smith }
16701a81e61SBarry Smith
PCSPAISetNBSteps_SPAI(PC pc,PetscInt nbsteps1)168d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc, PetscInt nbsteps1)
169d71ae5a4SJacob Faibussowitsch {
17001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
1715fd66863SKarl Rupp
17201a81e61SBarry Smith PetscFunctionBegin;
173b03515a0SUmberto Zerbinati ispai->nbsteps = (int)nbsteps1;
1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
17501a81e61SBarry Smith }
17601a81e61SBarry Smith
17701a81e61SBarry Smith /* added 1/7/99 g.h. */
PCSPAISetMax_SPAI(PC pc,PetscInt max1)178d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMax_SPAI(PC pc, PetscInt max1)
179d71ae5a4SJacob Faibussowitsch {
18001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
1815fd66863SKarl Rupp
18201a81e61SBarry Smith PetscFunctionBegin;
183b03515a0SUmberto Zerbinati ispai->max = (int)max1;
1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
18501a81e61SBarry Smith }
18601a81e61SBarry Smith
PCSPAISetMaxNew_SPAI(PC pc,PetscInt maxnew1)187d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc, PetscInt maxnew1)
188d71ae5a4SJacob Faibussowitsch {
18901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
1905fd66863SKarl Rupp
19101a81e61SBarry Smith PetscFunctionBegin;
192b03515a0SUmberto Zerbinati ispai->maxnew = (int)maxnew1;
1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
19401a81e61SBarry Smith }
19501a81e61SBarry Smith
PCSPAISetBlockSize_SPAI(PC pc,PetscInt block_size1)196d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc, PetscInt block_size1)
197d71ae5a4SJacob Faibussowitsch {
19801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
1995fd66863SKarl Rupp
20001a81e61SBarry Smith PetscFunctionBegin;
201b03515a0SUmberto Zerbinati ispai->block_size = (int)block_size1;
2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
20301a81e61SBarry Smith }
20401a81e61SBarry Smith
PCSPAISetCacheSize_SPAI(PC pc,PetscInt cache_size)205d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc, PetscInt cache_size)
206d71ae5a4SJacob Faibussowitsch {
20701a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
2085fd66863SKarl Rupp
20901a81e61SBarry Smith PetscFunctionBegin;
210b03515a0SUmberto Zerbinati ispai->cache_size = (int)cache_size;
2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
21201a81e61SBarry Smith }
21301a81e61SBarry Smith
PCSPAISetVerbose_SPAI(PC pc,PetscInt verbose)214d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc, PetscInt verbose)
215d71ae5a4SJacob Faibussowitsch {
21601a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
2175fd66863SKarl Rupp
21801a81e61SBarry Smith PetscFunctionBegin;
219b03515a0SUmberto Zerbinati ispai->verbose = (int)verbose;
2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
22101a81e61SBarry Smith }
22201a81e61SBarry Smith
PCSPAISetSp_SPAI(PC pc,PetscInt sp)223d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSPAISetSp_SPAI(PC pc, PetscInt sp)
224d71ae5a4SJacob Faibussowitsch {
22501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data;
2265fd66863SKarl Rupp
22701a81e61SBarry Smith PetscFunctionBegin;
228b03515a0SUmberto Zerbinati ispai->sp = (int)sp;
2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
23001a81e61SBarry Smith }
23101a81e61SBarry Smith
23201a81e61SBarry Smith /*@
233f1580f4eSBarry Smith PCSPAISetEpsilon -- Set the tolerance for the `PCSPAI` preconditioner
23401a81e61SBarry Smith
23501a81e61SBarry Smith Input Parameters:
23601a81e61SBarry Smith + pc - the preconditioner
2371d27aa22SBarry Smith - epsilon1 - the tolerance (default .4)
23801a81e61SBarry Smith
23901a81e61SBarry Smith Level: intermediate
24001a81e61SBarry Smith
2411d27aa22SBarry Smith Note:
2421d27aa22SBarry Smith `espilon1` must be between 0 and 1. It controls the
2431d27aa22SBarry Smith quality of the approximation of M to the inverse of
2441d27aa22SBarry Smith A. Higher values of `epsilon1` lead to more work, more
2451d27aa22SBarry Smith fill, and usually better preconditioners. In many
2461d27aa22SBarry Smith cases the best choice of `epsilon1` is the one that
2471d27aa22SBarry Smith divides the total solution time equally between the
2481d27aa22SBarry Smith preconditioner and the solver.
2491d27aa22SBarry Smith
250562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
25101a81e61SBarry Smith @*/
PCSPAISetEpsilon(PC pc,PetscReal epsilon1)252d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1)
253d71ae5a4SJacob Faibussowitsch {
25401a81e61SBarry Smith PetscFunctionBegin;
255b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1));
2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
25701a81e61SBarry Smith }
25801a81e61SBarry Smith
25901a81e61SBarry Smith /*@
26001a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in
261f1580f4eSBarry Smith the `PCSPAI` preconditioner
26201a81e61SBarry Smith
26301a81e61SBarry Smith Input Parameters:
26401a81e61SBarry Smith + pc - the preconditioner
265feefa0e1SJacob Faibussowitsch - nbsteps1 - number of steps (default 5)
26601a81e61SBarry Smith
267f1580f4eSBarry Smith Note:
268f1580f4eSBarry Smith `PCSPAI` constructs to approximation to every column of
26901a81e61SBarry Smith the exact inverse of A in a series of improvement
27001a81e61SBarry Smith steps. The quality of the approximation is determined
27101a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy
2721d27aa22SBarry Smith of epsilon is not obtained after `nbsteps1` steps, `PCSPAI` simply
27301a81e61SBarry Smith uses the best approximation constructed so far.
27401a81e61SBarry Smith
27501a81e61SBarry Smith Level: intermediate
27601a81e61SBarry Smith
277562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
27801a81e61SBarry Smith @*/
PCSPAISetNBSteps(PC pc,PetscInt nbsteps1)279d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1)
280d71ae5a4SJacob Faibussowitsch {
28101a81e61SBarry Smith PetscFunctionBegin;
282b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1));
2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
28401a81e61SBarry Smith }
28501a81e61SBarry Smith
28601a81e61SBarry Smith /* added 1/7/99 g.h. */
28701a81e61SBarry Smith /*@
2881d27aa22SBarry Smith PCSPAISetMax - set the size of various working buffers in the `PCSPAI` preconditioner
28901a81e61SBarry Smith
29001a81e61SBarry Smith Input Parameters:
29101a81e61SBarry Smith + pc - the preconditioner
292feefa0e1SJacob Faibussowitsch - max1 - size (default is 5000)
29301a81e61SBarry Smith
29401a81e61SBarry Smith Level: intermediate
29501a81e61SBarry Smith
296562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
29701a81e61SBarry Smith @*/
PCSPAISetMax(PC pc,PetscInt max1)298d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1)
299d71ae5a4SJacob Faibussowitsch {
30001a81e61SBarry Smith PetscFunctionBegin;
301b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1));
3023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30301a81e61SBarry Smith }
30401a81e61SBarry Smith
30501a81e61SBarry Smith /*@
3061d27aa22SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step in the `PCSPAI` preconditioner
30701a81e61SBarry Smith
30801a81e61SBarry Smith Input Parameters:
30901a81e61SBarry Smith + pc - the preconditioner
310feefa0e1SJacob Faibussowitsch - maxnew1 - maximum number (default 5)
31101a81e61SBarry Smith
31201a81e61SBarry Smith Level: intermediate
31301a81e61SBarry Smith
314562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
31501a81e61SBarry Smith @*/
PCSPAISetMaxNew(PC pc,PetscInt maxnew1)316d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1)
317d71ae5a4SJacob Faibussowitsch {
31801a81e61SBarry Smith PetscFunctionBegin;
319b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1));
3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
32101a81e61SBarry Smith }
32201a81e61SBarry Smith
32301a81e61SBarry Smith /*@
324f1580f4eSBarry Smith PCSPAISetBlockSize - set the block size for the `PCSPAI` preconditioner
32501a81e61SBarry Smith
32601a81e61SBarry Smith Input Parameters:
32701a81e61SBarry Smith + pc - the preconditioner
328feefa0e1SJacob Faibussowitsch - block_size1 - block size (default 1)
32901a81e61SBarry Smith
3301d27aa22SBarry Smith Level: intermediate
3311d27aa22SBarry Smith
33295452b02SPatrick Sanan Notes:
33395452b02SPatrick Sanan A block
33401a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A
33501a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs
33601a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with
33701a81e61SBarry Smith variable sized blocks, which are determined by
33801a81e61SBarry Smith searching for dense square diagonal blocks in A.
33901a81e61SBarry Smith This can be very effective for finite-element
34001a81e61SBarry Smith matrices.
34101a81e61SBarry Smith
34201a81e61SBarry Smith SPAI will convert A to block form, use a block
34301a81e61SBarry Smith version of the preconditioner algorithm, and then
34401a81e61SBarry Smith convert the result back to scalar form.
34501a81e61SBarry Smith
34601a81e61SBarry Smith In many cases the a block-size parameter other than 1
34701a81e61SBarry Smith can lead to very significant improvement in
34801a81e61SBarry Smith performance.
34901a81e61SBarry Smith
3501d27aa22SBarry Smith Developer Note:
3511d27aa22SBarry Smith This preconditioner could use the matrix block size as the default block size to use
35201a81e61SBarry Smith
353562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
35401a81e61SBarry Smith @*/
PCSPAISetBlockSize(PC pc,PetscInt block_size1)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));
3593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
367feefa0e1SJacob Faibussowitsch - cache_size - cache size {0,1,2,3,4,5} (default 5)
36801a81e61SBarry Smith
3691d27aa22SBarry Smith Level: intermediate
3701d27aa22SBarry 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
377562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
37801a81e61SBarry Smith @*/
PCSPAISetCacheSize(PC pc,PetscInt cache_size)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));
3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
391feefa0e1SJacob Faibussowitsch - verbose - level (default 1)
39201a81e61SBarry Smith
39301a81e61SBarry Smith Level: intermediate
39401a81e61SBarry Smith
3951d27aa22SBarry Smith Note:
3961d27aa22SBarry Smith Prints parameters, timings and matrix statistics
3971d27aa22SBarry Smith
398562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
39901a81e61SBarry Smith @*/
PCSPAISetVerbose(PC pc,PetscInt verbose)400d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose)
401d71ae5a4SJacob Faibussowitsch {
40201a81e61SBarry Smith PetscFunctionBegin;
403b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose));
4043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
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
412feefa0e1SJacob Faibussowitsch - sp - 0 or 1
41301a81e61SBarry Smith
4141d27aa22SBarry Smith Level: intermediate
4151d27aa22SBarry Smith
416f1580f4eSBarry Smith Note:
4171d27aa22SBarry Smith 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
4201d27aa22SBarry 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
424562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSPAI`, `PCSetType()`
42501a81e61SBarry Smith @*/
PCSPAISetSp(PC pc,PetscInt sp)426d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp)
427d71ae5a4SJacob Faibussowitsch {
42801a81e61SBarry Smith PetscFunctionBegin;
429b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp));
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
43101a81e61SBarry Smith }
43201a81e61SBarry Smith
PCSetFromOptions_SPAI(PC pc,PetscOptionItems PetscOptionsObject)433ce78bad3SBarry Smith 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();
4603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
46101a81e61SBarry Smith }
46201a81e61SBarry Smith
46301a81e61SBarry Smith /*MC
4641d27aa22SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method {cite}`gh97`
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
481562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
482db781477SPatrick Sanan `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`,
483f1580f4eSBarry Smith `PCSPAISetVerbose()`, `PCSPAISetSp()`, `PCSPAISetNBSteps()`, `PCSPAISetCacheSize()`
48401a81e61SBarry Smith M*/
48501a81e61SBarry Smith
PCCreate_SPAI(PC pc)486d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc)
487d71ae5a4SJacob Faibussowitsch {
48801a81e61SBarry Smith PC_SPAI *ispai;
48901a81e61SBarry Smith
49001a81e61SBarry Smith PetscFunctionBegin;
4914dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&ispai));
4922a4967abSBarry Smith pc->data = ispai;
49301a81e61SBarry Smith
49401a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI;
49501a81e61SBarry Smith pc->ops->apply = PCApply_SPAI;
49648e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_SPAI;
49701a81e61SBarry Smith pc->ops->applyrichardson = 0;
49801a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI;
49901a81e61SBarry Smith pc->ops->view = PCView_SPAI;
50001a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI;
50101a81e61SBarry Smith
50201a81e61SBarry Smith ispai->epsilon = .4;
50301a81e61SBarry Smith ispai->nbsteps = 5;
50401a81e61SBarry Smith ispai->max = 5000;
50501a81e61SBarry Smith ispai->maxnew = 5;
50601a81e61SBarry Smith ispai->block_size = 1;
50701a81e61SBarry Smith ispai->cache_size = 5;
50801a81e61SBarry Smith ispai->verbose = 0;
50901a81e61SBarry Smith
51001a81e61SBarry Smith ispai->sp = 1;
511f4f49eeaSPierre Jolivet PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &ispai->comm_spai));
51201a81e61SBarry Smith
5139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", PCSPAISetEpsilon_SPAI));
5149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", PCSPAISetNBSteps_SPAI));
5159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", PCSPAISetMax_SPAI));
5169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", PCSPAISetMaxNew_SPAI));
5179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", PCSPAISetBlockSize_SPAI));
5189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", PCSPAISetCacheSize_SPAI));
5199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", PCSPAISetVerbose_SPAI));
5209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", PCSPAISetSp_SPAI));
5213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
52201a81e61SBarry Smith }
52301a81e61SBarry Smith
52401a81e61SBarry Smith /*
52501a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix
52601a81e61SBarry Smith */
ConvertMatToMatrix(MPI_Comm comm,Mat A,Mat AT,matrix ** B)527ba38deedSJacob Faibussowitsch static PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A, Mat AT, matrix **B)
528d71ae5a4SJacob Faibussowitsch {
52901a81e61SBarry Smith matrix *M;
53001a81e61SBarry Smith int i, j, col;
53101a81e61SBarry Smith int row_indx;
53201a81e61SBarry Smith int len, pe, local_indx, start_indx;
53301a81e61SBarry Smith int *mapping;
53401a81e61SBarry Smith const int *cols;
53501a81e61SBarry Smith const double *vals;
5362122902bSSatish Balay int n, mnl, nnl, nz, rstart, rend;
5372a4967abSBarry Smith PetscMPIInt size, rank;
53801a81e61SBarry Smith struct compressed_lines *rows;
53901a81e61SBarry Smith
54001a81e61SBarry Smith PetscFunctionBegin;
5419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
5429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
5439566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &n, &n));
5449566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &mnl, &nnl));
54501a81e61SBarry Smith
54601a81e61SBarry Smith /*
54701a81e61SBarry Smith not sure why a barrier is required. commenting out
5489566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(comm));
54901a81e61SBarry Smith */
55001a81e61SBarry Smith
55129b2f704SSatish Balay M = new_matrix((SPAI_Comm)comm);
55201a81e61SBarry Smith
55301a81e61SBarry Smith M->n = n;
55401a81e61SBarry Smith M->bs = 1;
55501a81e61SBarry Smith M->max_block_size = 1;
55601a81e61SBarry Smith
55701a81e61SBarry Smith M->mnls = (int *)malloc(sizeof(int) * size);
55801a81e61SBarry Smith M->start_indices = (int *)malloc(sizeof(int) * size);
55901a81e61SBarry Smith M->pe = (int *)malloc(sizeof(int) * n);
56001a81e61SBarry Smith M->block_sizes = (int *)malloc(sizeof(int) * n);
56101a81e61SBarry Smith for (i = 0; i < n; i++) M->block_sizes[i] = 1;
56201a81e61SBarry Smith
5639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnl, 1, MPI_INT, M->mnls, 1, MPI_INT, comm));
56401a81e61SBarry Smith
56501a81e61SBarry Smith M->start_indices[0] = 0;
5662fa5cd67SKarl Rupp for (i = 1; i < size; i++) M->start_indices[i] = M->start_indices[i - 1] + M->mnls[i - 1];
56701a81e61SBarry Smith
56801a81e61SBarry Smith M->mnl = M->mnls[M->myid];
56901a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid];
57001a81e61SBarry Smith
57101a81e61SBarry Smith for (i = 0; i < size; i++) {
57201a81e61SBarry Smith start_indx = M->start_indices[i];
5732fa5cd67SKarl Rupp for (j = 0; j < M->mnls[i]; j++) M->pe[start_indx + j] = i;
57401a81e61SBarry Smith }
57501a81e61SBarry Smith
57601a81e61SBarry Smith if (AT) {
5772f613bf5SBarry Smith M->lines = new_compressed_lines(M->mnls[rank], 1);
57801a81e61SBarry Smith } else {
5792f613bf5SBarry Smith M->lines = new_compressed_lines(M->mnls[rank], 0);
58001a81e61SBarry Smith }
58101a81e61SBarry Smith
58201a81e61SBarry Smith rows = M->lines;
58301a81e61SBarry Smith
58401a81e61SBarry Smith /* Determine the mapping from global indices to pointers */
5859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M->n, &mapping));
58601a81e61SBarry Smith pe = 0;
58701a81e61SBarry Smith local_indx = 0;
58801a81e61SBarry Smith for (i = 0; i < M->n; i++) {
58901a81e61SBarry Smith if (local_indx >= M->mnls[pe]) {
59001a81e61SBarry Smith pe++;
59101a81e61SBarry Smith local_indx = 0;
59201a81e61SBarry Smith }
59301a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe];
59401a81e61SBarry Smith local_indx++;
59501a81e61SBarry Smith }
59601a81e61SBarry Smith
59701a81e61SBarry Smith /************** Set up the row structure *****************/
59801a81e61SBarry Smith
5999566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
60001a81e61SBarry Smith for (i = rstart; i < rend; i++) {
60101a81e61SBarry Smith row_indx = i - rstart;
6029566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nz, &cols, &vals));
6032122902bSSatish Balay /* allocate buffers */
6042122902bSSatish Balay rows->ptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6052122902bSSatish Balay rows->A[row_indx] = (double *)malloc(nz * sizeof(double));
6062122902bSSatish Balay /* copy the matrix */
60701a81e61SBarry Smith for (j = 0; j < nz; j++) {
60801a81e61SBarry Smith col = cols[j];
60901a81e61SBarry Smith len = rows->len[row_indx]++;
6102fa5cd67SKarl Rupp
61101a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col];
61201a81e61SBarry Smith rows->A[row_indx][len] = vals[j];
61301a81e61SBarry Smith }
61401a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx];
6152fa5cd67SKarl Rupp
6169566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nz, &cols, &vals));
61701a81e61SBarry Smith }
61801a81e61SBarry Smith
61901a81e61SBarry Smith /************** Set up the column structure *****************/
62001a81e61SBarry Smith
62101a81e61SBarry Smith if (AT) {
62201a81e61SBarry Smith for (i = rstart; i < rend; i++) {
62301a81e61SBarry Smith row_indx = i - rstart;
6249566063dSJacob Faibussowitsch PetscCall(MatGetRow(AT, i, &nz, &cols, &vals));
6252122902bSSatish Balay /* allocate buffers */
6262122902bSSatish Balay rows->rptrs[row_indx] = (int *)malloc(nz * sizeof(int));
6272122902bSSatish Balay /* copy the matrix (i.e., the structure) */
62801a81e61SBarry Smith for (j = 0; j < nz; j++) {
62901a81e61SBarry Smith col = cols[j];
63001a81e61SBarry Smith len = rows->rlen[row_indx]++;
6312fa5cd67SKarl Rupp
63201a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col];
63301a81e61SBarry Smith }
6349566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(AT, i, &nz, &cols, &vals));
63501a81e61SBarry Smith }
63601a81e61SBarry Smith }
63701a81e61SBarry Smith
6389566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping));
63901a81e61SBarry Smith
64001a81e61SBarry Smith order_pointers(M);
64101a81e61SBarry Smith M->maxnz = calc_maxnz(M);
64201a81e61SBarry Smith *B = M;
6433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
64401a81e61SBarry Smith }
64501a81e61SBarry Smith
64601a81e61SBarry Smith /*
64701a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB.
648df3898eeSBarry Smith This assumes that the SPAI matrix B is stored in
64901a81e61SBarry Smith COMPRESSED-ROW format.
65001a81e61SBarry Smith */
ConvertMatrixToMat(MPI_Comm comm,matrix * B,Mat * PB)651ba38deedSJacob Faibussowitsch static PetscErrorCode ConvertMatrixToMat(MPI_Comm comm, matrix *B, Mat *PB)
652d71ae5a4SJacob Faibussowitsch {
6534b6b5fe1SBarry Smith PetscMPIInt size, rank;
65401a81e61SBarry Smith int m, n, M, N;
65501a81e61SBarry Smith int d_nz, o_nz;
65601a81e61SBarry Smith int *d_nnz, *o_nnz;
65701a81e61SBarry Smith int i, k, global_row, global_col, first_diag_col, last_diag_col;
65801a81e61SBarry Smith PetscScalar val;
65901a81e61SBarry Smith
66001a81e61SBarry Smith PetscFunctionBegin;
6619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size));
6629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
66301a81e61SBarry Smith
66401a81e61SBarry Smith m = n = B->mnls[rank];
66501a81e61SBarry Smith d_nz = o_nz = 0;
66601a81e61SBarry Smith
66706946f3aSJose E. Roman /* Determine preallocation for MatCreateAIJ */
6689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &d_nnz));
6699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &o_nnz));
67001a81e61SBarry Smith for (i = 0; i < m; i++) d_nnz[i] = o_nnz[i] = 0;
67101a81e61SBarry Smith first_diag_col = B->start_indices[rank];
67201a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank];
67301a81e61SBarry Smith for (i = 0; i < B->mnls[rank]; i++) {
67401a81e61SBarry Smith for (k = 0; k < B->lines->len[i]; k++) {
67501a81e61SBarry Smith global_col = B->lines->ptrs[i][k];
676db4deed7SKarl Rupp if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
677db4deed7SKarl Rupp else o_nnz[i]++;
67801a81e61SBarry Smith }
67901a81e61SBarry Smith }
68001a81e61SBarry Smith
68101a81e61SBarry Smith M = N = B->n;
68201a81e61SBarry Smith /* Here we only know how to create AIJ format */
6839566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, PB));
6849566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*PB, m, n, M, N));
6859566063dSJacob Faibussowitsch PetscCall(MatSetType(*PB, MATAIJ));
6869566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*PB, d_nz, d_nnz));
6879566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*PB, d_nz, d_nnz, o_nz, o_nnz));
68801a81e61SBarry Smith
68901a81e61SBarry Smith for (i = 0; i < B->mnls[rank]; i++) {
69001a81e61SBarry Smith global_row = B->start_indices[rank] + i;
69101a81e61SBarry Smith for (k = 0; k < B->lines->len[i]; k++) {
69201a81e61SBarry Smith global_col = B->lines->ptrs[i][k];
6932fa5cd67SKarl Rupp
69401a81e61SBarry Smith val = B->lines->A[i][k];
6959566063dSJacob Faibussowitsch PetscCall(MatSetValues(*PB, 1, &global_row, 1, &global_col, &val, ADD_VALUES));
69601a81e61SBarry Smith }
69701a81e61SBarry Smith }
6989566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz));
6999566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nnz));
70001a81e61SBarry Smith
7019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*PB, MAT_FINAL_ASSEMBLY));
7029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*PB, MAT_FINAL_ASSEMBLY));
7033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
70401a81e61SBarry Smith }
705