xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision 2e956fe4fc852fabc23b437482e1fb7b77fddb0d)
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 
4101a81e61SBarry Smith   matrix *B;                /* matrix in SPAI format */
4201a81e61SBarry Smith   matrix *BT;               /* transpose of matrix in SPAI format */
4301a81e61SBarry Smith   matrix *M;                /* the approximate inverse in SPAI format */
4401a81e61SBarry Smith 
4501a81e61SBarry Smith   Mat PM;                   /* the approximate inverse PETSc format */
4601a81e61SBarry Smith 
4701a81e61SBarry Smith   double epsilon;           /* tolerance */
4801a81e61SBarry Smith   int    nbsteps;           /* max number of "improvement" steps per line */
4901a81e61SBarry Smith   int    max;               /* max dimensions of is_I, q, etc. */
5001a81e61SBarry Smith   int    maxnew;            /* max number of new entries per step */
5101a81e61SBarry Smith   int    block_size;        /* constant block size */
5201a81e61SBarry Smith   int    cache_size;        /* one of (1,2,3,4,5,6) indicting size of cache */
5301a81e61SBarry Smith   int    verbose;           /* SPAI prints timing and statistics */
5401a81e61SBarry Smith 
5501a81e61SBarry Smith   int      sp;              /* symmetric nonzero pattern */
5601a81e61SBarry Smith   MPI_Comm comm_spai;     /* communicator to be used with spai */
5701a81e61SBarry Smith } PC_SPAI;
5801a81e61SBarry Smith 
5901a81e61SBarry Smith /**********************************************************************/
6001a81e61SBarry Smith 
6101a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc)
6201a81e61SBarry Smith {
6301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
6401a81e61SBarry Smith   Mat      AT;
6501a81e61SBarry Smith 
6601a81e61SBarry Smith   PetscFunctionBegin;
6701a81e61SBarry Smith   init_SPAI();
6801a81e61SBarry Smith 
6901a81e61SBarry Smith   if (ispai->sp) {
709566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B));
7101a81e61SBarry Smith   } else {
7201a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
739566063dSJacob Faibussowitsch     PetscCall(MatTranspose(pc->pmat,MAT_INITIAL_MATRIX,&AT));
749566063dSJacob Faibussowitsch     PetscCall(ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B));
759566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&AT));
7601a81e61SBarry Smith   }
7701a81e61SBarry Smith 
7801a81e61SBarry Smith   /* Destroy the transpose */
7901a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
8001a81e61SBarry Smith 
8101a81e61SBarry Smith   /* construct SPAI preconditioner */
8201a81e61SBarry Smith   /* FILE *messages */     /* file for warning messages */
8301a81e61SBarry Smith   /* double epsilon */     /* tolerance */
8401a81e61SBarry Smith   /* int nbsteps */        /* max number of "improvement" steps per line */
8501a81e61SBarry Smith   /* int max */            /* max dimensions of is_I, q, etc. */
8601a81e61SBarry Smith   /* int maxnew */         /* max number of new entries per step */
8701a81e61SBarry Smith   /* int block_size */     /* block_size == 1 specifies scalar elments
8801a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
8901a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
902fa5cd67SKarl Rupp   /* int cache_size */     /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */
9101a81e61SBarry Smith   /* int    verbose    */  /* verbose == 0 specifies that SPAI is silent
9201a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9301a81e61SBarry Smith 
949566063dSJacob Faibussowitsch   PetscCall(bspai(ispai->B,&ispai->M,
9501a81e61SBarry Smith                 stdout,
9601a81e61SBarry Smith                 ispai->epsilon,
9701a81e61SBarry Smith                 ispai->nbsteps,
9801a81e61SBarry Smith                 ispai->max,
9901a81e61SBarry Smith                 ispai->maxnew,
10001a81e61SBarry Smith                 ispai->block_size,
10101a81e61SBarry Smith                 ispai->cache_size,
1025f80ce2aSJacob Faibussowitsch                 ispai->verbose));
10301a81e61SBarry Smith 
1049566063dSJacob Faibussowitsch   PetscCall(ConvertMatrixToMat(PetscObjectComm((PetscObject)pc),ispai->M,&ispai->PM));
10501a81e61SBarry Smith 
10601a81e61SBarry Smith   /* free the SPAI matrices */
10701a81e61SBarry Smith   sp_free_matrix(ispai->B);
10801a81e61SBarry Smith   sp_free_matrix(ispai->M);
10901a81e61SBarry Smith   PetscFunctionReturn(0);
11001a81e61SBarry Smith }
11101a81e61SBarry Smith 
11201a81e61SBarry Smith /**********************************************************************/
11301a81e61SBarry Smith 
11401a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y)
11501a81e61SBarry Smith {
11601a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
11701a81e61SBarry Smith 
11801a81e61SBarry Smith   PetscFunctionBegin;
11901a81e61SBarry Smith   /* Now using PETSc's multiply */
1209566063dSJacob Faibussowitsch   PetscCall(MatMult(ispai->PM,xx,y));
12101a81e61SBarry Smith   PetscFunctionReturn(0);
12201a81e61SBarry Smith }
12301a81e61SBarry Smith 
12448e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_SPAI(PC pc,Mat X,Mat Y)
12548e38a8aSPierre Jolivet {
12648e38a8aSPierre Jolivet   PC_SPAI *ispai = (PC_SPAI*)pc->data;
12748e38a8aSPierre Jolivet 
12848e38a8aSPierre Jolivet   PetscFunctionBegin;
12948e38a8aSPierre Jolivet   /* Now using PETSc's multiply */
1309566063dSJacob Faibussowitsch   PetscCall(MatMatMult(ispai->PM,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y));
13148e38a8aSPierre Jolivet   PetscFunctionReturn(0);
13248e38a8aSPierre Jolivet }
13348e38a8aSPierre Jolivet 
13401a81e61SBarry Smith /**********************************************************************/
13501a81e61SBarry Smith 
13601a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc)
13701a81e61SBarry Smith {
13801a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
13901a81e61SBarry Smith 
14001a81e61SBarry Smith   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ispai->PM));
1429566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_free(&(ispai->comm_spai)));
1439566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
144*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",NULL));
145*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",NULL));
146*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",NULL));
147*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",NULL));
148*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",NULL));
149*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",NULL));
150*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",NULL));
151*2e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",NULL));
15201a81e61SBarry Smith   PetscFunctionReturn(0);
15301a81e61SBarry Smith }
15401a81e61SBarry Smith 
15501a81e61SBarry Smith /**********************************************************************/
15601a81e61SBarry Smith 
15701a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer)
15801a81e61SBarry Smith {
15901a81e61SBarry Smith   PC_SPAI   *ispai = (PC_SPAI*)pc->data;
160ace3abfcSBarry Smith   PetscBool  iascii;
16101a81e61SBarry Smith 
16201a81e61SBarry Smith   PetscFunctionBegin;
1639566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
16401a81e61SBarry Smith   if (iascii) {
165b03515a0SUmberto Zerbinati     PetscCall(PetscViewerASCIIPrintf(viewer,"    epsilon %g\n",   ispai->epsilon));
1669566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps));
1679566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max));
1689566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew));
1699566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size));
1709566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size));
1719566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose));
1729566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp));
17301a81e61SBarry Smith   }
17401a81e61SBarry Smith   PetscFunctionReturn(0);
17501a81e61SBarry Smith }
17601a81e61SBarry Smith 
177b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1)
17801a81e61SBarry Smith {
17901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
1805fd66863SKarl Rupp 
18101a81e61SBarry Smith   PetscFunctionBegin;
182b03515a0SUmberto Zerbinati   ispai->epsilon = (double)epsilon1;
18301a81e61SBarry Smith   PetscFunctionReturn(0);
18401a81e61SBarry Smith }
18501a81e61SBarry Smith 
18601a81e61SBarry Smith /**********************************************************************/
18701a81e61SBarry Smith 
188b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetNBSteps_SPAI(PC pc,PetscInt nbsteps1)
18901a81e61SBarry Smith {
19001a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
1915fd66863SKarl Rupp 
19201a81e61SBarry Smith   PetscFunctionBegin;
193b03515a0SUmberto Zerbinati   ispai->nbsteps = (int)nbsteps1;
19401a81e61SBarry Smith   PetscFunctionReturn(0);
19501a81e61SBarry Smith }
19601a81e61SBarry Smith 
19701a81e61SBarry Smith /**********************************************************************/
19801a81e61SBarry Smith 
19901a81e61SBarry Smith /* added 1/7/99 g.h. */
200b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetMax_SPAI(PC pc,PetscInt max1)
20101a81e61SBarry Smith {
20201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2035fd66863SKarl Rupp 
20401a81e61SBarry Smith   PetscFunctionBegin;
205b03515a0SUmberto Zerbinati   ispai->max = (int)max1;
20601a81e61SBarry Smith   PetscFunctionReturn(0);
20701a81e61SBarry Smith }
20801a81e61SBarry Smith 
20901a81e61SBarry Smith /**********************************************************************/
21001a81e61SBarry Smith 
211b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetMaxNew_SPAI(PC pc,PetscInt maxnew1)
21201a81e61SBarry Smith {
21301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2145fd66863SKarl Rupp 
21501a81e61SBarry Smith   PetscFunctionBegin;
216b03515a0SUmberto Zerbinati   ispai->maxnew = (int)maxnew1;
21701a81e61SBarry Smith   PetscFunctionReturn(0);
21801a81e61SBarry Smith }
21901a81e61SBarry Smith 
22001a81e61SBarry Smith /**********************************************************************/
22101a81e61SBarry Smith 
222b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetBlockSize_SPAI(PC pc,PetscInt block_size1)
22301a81e61SBarry Smith {
22401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2255fd66863SKarl Rupp 
22601a81e61SBarry Smith   PetscFunctionBegin;
227b03515a0SUmberto Zerbinati   ispai->block_size = (int)block_size1;
22801a81e61SBarry Smith   PetscFunctionReturn(0);
22901a81e61SBarry Smith }
23001a81e61SBarry Smith 
23101a81e61SBarry Smith /**********************************************************************/
23201a81e61SBarry Smith 
233b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetCacheSize_SPAI(PC pc,PetscInt cache_size)
23401a81e61SBarry Smith {
23501a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2365fd66863SKarl Rupp 
23701a81e61SBarry Smith   PetscFunctionBegin;
238b03515a0SUmberto Zerbinati   ispai->cache_size = (int)cache_size;
23901a81e61SBarry Smith   PetscFunctionReturn(0);
24001a81e61SBarry Smith }
24101a81e61SBarry Smith 
24201a81e61SBarry Smith /**********************************************************************/
24301a81e61SBarry Smith 
244b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetVerbose_SPAI(PC pc,PetscInt verbose)
24501a81e61SBarry Smith {
24601a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2475fd66863SKarl Rupp 
24801a81e61SBarry Smith   PetscFunctionBegin;
249b03515a0SUmberto Zerbinati   ispai->verbose = (int)verbose;
25001a81e61SBarry Smith   PetscFunctionReturn(0);
25101a81e61SBarry Smith }
25201a81e61SBarry Smith 
25301a81e61SBarry Smith /**********************************************************************/
25401a81e61SBarry Smith 
255b03515a0SUmberto Zerbinati static PetscErrorCode  PCSPAISetSp_SPAI(PC pc,PetscInt sp)
25601a81e61SBarry Smith {
25701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2585fd66863SKarl Rupp 
25901a81e61SBarry Smith   PetscFunctionBegin;
260b03515a0SUmberto Zerbinati   ispai->sp = (int)sp;
26101a81e61SBarry Smith   PetscFunctionReturn(0);
26201a81e61SBarry Smith }
26301a81e61SBarry Smith 
26401a81e61SBarry Smith /* -------------------------------------------------------------------*/
26501a81e61SBarry Smith 
26601a81e61SBarry Smith /*@
26701a81e61SBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner
26801a81e61SBarry Smith 
26901a81e61SBarry Smith   Input Parameters:
27001a81e61SBarry Smith + pc - the preconditioner
27101a81e61SBarry Smith - eps - epsilon (default .4)
27201a81e61SBarry Smith 
27395452b02SPatrick Sanan   Notes:
27495452b02SPatrick Sanan     Espilon must be between 0 and 1. It controls the
27501a81e61SBarry Smith                  quality of the approximation of M to the inverse of
27601a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
27701a81e61SBarry Smith                  fill, and usually better preconditioners. In many
27801a81e61SBarry Smith                  cases the best choice of epsilon is the one that
27901a81e61SBarry Smith                  divides the total solution time equally between the
28001a81e61SBarry Smith                  preconditioner and the solver.
28101a81e61SBarry Smith 
28201a81e61SBarry Smith   Level: intermediate
28301a81e61SBarry Smith 
284db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
28501a81e61SBarry Smith   @*/
286b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetEpsilon(PC pc,PetscReal epsilon1)
28701a81e61SBarry Smith {
28801a81e61SBarry Smith   PetscFunctionBegin;
289b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,PetscReal),(pc,epsilon1));
29001a81e61SBarry Smith   PetscFunctionReturn(0);
29101a81e61SBarry Smith }
29201a81e61SBarry Smith 
29301a81e61SBarry Smith /**********************************************************************/
29401a81e61SBarry Smith 
29501a81e61SBarry Smith /*@
29601a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
29701a81e61SBarry Smith         the SPAI preconditioner
29801a81e61SBarry Smith 
29901a81e61SBarry Smith   Input Parameters:
30001a81e61SBarry Smith + pc - the preconditioner
30101a81e61SBarry Smith - n - number of steps (default 5)
30201a81e61SBarry Smith 
30395452b02SPatrick Sanan   Notes:
30495452b02SPatrick Sanan     SPAI constructs to approximation to every column of
30501a81e61SBarry Smith                  the exact inverse of A in a series of improvement
30601a81e61SBarry Smith                  steps. The quality of the approximation is determined
30701a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
30801a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
30901a81e61SBarry Smith                  uses the best approximation constructed so far.
31001a81e61SBarry Smith 
31101a81e61SBarry Smith   Level: intermediate
31201a81e61SBarry Smith 
313db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()`
31401a81e61SBarry Smith @*/
315b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetNBSteps(PC pc,PetscInt nbsteps1)
31601a81e61SBarry Smith {
31701a81e61SBarry Smith   PetscFunctionBegin;
318b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,PetscInt),(pc,nbsteps1));
31901a81e61SBarry Smith   PetscFunctionReturn(0);
32001a81e61SBarry Smith }
32101a81e61SBarry Smith 
32201a81e61SBarry Smith /**********************************************************************/
32301a81e61SBarry Smith 
32401a81e61SBarry Smith /* added 1/7/99 g.h. */
32501a81e61SBarry Smith /*@
32601a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
32701a81e61SBarry Smith         the SPAI preconditioner
32801a81e61SBarry Smith 
32901a81e61SBarry Smith   Input Parameters:
33001a81e61SBarry Smith + pc - the preconditioner
33101a81e61SBarry Smith - n - size (default is 5000)
33201a81e61SBarry Smith 
33301a81e61SBarry Smith   Level: intermediate
33401a81e61SBarry Smith 
335db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
33601a81e61SBarry Smith @*/
337b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetMax(PC pc,PetscInt max1)
33801a81e61SBarry Smith {
33901a81e61SBarry Smith   PetscFunctionBegin;
340b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetMax_C",(PC,PetscInt),(pc,max1));
34101a81e61SBarry Smith   PetscFunctionReturn(0);
34201a81e61SBarry Smith }
34301a81e61SBarry Smith 
34401a81e61SBarry Smith /**********************************************************************/
34501a81e61SBarry Smith 
34601a81e61SBarry Smith /*@
34701a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
34801a81e61SBarry Smith    in SPAI preconditioner
34901a81e61SBarry Smith 
35001a81e61SBarry Smith   Input Parameters:
35101a81e61SBarry Smith + pc - the preconditioner
35201a81e61SBarry Smith - n - maximum number (default 5)
35301a81e61SBarry Smith 
35401a81e61SBarry Smith   Level: intermediate
35501a81e61SBarry Smith 
356db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()`
35701a81e61SBarry Smith @*/
358b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetMaxNew(PC pc,PetscInt maxnew1)
35901a81e61SBarry Smith {
36001a81e61SBarry Smith   PetscFunctionBegin;
361b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,PetscInt),(pc,maxnew1));
36201a81e61SBarry Smith   PetscFunctionReturn(0);
36301a81e61SBarry Smith }
36401a81e61SBarry Smith 
36501a81e61SBarry Smith /**********************************************************************/
36601a81e61SBarry Smith 
36701a81e61SBarry Smith /*@
36801a81e61SBarry Smith   PCSPAISetBlockSize - set the block size for the SPAI preconditioner
36901a81e61SBarry Smith 
37001a81e61SBarry Smith   Input Parameters:
37101a81e61SBarry Smith + pc - the preconditioner
37201a81e61SBarry Smith - n - block size (default 1)
37301a81e61SBarry Smith 
37495452b02SPatrick Sanan   Notes:
37595452b02SPatrick Sanan     A block
37601a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
37701a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
37801a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
37901a81e61SBarry Smith                  variable sized blocks, which are determined by
38001a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
38101a81e61SBarry Smith                  This can be very effective for finite-element
38201a81e61SBarry Smith                  matrices.
38301a81e61SBarry Smith 
38401a81e61SBarry Smith                  SPAI will convert A to block form, use a block
38501a81e61SBarry Smith                  version of the preconditioner algorithm, and then
38601a81e61SBarry Smith                  convert the result back to scalar form.
38701a81e61SBarry Smith 
38801a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
38901a81e61SBarry Smith                  can lead to very significant improvement in
39001a81e61SBarry Smith                  performance.
39101a81e61SBarry Smith 
39201a81e61SBarry Smith   Level: intermediate
39301a81e61SBarry Smith 
394db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
39501a81e61SBarry Smith @*/
396b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetBlockSize(PC pc,PetscInt block_size1)
39701a81e61SBarry Smith {
39801a81e61SBarry Smith   PetscFunctionBegin;
399b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,PetscInt),(pc,block_size1));
40001a81e61SBarry Smith   PetscFunctionReturn(0);
40101a81e61SBarry Smith }
40201a81e61SBarry Smith 
40301a81e61SBarry Smith /**********************************************************************/
40401a81e61SBarry Smith 
40501a81e61SBarry Smith /*@
40601a81e61SBarry Smith   PCSPAISetCacheSize - specify cache size in the SPAI preconditioner
40701a81e61SBarry Smith 
40801a81e61SBarry Smith   Input Parameters:
40901a81e61SBarry Smith + pc - the preconditioner
41001a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
41101a81e61SBarry Smith 
41295452b02SPatrick Sanan   Notes:
41395452b02SPatrick Sanan     SPAI uses a hash table to cache messages and avoid
41401a81e61SBarry Smith                  redundant communication. If suggest always using
41501a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
41601a81e61SBarry Smith                  version.
41701a81e61SBarry Smith 
41801a81e61SBarry Smith   Level: intermediate
41901a81e61SBarry Smith 
420db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
42101a81e61SBarry Smith @*/
422b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetCacheSize(PC pc,PetscInt cache_size)
42301a81e61SBarry Smith {
42401a81e61SBarry Smith   PetscFunctionBegin;
425b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,PetscInt),(pc,cache_size));
42601a81e61SBarry Smith   PetscFunctionReturn(0);
42701a81e61SBarry Smith }
42801a81e61SBarry Smith 
42901a81e61SBarry Smith /**********************************************************************/
43001a81e61SBarry Smith 
43101a81e61SBarry Smith /*@
43201a81e61SBarry Smith   PCSPAISetVerbose - verbosity level for the SPAI preconditioner
43301a81e61SBarry Smith 
43401a81e61SBarry Smith   Input Parameters:
43501a81e61SBarry Smith + pc - the preconditioner
43601a81e61SBarry Smith - n - level (default 1)
43701a81e61SBarry Smith 
43895452b02SPatrick Sanan   Notes:
43995452b02SPatrick Sanan     print parameters, timings and matrix statistics
44001a81e61SBarry Smith 
44101a81e61SBarry Smith   Level: intermediate
44201a81e61SBarry Smith 
443db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
44401a81e61SBarry Smith @*/
445b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetVerbose(PC pc,PetscInt verbose)
44601a81e61SBarry Smith {
44701a81e61SBarry Smith   PetscFunctionBegin;
448b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,PetscInt),(pc,verbose));
44901a81e61SBarry Smith   PetscFunctionReturn(0);
45001a81e61SBarry Smith }
45101a81e61SBarry Smith 
45201a81e61SBarry Smith /**********************************************************************/
45301a81e61SBarry Smith 
45401a81e61SBarry Smith /*@
45501a81e61SBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner
45601a81e61SBarry Smith 
45701a81e61SBarry Smith   Input Parameters:
45801a81e61SBarry Smith + pc - the preconditioner
45901a81e61SBarry Smith - n - 0 or 1
46001a81e61SBarry Smith 
46195452b02SPatrick Sanan   Notes:
46295452b02SPatrick Sanan     If A has a symmetric nonzero pattern use -sp 1 to
46301a81e61SBarry Smith                  improve performance by eliminating some communication
46401a81e61SBarry Smith                  in the parallel version. Even if A does not have a
46501a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
46601a81e61SBarry Smith                  results, but the code will not follow the published
46701a81e61SBarry Smith                  SPAI algorithm exactly.
46801a81e61SBarry Smith 
46901a81e61SBarry Smith   Level: intermediate
47001a81e61SBarry Smith 
471db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`
47201a81e61SBarry Smith @*/
473b03515a0SUmberto Zerbinati PetscErrorCode  PCSPAISetSp(PC pc,PetscInt sp)
47401a81e61SBarry Smith {
47501a81e61SBarry Smith   PetscFunctionBegin;
476b03515a0SUmberto Zerbinati   PetscTryMethod(pc,"PCSPAISetSp_C",(PC,PetscInt),(pc,sp));
47701a81e61SBarry Smith   PetscFunctionReturn(0);
47801a81e61SBarry Smith }
47901a81e61SBarry Smith 
48001a81e61SBarry Smith /**********************************************************************/
48101a81e61SBarry Smith 
48201a81e61SBarry Smith /**********************************************************************/
48301a81e61SBarry Smith 
4844416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PetscOptionItems *PetscOptionsObject,PC pc)
48501a81e61SBarry Smith {
48601a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
48701a81e61SBarry Smith   int            nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
48801a81e61SBarry Smith   double         epsilon1;
489ace3abfcSBarry Smith   PetscBool      flg;
49001a81e61SBarry Smith 
49101a81e61SBarry Smith   PetscFunctionBegin;
492d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"SPAI options");
4939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg));
49401a81e61SBarry Smith   if (flg) {
4959566063dSJacob Faibussowitsch     PetscCall(PCSPAISetEpsilon(pc,epsilon1));
49601a81e61SBarry Smith   }
4979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg));
49801a81e61SBarry Smith   if (flg) {
4999566063dSJacob Faibussowitsch     PetscCall(PCSPAISetNBSteps(pc,nbsteps1));
50001a81e61SBarry Smith   }
50101a81e61SBarry Smith   /* added 1/7/99 g.h. */
5029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg));
50301a81e61SBarry Smith   if (flg) {
5049566063dSJacob Faibussowitsch     PetscCall(PCSPAISetMax(pc,max1));
50501a81e61SBarry Smith   }
5069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg));
50701a81e61SBarry Smith   if (flg) {
5089566063dSJacob Faibussowitsch     PetscCall(PCSPAISetMaxNew(pc,maxnew1));
50901a81e61SBarry Smith   }
5109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg));
51101a81e61SBarry Smith   if (flg) {
5129566063dSJacob Faibussowitsch     PetscCall(PCSPAISetBlockSize(pc,block_size1));
51301a81e61SBarry Smith   }
5149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg));
51501a81e61SBarry Smith   if (flg) {
5169566063dSJacob Faibussowitsch     PetscCall(PCSPAISetCacheSize(pc,cache_size));
51701a81e61SBarry Smith   }
5189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg));
51901a81e61SBarry Smith   if (flg) {
5209566063dSJacob Faibussowitsch     PetscCall(PCSPAISetVerbose(pc,verbose));
52101a81e61SBarry Smith   }
5229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg));
52301a81e61SBarry Smith   if (flg) {
5249566063dSJacob Faibussowitsch     PetscCall(PCSPAISetSp(pc,sp));
52501a81e61SBarry Smith   }
526d0609cedSBarry Smith   PetscOptionsHeadEnd();
52701a81e61SBarry Smith   PetscFunctionReturn(0);
52801a81e61SBarry Smith }
52901a81e61SBarry Smith 
53001a81e61SBarry Smith /**********************************************************************/
53101a81e61SBarry Smith 
53201a81e61SBarry Smith /*MC
53301a81e61SBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard
53401a81e61SBarry Smith      as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)
53501a81e61SBarry Smith 
53601a81e61SBarry Smith    Options Database Keys:
537c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
53801a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
53901a81e61SBarry Smith .  -pc_spai_max <m> - set max
54001a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
54101a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
54201a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
54301a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
54401a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
54501a81e61SBarry Smith 
54695452b02SPatrick Sanan    Notes:
54795452b02SPatrick Sanan     This only works with AIJ matrices.
54801a81e61SBarry Smith 
54901a81e61SBarry Smith    Level: beginner
55001a81e61SBarry Smith 
551db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
552db781477SPatrick Sanan           `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`,
553db781477SPatrick Sanan           `PCSPAISetVerbose()`, `PCSPAISetSp()`
55401a81e61SBarry Smith M*/
55501a81e61SBarry Smith 
5568cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc)
55701a81e61SBarry Smith {
55801a81e61SBarry Smith   PC_SPAI *ispai;
55901a81e61SBarry Smith 
56001a81e61SBarry Smith   PetscFunctionBegin;
5619566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc,&ispai));
5622a4967abSBarry Smith   pc->data = ispai;
56301a81e61SBarry Smith 
56401a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
56501a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
56648e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_SPAI;
56701a81e61SBarry Smith   pc->ops->applyrichardson = 0;
56801a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
56901a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
57001a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
57101a81e61SBarry Smith 
57201a81e61SBarry Smith   ispai->epsilon    = .4;
57301a81e61SBarry Smith   ispai->nbsteps    = 5;
57401a81e61SBarry Smith   ispai->max        = 5000;
57501a81e61SBarry Smith   ispai->maxnew     = 5;
57601a81e61SBarry Smith   ispai->block_size = 1;
57701a81e61SBarry Smith   ispai->cache_size = 5;
57801a81e61SBarry Smith   ispai->verbose    = 0;
57901a81e61SBarry Smith 
58001a81e61SBarry Smith   ispai->sp = 1;
5819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai)));
58201a81e61SBarry Smith 
5839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",PCSPAISetEpsilon_SPAI));
5849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",PCSPAISetNBSteps_SPAI));
5859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",PCSPAISetMax_SPAI));
5869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",PCSPAISetMaxNew_SPAI));
5879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",PCSPAISetBlockSize_SPAI));
5889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",PCSPAISetCacheSize_SPAI));
5899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",PCSPAISetVerbose_SPAI));
5909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",PCSPAISetSp_SPAI));
59101a81e61SBarry Smith   PetscFunctionReturn(0);
59201a81e61SBarry Smith }
59301a81e61SBarry Smith 
59401a81e61SBarry Smith /**********************************************************************/
59501a81e61SBarry Smith 
59601a81e61SBarry Smith /*
59701a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
59801a81e61SBarry Smith */
59901a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B)
60001a81e61SBarry Smith {
60101a81e61SBarry Smith   matrix                  *M;
60201a81e61SBarry Smith   int                     i,j,col;
60301a81e61SBarry Smith   int                     row_indx;
60401a81e61SBarry Smith   int                     len,pe,local_indx,start_indx;
60501a81e61SBarry Smith   int                     *mapping;
60601a81e61SBarry Smith   const int               *cols;
60701a81e61SBarry Smith   const double            *vals;
6082122902bSSatish Balay   int                     n,mnl,nnl,nz,rstart,rend;
6092a4967abSBarry Smith   PetscMPIInt             size,rank;
61001a81e61SBarry Smith   struct compressed_lines *rows;
61101a81e61SBarry Smith 
61201a81e61SBarry Smith   PetscFunctionBegin;
6139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
6149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
6159566063dSJacob Faibussowitsch   PetscCall(MatGetSize(A,&n,&n));
6169566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&mnl,&nnl));
61701a81e61SBarry Smith 
61801a81e61SBarry Smith   /*
61901a81e61SBarry Smith     not sure why a barrier is required. commenting out
6209566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Barrier(comm));
62101a81e61SBarry Smith   */
62201a81e61SBarry Smith 
62329b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
62401a81e61SBarry Smith 
62501a81e61SBarry Smith   M->n              = n;
62601a81e61SBarry Smith   M->bs             = 1;
62701a81e61SBarry Smith   M->max_block_size = 1;
62801a81e61SBarry Smith 
62901a81e61SBarry Smith   M->mnls          = (int*)malloc(sizeof(int)*size);
63001a81e61SBarry Smith   M->start_indices = (int*)malloc(sizeof(int)*size);
63101a81e61SBarry Smith   M->pe            = (int*)malloc(sizeof(int)*n);
63201a81e61SBarry Smith   M->block_sizes   = (int*)malloc(sizeof(int)*n);
63301a81e61SBarry Smith   for (i=0; i<n; i++) M->block_sizes[i] = 1;
63401a81e61SBarry Smith 
6359566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm));
63601a81e61SBarry Smith 
63701a81e61SBarry Smith   M->start_indices[0] = 0;
6382fa5cd67SKarl Rupp   for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
63901a81e61SBarry Smith 
64001a81e61SBarry Smith   M->mnl            = M->mnls[M->myid];
64101a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
64201a81e61SBarry Smith 
64301a81e61SBarry Smith   for (i=0; i<size; i++) {
64401a81e61SBarry Smith     start_indx = M->start_indices[i];
6452fa5cd67SKarl Rupp     for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i;
64601a81e61SBarry Smith   }
64701a81e61SBarry Smith 
64801a81e61SBarry Smith   if (AT) {
6492f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],1);
65001a81e61SBarry Smith   } else {
6512f613bf5SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],0);
65201a81e61SBarry Smith   }
65301a81e61SBarry Smith 
65401a81e61SBarry Smith   rows = M->lines;
65501a81e61SBarry Smith 
65601a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
6579566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(M->n,&mapping));
65801a81e61SBarry Smith   pe         = 0;
65901a81e61SBarry Smith   local_indx = 0;
66001a81e61SBarry Smith   for (i=0; i<M->n; i++) {
66101a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
66201a81e61SBarry Smith       pe++;
66301a81e61SBarry Smith       local_indx = 0;
66401a81e61SBarry Smith     }
66501a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
66601a81e61SBarry Smith     local_indx++;
66701a81e61SBarry Smith   }
66801a81e61SBarry Smith 
66901a81e61SBarry Smith   /*********************************************************/
67001a81e61SBarry Smith   /************** Set up the row structure *****************/
67101a81e61SBarry Smith   /*********************************************************/
67201a81e61SBarry Smith 
6739566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A,&rstart,&rend));
67401a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
67501a81e61SBarry Smith     row_indx = i - rstart;
6769566063dSJacob Faibussowitsch     PetscCall(MatGetRow(A,i,&nz,&cols,&vals));
6772122902bSSatish Balay     /* allocate buffers */
6782122902bSSatish Balay     rows->ptrs[row_indx] = (int*)malloc(nz*sizeof(int));
6792122902bSSatish Balay     rows->A[row_indx]    = (double*)malloc(nz*sizeof(double));
6802122902bSSatish Balay     /* copy the matrix */
68101a81e61SBarry Smith     for (j=0; j<nz; j++) {
68201a81e61SBarry Smith       col = cols[j];
68301a81e61SBarry Smith       len = rows->len[row_indx]++;
6842fa5cd67SKarl Rupp 
68501a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
68601a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
68701a81e61SBarry Smith     }
68801a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
6892fa5cd67SKarl Rupp 
6909566063dSJacob Faibussowitsch     PetscCall(MatRestoreRow(A,i,&nz,&cols,&vals));
69101a81e61SBarry Smith   }
69201a81e61SBarry Smith 
69301a81e61SBarry Smith   /************************************************************/
69401a81e61SBarry Smith   /************** Set up the column structure *****************/
69501a81e61SBarry Smith   /*********************************************************/
69601a81e61SBarry Smith 
69701a81e61SBarry Smith   if (AT) {
69801a81e61SBarry Smith 
69901a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
70001a81e61SBarry Smith       row_indx = i - rstart;
7019566063dSJacob Faibussowitsch       PetscCall(MatGetRow(AT,i,&nz,&cols,&vals));
7022122902bSSatish Balay       /* allocate buffers */
7032122902bSSatish Balay       rows->rptrs[row_indx] = (int*)malloc(nz*sizeof(int));
7042122902bSSatish Balay       /* copy the matrix (i.e., the structure) */
70501a81e61SBarry Smith       for (j=0; j<nz; j++) {
70601a81e61SBarry Smith         col = cols[j];
70701a81e61SBarry Smith         len = rows->rlen[row_indx]++;
7082fa5cd67SKarl Rupp 
70901a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
71001a81e61SBarry Smith       }
7119566063dSJacob Faibussowitsch       PetscCall(MatRestoreRow(AT,i,&nz,&cols,&vals));
71201a81e61SBarry Smith     }
71301a81e61SBarry Smith   }
71401a81e61SBarry Smith 
7159566063dSJacob Faibussowitsch   PetscCall(PetscFree(mapping));
71601a81e61SBarry Smith 
71701a81e61SBarry Smith   order_pointers(M);
71801a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
71901a81e61SBarry Smith   *B       = M;
72001a81e61SBarry Smith   PetscFunctionReturn(0);
72101a81e61SBarry Smith }
72201a81e61SBarry Smith 
72301a81e61SBarry Smith /**********************************************************************/
72401a81e61SBarry Smith 
72501a81e61SBarry Smith /*
72601a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
727df3898eeSBarry Smith    This assumes that the SPAI matrix B is stored in
72801a81e61SBarry Smith    COMPRESSED-ROW format.
72901a81e61SBarry Smith */
73001a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB)
73101a81e61SBarry Smith {
7324b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
73301a81e61SBarry Smith   int            m,n,M,N;
73401a81e61SBarry Smith   int            d_nz,o_nz;
73501a81e61SBarry Smith   int            *d_nnz,*o_nnz;
73601a81e61SBarry Smith   int            i,k,global_row,global_col,first_diag_col,last_diag_col;
73701a81e61SBarry Smith   PetscScalar    val;
73801a81e61SBarry Smith 
73901a81e61SBarry Smith   PetscFunctionBegin;
7409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
7419566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
74201a81e61SBarry Smith 
74301a81e61SBarry Smith   m    = n = B->mnls[rank];
74401a81e61SBarry Smith   d_nz = o_nz = 0;
74501a81e61SBarry Smith 
74606946f3aSJose E. Roman   /* Determine preallocation for MatCreateAIJ */
7479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&d_nnz));
7489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(m,&o_nnz));
74901a81e61SBarry Smith   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
75001a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
75101a81e61SBarry Smith   last_diag_col  = first_diag_col + B->mnls[rank];
75201a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
75301a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
75401a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
755db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
756db4deed7SKarl Rupp       else o_nnz[i]++;
75701a81e61SBarry Smith     }
75801a81e61SBarry Smith   }
75901a81e61SBarry Smith 
76001a81e61SBarry Smith   M = N = B->n;
76101a81e61SBarry Smith   /* Here we only know how to create AIJ format */
7629566063dSJacob Faibussowitsch   PetscCall(MatCreate(comm,PB));
7639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*PB,m,n,M,N));
7649566063dSJacob Faibussowitsch   PetscCall(MatSetType(*PB,MATAIJ));
7659566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz));
7669566063dSJacob Faibussowitsch   PetscCall(MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz));
76701a81e61SBarry Smith 
76801a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
76901a81e61SBarry Smith     global_row = B->start_indices[rank]+i;
77001a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
77101a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
7722fa5cd67SKarl Rupp 
77301a81e61SBarry Smith       val  = B->lines->A[i][k];
7749566063dSJacob Faibussowitsch       PetscCall(MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES));
77501a81e61SBarry Smith     }
77601a81e61SBarry Smith   }
77701a81e61SBarry Smith 
7789566063dSJacob Faibussowitsch   PetscCall(PetscFree(d_nnz));
7799566063dSJacob Faibussowitsch   PetscCall(PetscFree(o_nnz));
78001a81e61SBarry Smith 
7819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY));
7829566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY));
78301a81e61SBarry Smith   PetscFunctionReturn(0);
78401a81e61SBarry Smith }
78501a81e61SBarry Smith 
78601a81e61SBarry Smith /**********************************************************************/
78701a81e61SBarry Smith 
78801a81e61SBarry Smith /*
78901a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
79001a81e61SBarry Smith */
79101a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv)
79201a81e61SBarry Smith {
7934b6b5fe1SBarry Smith   PetscMPIInt size,rank;
7944b6b5fe1SBarry Smith   int         m,M,i,*mnls,*start_indices,*global_indices;
79501a81e61SBarry Smith 
79601a81e61SBarry Smith   PetscFunctionBegin;
7979566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
7989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
79901a81e61SBarry Smith 
80001a81e61SBarry Smith   m = v->mnl;
80101a81e61SBarry Smith   M = v->n;
80201a81e61SBarry Smith 
8039566063dSJacob Faibussowitsch   PetscCall(VecCreateMPI(comm,m,M,Pv));
80401a81e61SBarry Smith 
8059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&mnls));
8069566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm));
80701a81e61SBarry Smith 
8089566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(size,&start_indices));
8092fa5cd67SKarl Rupp 
81001a81e61SBarry Smith   start_indices[0] = 0;
8112fa5cd67SKarl Rupp   for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1];
81201a81e61SBarry Smith 
8139566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(v->mnl,&global_indices));
8142fa5cd67SKarl Rupp   for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i;
81501a81e61SBarry Smith 
8169566063dSJacob Faibussowitsch   PetscCall(PetscFree(mnls));
8179566063dSJacob Faibussowitsch   PetscCall(PetscFree(start_indices));
81801a81e61SBarry Smith 
8199566063dSJacob Faibussowitsch   PetscCall(VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES));
8209566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(*Pv));
8219566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(*Pv));
82201a81e61SBarry Smith 
8239566063dSJacob Faibussowitsch   PetscCall(PetscFree(global_indices));
82401a81e61SBarry Smith   PetscFunctionReturn(0);
82501a81e61SBarry Smith }
826