xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision 2a4967ab70d56c5b59dad8fe3ce4d351b67418cd)
101a81e61SBarry Smith #define PETSCKSP_DLL
201a81e61SBarry Smith 
301a81e61SBarry Smith /*
401a81e61SBarry Smith    3/99 Modified by Stephen Barnard to support SPAI version 3.0
501a81e61SBarry Smith */
601a81e61SBarry Smith 
701a81e61SBarry Smith /*
801a81e61SBarry Smith       Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner
901a81e61SBarry Smith    Code written by Stephen Barnard.
1001a81e61SBarry Smith 
1101a81e61SBarry Smith       Note: there is some BAD memory bleeding below!
1201a81e61SBarry Smith 
1301a81e61SBarry Smith       This code needs work
1401a81e61SBarry Smith 
1501a81e61SBarry Smith    1) get rid of all memory bleeding
1601a81e61SBarry Smith    2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix
1701a81e61SBarry Smith       rather than having the sp flag for PC_SPAI
18*2a4967abSBarry Smith    3) fix to set the block size based on the matrix block size
1901a81e61SBarry Smith 
2001a81e61SBarry Smith */
2101a81e61SBarry Smith 
226356e834SBarry Smith #include "private/pcimpl.h"        /*I "petscpc.h" I*/
2375031707SSatish Balay #include "petscspai.h"
2401a81e61SBarry Smith 
2501a81e61SBarry Smith /*
2601a81e61SBarry Smith     These are the SPAI include files
2701a81e61SBarry Smith */
2801a81e61SBarry Smith EXTERN_C_BEGIN
29b7f9c950SSatish Balay #define MPI /* required for setting SPAI_Comm correctly in basics.h */
3001a81e61SBarry Smith #include "spai.h"
3101a81e61SBarry Smith #include "matrix.h"
3201a81e61SBarry Smith EXTERN_C_END
3301a81e61SBarry Smith 
3401a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**);
3501a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix *,Mat *);
3601a81e61SBarry Smith EXTERN PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv);
3701a81e61SBarry 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 #undef __FUNCT__
6201a81e61SBarry Smith #define __FUNCT__ "PCSetUp_SPAI"
6301a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc)
6401a81e61SBarry Smith {
6501a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
6601a81e61SBarry Smith   PetscErrorCode ierr;
6701a81e61SBarry Smith   Mat            AT;
6801a81e61SBarry Smith 
6901a81e61SBarry Smith   PetscFunctionBegin;
7001a81e61SBarry Smith 
7101a81e61SBarry Smith   init_SPAI();
7201a81e61SBarry Smith 
7301a81e61SBarry Smith   if (ispai->sp) {
7401a81e61SBarry Smith     ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr);
7501a81e61SBarry Smith   } else {
7601a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
7701a81e61SBarry Smith     ierr = MatTranspose(pc->pmat,&AT);CHKERRQ(ierr);
7801a81e61SBarry Smith     ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr);
7901a81e61SBarry Smith     ierr = MatDestroy(AT);CHKERRQ(ierr);
8001a81e61SBarry Smith   }
8101a81e61SBarry Smith 
8201a81e61SBarry Smith   /* Destroy the transpose */
8301a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
8401a81e61SBarry Smith 
8501a81e61SBarry Smith   /* construct SPAI preconditioner */
8601a81e61SBarry Smith   /* FILE *messages */     /* file for warning messages */
8701a81e61SBarry Smith   /* double epsilon */     /* tolerance */
8801a81e61SBarry Smith   /* int nbsteps */        /* max number of "improvement" steps per line */
8901a81e61SBarry Smith   /* int max */            /* max dimensions of is_I, q, etc. */
9001a81e61SBarry Smith   /* int maxnew */         /* max number of new entries per step */
9101a81e61SBarry Smith   /* int block_size */     /* block_size == 1 specifies scalar elments
9201a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
9301a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
9401a81e61SBarry Smith   /* int cache_size */     /* one of (1,2,3,4,5,6) indicting size of cache */
9501a81e61SBarry Smith                            /* cache_size == 0 indicates no caching */
9601a81e61SBarry Smith   /* int    verbose    */  /* verbose == 0 specifies that SPAI is silent
9701a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9801a81e61SBarry Smith 
9901a81e61SBarry Smith   ierr = bspai(ispai->B,&ispai->M,
10001a81e61SBarry Smith 		   stdout,
10101a81e61SBarry Smith 		   ispai->epsilon,
10201a81e61SBarry Smith 		   ispai->nbsteps,
10301a81e61SBarry Smith 		   ispai->max,
10401a81e61SBarry Smith 		   ispai->maxnew,
10501a81e61SBarry Smith 		   ispai->block_size,
10601a81e61SBarry Smith 		   ispai->cache_size,
10701a81e61SBarry Smith 	       ispai->verbose); CHKERRQ(ierr);
10801a81e61SBarry Smith 
1097adad957SLisandro Dalcin   ierr = ConvertMatrixToMat(((PetscObject)pc)->comm,ispai->M,&ispai->PM);CHKERRQ(ierr);
11001a81e61SBarry Smith 
11101a81e61SBarry Smith   /* free the SPAI matrices */
11201a81e61SBarry Smith   sp_free_matrix(ispai->B);
11301a81e61SBarry Smith   sp_free_matrix(ispai->M);
11401a81e61SBarry Smith 
11501a81e61SBarry Smith   PetscFunctionReturn(0);
11601a81e61SBarry Smith }
11701a81e61SBarry Smith 
11801a81e61SBarry Smith /**********************************************************************/
11901a81e61SBarry Smith 
12001a81e61SBarry Smith #undef __FUNCT__
12101a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI"
12201a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y)
12301a81e61SBarry Smith {
12401a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
12501a81e61SBarry Smith   PetscErrorCode ierr;
12601a81e61SBarry Smith 
12701a81e61SBarry Smith   PetscFunctionBegin;
12801a81e61SBarry Smith   /* Now using PETSc's multiply */
12901a81e61SBarry Smith   ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr);
13001a81e61SBarry Smith   PetscFunctionReturn(0);
13101a81e61SBarry Smith }
13201a81e61SBarry Smith 
13301a81e61SBarry Smith /**********************************************************************/
13401a81e61SBarry Smith 
13501a81e61SBarry Smith #undef __FUNCT__
13601a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI"
13701a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc)
13801a81e61SBarry Smith {
13901a81e61SBarry Smith   PetscErrorCode ierr;
14001a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
14101a81e61SBarry Smith 
14201a81e61SBarry Smith   PetscFunctionBegin;
14301a81e61SBarry Smith   if (ispai->PM) {ierr = MatDestroy(ispai->PM);CHKERRQ(ierr);}
14401a81e61SBarry Smith   ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr);
14501a81e61SBarry Smith   ierr = PetscFree(ispai);CHKERRQ(ierr);
14601a81e61SBarry Smith   PetscFunctionReturn(0);
14701a81e61SBarry Smith }
14801a81e61SBarry Smith 
14901a81e61SBarry Smith /**********************************************************************/
15001a81e61SBarry Smith 
15101a81e61SBarry Smith #undef __FUNCT__
15201a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI"
15301a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer)
15401a81e61SBarry Smith {
15501a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
15601a81e61SBarry Smith   PetscErrorCode ierr;
15701a81e61SBarry Smith   PetscTruth     iascii;
15801a81e61SBarry Smith 
15901a81e61SBarry Smith   PetscFunctionBegin;
16001a81e61SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
16101a81e61SBarry Smith   if (iascii) {
16201a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    SPAI preconditioner\n");CHKERRQ(ierr);
163a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    epsilon %G\n",   ispai->epsilon);CHKERRQ(ierr);
16401a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps);CHKERRQ(ierr);
16501a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max);CHKERRQ(ierr);
16601a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew);CHKERRQ(ierr);
16701a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size);CHKERRQ(ierr);
16801a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size);CHKERRQ(ierr);
16901a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose);CHKERRQ(ierr);
17001a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp);CHKERRQ(ierr);
17101a81e61SBarry Smith   }
17201a81e61SBarry Smith   PetscFunctionReturn(0);
17301a81e61SBarry Smith }
17401a81e61SBarry Smith 
17501a81e61SBarry Smith EXTERN_C_BEGIN
17601a81e61SBarry Smith #undef __FUNCT__
17701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI"
17801a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon_SPAI(PC pc,double epsilon1)
17901a81e61SBarry Smith {
18001a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
18101a81e61SBarry Smith   PetscFunctionBegin;
18201a81e61SBarry Smith   ispai->epsilon = epsilon1;
18301a81e61SBarry Smith   PetscFunctionReturn(0);
18401a81e61SBarry Smith }
18501a81e61SBarry Smith EXTERN_C_END
18601a81e61SBarry Smith 
18701a81e61SBarry Smith /**********************************************************************/
18801a81e61SBarry Smith 
18901a81e61SBarry Smith EXTERN_C_BEGIN
19001a81e61SBarry Smith #undef __FUNCT__
19101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI"
19201a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1)
19301a81e61SBarry Smith {
19401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
19501a81e61SBarry Smith   PetscFunctionBegin;
19601a81e61SBarry Smith   ispai->nbsteps = nbsteps1;
19701a81e61SBarry Smith   PetscFunctionReturn(0);
19801a81e61SBarry Smith }
19901a81e61SBarry Smith EXTERN_C_END
20001a81e61SBarry Smith 
20101a81e61SBarry Smith /**********************************************************************/
20201a81e61SBarry Smith 
20301a81e61SBarry Smith /* added 1/7/99 g.h. */
20401a81e61SBarry Smith EXTERN_C_BEGIN
20501a81e61SBarry Smith #undef __FUNCT__
20601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax_SPAI"
20701a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax_SPAI(PC pc,int max1)
20801a81e61SBarry Smith {
20901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
21001a81e61SBarry Smith   PetscFunctionBegin;
21101a81e61SBarry Smith   ispai->max = max1;
21201a81e61SBarry Smith   PetscFunctionReturn(0);
21301a81e61SBarry Smith }
21401a81e61SBarry Smith EXTERN_C_END
21501a81e61SBarry Smith 
21601a81e61SBarry Smith /**********************************************************************/
21701a81e61SBarry Smith 
21801a81e61SBarry Smith EXTERN_C_BEGIN
21901a81e61SBarry Smith #undef __FUNCT__
22001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI"
22101a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew_SPAI(PC pc,int maxnew1)
22201a81e61SBarry Smith {
22301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
22401a81e61SBarry Smith   PetscFunctionBegin;
22501a81e61SBarry Smith   ispai->maxnew = maxnew1;
22601a81e61SBarry Smith   PetscFunctionReturn(0);
22701a81e61SBarry Smith }
22801a81e61SBarry Smith EXTERN_C_END
22901a81e61SBarry Smith 
23001a81e61SBarry Smith /**********************************************************************/
23101a81e61SBarry Smith 
23201a81e61SBarry Smith EXTERN_C_BEGIN
23301a81e61SBarry Smith #undef __FUNCT__
23401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize_SPAI"
23501a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize_SPAI(PC pc,int block_size1)
23601a81e61SBarry Smith {
23701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
23801a81e61SBarry Smith   PetscFunctionBegin;
23901a81e61SBarry Smith   ispai->block_size = block_size1;
24001a81e61SBarry Smith   PetscFunctionReturn(0);
24101a81e61SBarry Smith }
24201a81e61SBarry Smith EXTERN_C_END
24301a81e61SBarry Smith 
24401a81e61SBarry Smith /**********************************************************************/
24501a81e61SBarry Smith 
24601a81e61SBarry Smith EXTERN_C_BEGIN
24701a81e61SBarry Smith #undef __FUNCT__
24801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI"
24901a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
25001a81e61SBarry Smith {
25101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
25201a81e61SBarry Smith   PetscFunctionBegin;
25301a81e61SBarry Smith   ispai->cache_size = cache_size;
25401a81e61SBarry Smith   PetscFunctionReturn(0);
25501a81e61SBarry Smith }
25601a81e61SBarry Smith EXTERN_C_END
25701a81e61SBarry Smith 
25801a81e61SBarry Smith /**********************************************************************/
25901a81e61SBarry Smith 
26001a81e61SBarry Smith EXTERN_C_BEGIN
26101a81e61SBarry Smith #undef __FUNCT__
26201a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI"
26301a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose_SPAI(PC pc,int verbose)
26401a81e61SBarry Smith {
26501a81e61SBarry Smith   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
26601a81e61SBarry Smith   PetscFunctionBegin;
26701a81e61SBarry Smith   ispai->verbose = verbose;
26801a81e61SBarry Smith   PetscFunctionReturn(0);
26901a81e61SBarry Smith }
27001a81e61SBarry Smith EXTERN_C_END
27101a81e61SBarry Smith 
27201a81e61SBarry Smith /**********************************************************************/
27301a81e61SBarry Smith 
27401a81e61SBarry Smith EXTERN_C_BEGIN
27501a81e61SBarry Smith #undef __FUNCT__
27601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI"
27701a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp_SPAI(PC pc,int sp)
27801a81e61SBarry Smith {
27901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
28001a81e61SBarry Smith   PetscFunctionBegin;
28101a81e61SBarry Smith   ispai->sp = sp;
28201a81e61SBarry Smith   PetscFunctionReturn(0);
28301a81e61SBarry Smith }
28401a81e61SBarry Smith EXTERN_C_END
28501a81e61SBarry Smith 
28601a81e61SBarry Smith /* -------------------------------------------------------------------*/
28701a81e61SBarry Smith 
28801a81e61SBarry Smith #undef __FUNCT__
28901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon"
29001a81e61SBarry Smith /*@
29101a81e61SBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner
29201a81e61SBarry Smith 
29301a81e61SBarry Smith   Input Parameters:
29401a81e61SBarry Smith + pc - the preconditioner
29501a81e61SBarry Smith - eps - epsilon (default .4)
29601a81e61SBarry Smith 
29701a81e61SBarry Smith   Notes:  Espilon must be between 0 and 1. It controls the
29801a81e61SBarry Smith                  quality of the approximation of M to the inverse of
29901a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
30001a81e61SBarry Smith                  fill, and usually better preconditioners. In many
30101a81e61SBarry Smith                  cases the best choice of epsilon is the one that
30201a81e61SBarry Smith                  divides the total solution time equally between the
30301a81e61SBarry Smith                  preconditioner and the solver.
30401a81e61SBarry Smith 
30501a81e61SBarry Smith   Level: intermediate
30601a81e61SBarry Smith 
30701a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
30801a81e61SBarry Smith   @*/
30901a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC pc,double epsilon1)
31001a81e61SBarry Smith {
31101a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,double);
31201a81e61SBarry Smith   PetscFunctionBegin;
31301a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetEpsilon_C",(void (**)(void))&f);CHKERRQ(ierr);
31401a81e61SBarry Smith   if (f) {
31501a81e61SBarry Smith     ierr = (*f)(pc,epsilon1);CHKERRQ(ierr);
31601a81e61SBarry Smith   }
31701a81e61SBarry Smith   PetscFunctionReturn(0);
31801a81e61SBarry Smith }
31901a81e61SBarry Smith 
32001a81e61SBarry Smith /**********************************************************************/
32101a81e61SBarry Smith 
32201a81e61SBarry Smith #undef __FUNCT__
32301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps"
32401a81e61SBarry Smith /*@
32501a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
32601a81e61SBarry Smith         the SPAI preconditioner
32701a81e61SBarry Smith 
32801a81e61SBarry Smith   Input Parameters:
32901a81e61SBarry Smith + pc - the preconditioner
33001a81e61SBarry Smith - n - number of steps (default 5)
33101a81e61SBarry Smith 
33201a81e61SBarry Smith   Notes:  SPAI constructs to approximation to every column of
33301a81e61SBarry Smith                  the exact inverse of A in a series of improvement
33401a81e61SBarry Smith                  steps. The quality of the approximation is determined
33501a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
33601a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
33701a81e61SBarry Smith                  uses the best approximation constructed so far.
33801a81e61SBarry Smith 
33901a81e61SBarry Smith   Level: intermediate
34001a81e61SBarry Smith 
34101a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew()
34201a81e61SBarry Smith @*/
34301a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC pc,int nbsteps1)
34401a81e61SBarry Smith {
34501a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
34601a81e61SBarry Smith   PetscFunctionBegin;
34701a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetNBSteps_C",(void (**)(void))&f);CHKERRQ(ierr);
34801a81e61SBarry Smith   if (f) {
34901a81e61SBarry Smith     ierr = (*f)(pc,nbsteps1);CHKERRQ(ierr);
35001a81e61SBarry Smith   }
35101a81e61SBarry Smith   PetscFunctionReturn(0);
35201a81e61SBarry Smith }
35301a81e61SBarry Smith 
35401a81e61SBarry Smith /**********************************************************************/
35501a81e61SBarry Smith 
35601a81e61SBarry Smith /* added 1/7/99 g.h. */
35701a81e61SBarry Smith #undef __FUNCT__
35801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax"
35901a81e61SBarry Smith /*@
36001a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
36101a81e61SBarry Smith         the SPAI preconditioner
36201a81e61SBarry Smith 
36301a81e61SBarry Smith   Input Parameters:
36401a81e61SBarry Smith + pc - the preconditioner
36501a81e61SBarry Smith - n - size (default is 5000)
36601a81e61SBarry Smith 
36701a81e61SBarry Smith   Level: intermediate
36801a81e61SBarry Smith 
36901a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
37001a81e61SBarry Smith @*/
37101a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC pc,int max1)
37201a81e61SBarry Smith {
37301a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
37401a81e61SBarry Smith   PetscFunctionBegin;
37501a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMax_C",(void (**)(void))&f);CHKERRQ(ierr);
37601a81e61SBarry Smith   if (f) {
37701a81e61SBarry Smith     ierr = (*f)(pc,max1);CHKERRQ(ierr);
37801a81e61SBarry Smith   }
37901a81e61SBarry Smith   PetscFunctionReturn(0);
38001a81e61SBarry Smith }
38101a81e61SBarry Smith 
38201a81e61SBarry Smith /**********************************************************************/
38301a81e61SBarry Smith 
38401a81e61SBarry Smith #undef __FUNCT__
38501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew"
38601a81e61SBarry Smith /*@
38701a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
38801a81e61SBarry Smith    in SPAI preconditioner
38901a81e61SBarry Smith 
39001a81e61SBarry Smith   Input Parameters:
39101a81e61SBarry Smith + pc - the preconditioner
39201a81e61SBarry Smith - n - maximum number (default 5)
39301a81e61SBarry Smith 
39401a81e61SBarry Smith   Level: intermediate
39501a81e61SBarry Smith 
39601a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps()
39701a81e61SBarry Smith @*/
39801a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC pc,int maxnew1)
39901a81e61SBarry Smith {
40001a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
40101a81e61SBarry Smith   PetscFunctionBegin;
40201a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMaxNew_C",(void (**)(void))&f);CHKERRQ(ierr);
40301a81e61SBarry Smith   if (f) {
40401a81e61SBarry Smith     ierr = (*f)(pc,maxnew1);CHKERRQ(ierr);
40501a81e61SBarry Smith   }
40601a81e61SBarry Smith   PetscFunctionReturn(0);
40701a81e61SBarry Smith }
40801a81e61SBarry Smith 
40901a81e61SBarry Smith /**********************************************************************/
41001a81e61SBarry Smith 
41101a81e61SBarry Smith #undef __FUNCT__
41201a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize"
41301a81e61SBarry Smith /*@
41401a81e61SBarry Smith   PCSPAISetBlockSize - set the block size for the SPAI preconditioner
41501a81e61SBarry Smith 
41601a81e61SBarry Smith   Input Parameters:
41701a81e61SBarry Smith + pc - the preconditioner
41801a81e61SBarry Smith - n - block size (default 1)
41901a81e61SBarry Smith 
42001a81e61SBarry Smith   Notes: A block
42101a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
42201a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
42301a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
42401a81e61SBarry Smith                  variable sized blocks, which are determined by
42501a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
42601a81e61SBarry Smith                  This can be very effective for finite-element
42701a81e61SBarry Smith                  matrices.
42801a81e61SBarry Smith 
42901a81e61SBarry Smith                  SPAI will convert A to block form, use a block
43001a81e61SBarry Smith                  version of the preconditioner algorithm, and then
43101a81e61SBarry Smith                  convert the result back to scalar form.
43201a81e61SBarry Smith 
43301a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
43401a81e61SBarry Smith                  can lead to very significant improvement in
43501a81e61SBarry Smith                  performance.
43601a81e61SBarry Smith 
43701a81e61SBarry Smith 
43801a81e61SBarry Smith   Level: intermediate
43901a81e61SBarry Smith 
44001a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
44101a81e61SBarry Smith @*/
44201a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC pc,int block_size1)
44301a81e61SBarry Smith {
44401a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
44501a81e61SBarry Smith   PetscFunctionBegin;
44601a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr);
44701a81e61SBarry Smith   if (f) {
44801a81e61SBarry Smith     ierr = (*f)(pc,block_size1);CHKERRQ(ierr);
44901a81e61SBarry Smith   }
45001a81e61SBarry Smith   PetscFunctionReturn(0);
45101a81e61SBarry Smith }
45201a81e61SBarry Smith 
45301a81e61SBarry Smith /**********************************************************************/
45401a81e61SBarry Smith 
45501a81e61SBarry Smith #undef __FUNCT__
45601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize"
45701a81e61SBarry Smith /*@
45801a81e61SBarry Smith   PCSPAISetCacheSize - specify cache size in the SPAI preconditioner
45901a81e61SBarry Smith 
46001a81e61SBarry Smith   Input Parameters:
46101a81e61SBarry Smith + pc - the preconditioner
46201a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
46301a81e61SBarry Smith 
46401a81e61SBarry Smith   Notes:    SPAI uses a hash table to cache messages and avoid
46501a81e61SBarry Smith                  redundant communication. If suggest always using
46601a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
46701a81e61SBarry Smith                  version.
46801a81e61SBarry Smith 
46901a81e61SBarry Smith   Level: intermediate
47001a81e61SBarry Smith 
47101a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
47201a81e61SBarry Smith @*/
47301a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC pc,int cache_size)
47401a81e61SBarry Smith {
47501a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
47601a81e61SBarry Smith   PetscFunctionBegin;
47701a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetCacheSize_C",(void (**)(void))&f);CHKERRQ(ierr);
47801a81e61SBarry Smith   if (f) {
47901a81e61SBarry Smith     ierr = (*f)(pc,cache_size);CHKERRQ(ierr);
48001a81e61SBarry Smith   }
48101a81e61SBarry Smith   PetscFunctionReturn(0);
48201a81e61SBarry Smith }
48301a81e61SBarry Smith 
48401a81e61SBarry Smith /**********************************************************************/
48501a81e61SBarry Smith 
48601a81e61SBarry Smith #undef __FUNCT__
48701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose"
48801a81e61SBarry Smith /*@
48901a81e61SBarry Smith   PCSPAISetVerbose - verbosity level for the SPAI preconditioner
49001a81e61SBarry Smith 
49101a81e61SBarry Smith   Input Parameters:
49201a81e61SBarry Smith + pc - the preconditioner
49301a81e61SBarry Smith - n - level (default 1)
49401a81e61SBarry Smith 
49501a81e61SBarry Smith   Notes: print parameters, timings and matrix statistics
49601a81e61SBarry Smith 
49701a81e61SBarry Smith   Level: intermediate
49801a81e61SBarry Smith 
49901a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
50001a81e61SBarry Smith @*/
50101a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC pc,int verbose)
50201a81e61SBarry Smith {
50301a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
50401a81e61SBarry Smith   PetscFunctionBegin;
50501a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetVerbose_C",(void (**)(void))&f);CHKERRQ(ierr);
50601a81e61SBarry Smith   if (f) {
50701a81e61SBarry Smith     ierr = (*f)(pc,verbose);CHKERRQ(ierr);
50801a81e61SBarry Smith   }
50901a81e61SBarry Smith   PetscFunctionReturn(0);
51001a81e61SBarry Smith }
51101a81e61SBarry Smith 
51201a81e61SBarry Smith /**********************************************************************/
51301a81e61SBarry Smith 
51401a81e61SBarry Smith #undef __FUNCT__
51501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp"
51601a81e61SBarry Smith /*@
51701a81e61SBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner
51801a81e61SBarry Smith 
51901a81e61SBarry Smith   Input Parameters:
52001a81e61SBarry Smith + pc - the preconditioner
52101a81e61SBarry Smith - n - 0 or 1
52201a81e61SBarry Smith 
52301a81e61SBarry Smith   Notes: If A has a symmetric nonzero pattern use -sp 1 to
52401a81e61SBarry Smith                  improve performance by eliminating some communication
52501a81e61SBarry Smith                  in the parallel version. Even if A does not have a
52601a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
52701a81e61SBarry Smith                  results, but the code will not follow the published
52801a81e61SBarry Smith                  SPAI algorithm exactly.
52901a81e61SBarry Smith 
53001a81e61SBarry Smith 
53101a81e61SBarry Smith   Level: intermediate
53201a81e61SBarry Smith 
53301a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
53401a81e61SBarry Smith @*/
53501a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC pc,int sp)
53601a81e61SBarry Smith {
53701a81e61SBarry Smith   PetscErrorCode ierr,(*f)(PC,int);
53801a81e61SBarry Smith   PetscFunctionBegin;
53901a81e61SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetSp_C",(void (**)(void))&f);CHKERRQ(ierr);
54001a81e61SBarry Smith   if (f) {
54101a81e61SBarry Smith     ierr = (*f)(pc,sp);CHKERRQ(ierr);
54201a81e61SBarry Smith   }
54301a81e61SBarry Smith   PetscFunctionReturn(0);
54401a81e61SBarry Smith }
54501a81e61SBarry Smith 
54601a81e61SBarry Smith /**********************************************************************/
54701a81e61SBarry Smith 
54801a81e61SBarry Smith /**********************************************************************/
54901a81e61SBarry Smith 
55001a81e61SBarry Smith #undef __FUNCT__
55101a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI"
55201a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc)
55301a81e61SBarry Smith {
55401a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
55501a81e61SBarry Smith   PetscErrorCode ierr;
55601a81e61SBarry Smith   int            nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
55701a81e61SBarry Smith   double         epsilon1;
55801a81e61SBarry Smith   PetscTruth     flg;
55901a81e61SBarry Smith 
56001a81e61SBarry Smith   PetscFunctionBegin;
56101a81e61SBarry Smith   ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr);
56201a81e61SBarry Smith     ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr);
56301a81e61SBarry Smith     if (flg) {
56401a81e61SBarry Smith       ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr);
56501a81e61SBarry Smith     }
56601a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr);
56701a81e61SBarry Smith     if (flg) {
56801a81e61SBarry Smith       ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr);
56901a81e61SBarry Smith     }
57001a81e61SBarry Smith     /* added 1/7/99 g.h. */
57101a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr);
57201a81e61SBarry Smith     if (flg) {
57301a81e61SBarry Smith       ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr);
57401a81e61SBarry Smith     }
57501a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr);
57601a81e61SBarry Smith     if (flg) {
57701a81e61SBarry Smith       ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr);
57801a81e61SBarry Smith     }
57901a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr);
58001a81e61SBarry Smith     if (flg) {
58101a81e61SBarry Smith       ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr);
58201a81e61SBarry Smith     }
58301a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr);
58401a81e61SBarry Smith     if (flg) {
58501a81e61SBarry Smith       ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr);
58601a81e61SBarry Smith     }
58701a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr);
58801a81e61SBarry Smith     if (flg) {
58901a81e61SBarry Smith       ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr);
59001a81e61SBarry Smith     }
59101a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr);
59201a81e61SBarry Smith     if (flg) {
59301a81e61SBarry Smith       ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr);
59401a81e61SBarry Smith     }
59501a81e61SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
59601a81e61SBarry Smith   PetscFunctionReturn(0);
59701a81e61SBarry Smith }
59801a81e61SBarry Smith 
59901a81e61SBarry Smith /**********************************************************************/
60001a81e61SBarry Smith 
60101a81e61SBarry Smith /*MC
60201a81e61SBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard
60301a81e61SBarry Smith      as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)
60401a81e61SBarry Smith 
60501a81e61SBarry Smith    Options Database Keys:
60601a81e61SBarry Smith +  -pc_spai_set_epsilon <eps> - set tolerance
60701a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
60801a81e61SBarry Smith .  -pc_spai_max <m> - set max
60901a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
61001a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
61101a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
61201a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
61301a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
61401a81e61SBarry Smith 
61501a81e61SBarry Smith    Notes: This only works with AIJ matrices.
61601a81e61SBarry Smith 
61701a81e61SBarry Smith    Level: beginner
61801a81e61SBarry Smith 
61901a81e61SBarry Smith    Concepts: approximate inverse
62001a81e61SBarry Smith 
62101a81e61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
62201a81e61SBarry Smith     PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(),
62301a81e61SBarry Smith     PCSPAISetVerbose(), PCSPAISetSp()
62401a81e61SBarry Smith M*/
62501a81e61SBarry Smith 
62601a81e61SBarry Smith EXTERN_C_BEGIN
62701a81e61SBarry Smith #undef __FUNCT__
62801a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI"
62901a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SPAI(PC pc)
63001a81e61SBarry Smith {
63101a81e61SBarry Smith   PC_SPAI        *ispai;
63201a81e61SBarry Smith   PetscErrorCode ierr;
63301a81e61SBarry Smith 
63401a81e61SBarry Smith   PetscFunctionBegin;
63538f2d2fdSLisandro Dalcin   ierr               = PetscNewLog(pc,PC_SPAI,&ispai);CHKERRQ(ierr);
636*2a4967abSBarry Smith   pc->data           = ispai;
63701a81e61SBarry Smith 
63801a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
63901a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
64001a81e61SBarry Smith   pc->ops->applyrichardson = 0;
64101a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
64201a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
64301a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
64401a81e61SBarry Smith 
64501a81e61SBarry Smith   ispai->epsilon    = .4;
64601a81e61SBarry Smith   ispai->nbsteps    = 5;
64701a81e61SBarry Smith   ispai->max        = 5000;
64801a81e61SBarry Smith   ispai->maxnew     = 5;
64901a81e61SBarry Smith   ispai->block_size = 1;
65001a81e61SBarry Smith   ispai->cache_size = 5;
65101a81e61SBarry Smith   ispai->verbose    = 0;
65201a81e61SBarry Smith 
65301a81e61SBarry Smith   ispai->sp         = 1;
6547adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ispai->comm_spai));CHKERRQ(ierr);
65501a81e61SBarry Smith 
65601a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C",
65701a81e61SBarry Smith                     "PCSPAISetEpsilon_SPAI",
65801a81e61SBarry Smith                      PCSPAISetEpsilon_SPAI);CHKERRQ(ierr);
65901a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C",
66001a81e61SBarry Smith                     "PCSPAISetNBSteps_SPAI",
66101a81e61SBarry Smith                      PCSPAISetNBSteps_SPAI);CHKERRQ(ierr);
66201a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C",
66301a81e61SBarry Smith                     "PCSPAISetMax_SPAI",
66401a81e61SBarry Smith                      PCSPAISetMax_SPAI);CHKERRQ(ierr);
66501a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC",
66601a81e61SBarry Smith                     "PCSPAISetMaxNew_SPAI",
66701a81e61SBarry Smith                      PCSPAISetMaxNew_SPAI);CHKERRQ(ierr);
66801a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C",
66901a81e61SBarry Smith                     "PCSPAISetBlockSize_SPAI",
67001a81e61SBarry Smith                      PCSPAISetBlockSize_SPAI);CHKERRQ(ierr);
67101a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C",
67201a81e61SBarry Smith                     "PCSPAISetCacheSize_SPAI",
67301a81e61SBarry Smith                      PCSPAISetCacheSize_SPAI);CHKERRQ(ierr);
67401a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C",
67501a81e61SBarry Smith                     "PCSPAISetVerbose_SPAI",
67601a81e61SBarry Smith                      PCSPAISetVerbose_SPAI);CHKERRQ(ierr);
67701a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C",
67801a81e61SBarry Smith                     "PCSPAISetSp_SPAI",
67901a81e61SBarry Smith                      PCSPAISetSp_SPAI);CHKERRQ(ierr);
68001a81e61SBarry Smith 
68101a81e61SBarry Smith   PetscFunctionReturn(0);
68201a81e61SBarry Smith }
68301a81e61SBarry Smith EXTERN_C_END
68401a81e61SBarry Smith 
68501a81e61SBarry Smith /**********************************************************************/
68601a81e61SBarry Smith 
68701a81e61SBarry Smith /*
68801a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
68901a81e61SBarry Smith */
69001a81e61SBarry Smith #undef __FUNCT__
69101a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix"
69201a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B)
69301a81e61SBarry Smith {
69401a81e61SBarry Smith   matrix                  *M;
69501a81e61SBarry Smith   int                     i,j,col;
69601a81e61SBarry Smith   int                     row_indx;
69701a81e61SBarry Smith   int                     len,pe,local_indx,start_indx;
69801a81e61SBarry Smith   int                     *mapping;
69901a81e61SBarry Smith   PetscErrorCode          ierr;
70001a81e61SBarry Smith   const int               *cols;
70101a81e61SBarry Smith   const double            *vals;
702*2a4967abSBarry Smith   int                     *num_ptr,n,mnl,nnl,nz,rstart,rend;
703*2a4967abSBarry Smith   PetscMPIInt             size,rank;
70401a81e61SBarry Smith   struct compressed_lines *rows;
70501a81e61SBarry Smith 
70601a81e61SBarry Smith   PetscFunctionBegin;
70701a81e61SBarry Smith 
70801a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
70901a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
71001a81e61SBarry Smith   ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr);
71101a81e61SBarry Smith   ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr);
71201a81e61SBarry Smith 
71301a81e61SBarry Smith   /*
71401a81e61SBarry Smith     not sure why a barrier is required. commenting out
71501a81e61SBarry Smith   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
71601a81e61SBarry Smith   */
71701a81e61SBarry Smith 
71829b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
71901a81e61SBarry Smith 
72001a81e61SBarry Smith   M->n = n;
72101a81e61SBarry Smith   M->bs = 1;
72201a81e61SBarry Smith   M->max_block_size = 1;
72301a81e61SBarry Smith 
72401a81e61SBarry Smith   M->mnls          = (int*)malloc(sizeof(int)*size);
72501a81e61SBarry Smith   M->start_indices = (int*)malloc(sizeof(int)*size);
72601a81e61SBarry Smith   M->pe            = (int*)malloc(sizeof(int)*n);
72701a81e61SBarry Smith   M->block_sizes   = (int*)malloc(sizeof(int)*n);
72801a81e61SBarry Smith   for (i=0; i<n; i++) M->block_sizes[i] = 1;
72901a81e61SBarry Smith 
73001a81e61SBarry Smith   ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr);
73101a81e61SBarry Smith 
73201a81e61SBarry Smith   M->start_indices[0] = 0;
73301a81e61SBarry Smith   for (i=1; i<size; i++) {
73401a81e61SBarry Smith     M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
73501a81e61SBarry Smith   }
73601a81e61SBarry Smith 
73701a81e61SBarry Smith   M->mnl = M->mnls[M->myid];
73801a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
73901a81e61SBarry Smith 
74001a81e61SBarry Smith   for (i=0; i<size; i++) {
74101a81e61SBarry Smith     start_indx = M->start_indices[i];
74201a81e61SBarry Smith     for (j=0; j<M->mnls[i]; j++)
74301a81e61SBarry Smith       M->pe[start_indx+j] = i;
74401a81e61SBarry Smith   }
74501a81e61SBarry Smith 
74601a81e61SBarry Smith   if (AT) {
74701a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr);
74801a81e61SBarry Smith   } else {
74901a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr);
75001a81e61SBarry Smith   }
75101a81e61SBarry Smith 
75201a81e61SBarry Smith   rows = M->lines;
75301a81e61SBarry Smith 
75401a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
75501a81e61SBarry Smith   ierr       = PetscMalloc(M->n*sizeof(int),&mapping);CHKERRQ(ierr);
75601a81e61SBarry Smith   pe         = 0;
75701a81e61SBarry Smith   local_indx = 0;
75801a81e61SBarry Smith   for (i=0; i<M->n; i++) {
75901a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
76001a81e61SBarry Smith       pe++;
76101a81e61SBarry Smith       local_indx = 0;
76201a81e61SBarry Smith     }
76301a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
76401a81e61SBarry Smith     local_indx++;
76501a81e61SBarry Smith   }
76601a81e61SBarry Smith 
76701a81e61SBarry Smith 
76801a81e61SBarry Smith   ierr = PetscMalloc(mnl*sizeof(int),&num_ptr);CHKERRQ(ierr);
76901a81e61SBarry Smith 
77001a81e61SBarry Smith   /*********************************************************/
77101a81e61SBarry Smith   /************** Set up the row structure *****************/
77201a81e61SBarry Smith   /*********************************************************/
77301a81e61SBarry Smith 
77401a81e61SBarry Smith   /* count number of nonzeros in every row */
77501a81e61SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
77601a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
77701a81e61SBarry Smith     ierr = MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
77801a81e61SBarry Smith     ierr = MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
77901a81e61SBarry Smith   }
78001a81e61SBarry Smith 
78101a81e61SBarry Smith   /* allocate buffers */
78201a81e61SBarry Smith   len = 0;
78301a81e61SBarry Smith   for (i=0; i<mnl; i++) {
78401a81e61SBarry Smith     if (len < num_ptr[i]) len = num_ptr[i];
78501a81e61SBarry Smith   }
78601a81e61SBarry Smith 
78701a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
78801a81e61SBarry Smith     row_indx             = i-rstart;
78901a81e61SBarry Smith     len                  = num_ptr[row_indx];
79001a81e61SBarry Smith     rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int));
79101a81e61SBarry Smith     rows->A[row_indx]    = (double*)malloc(len*sizeof(double));
79201a81e61SBarry Smith   }
79301a81e61SBarry Smith 
79401a81e61SBarry Smith   /* copy the matrix */
79501a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
79601a81e61SBarry Smith     row_indx = i - rstart;
79701a81e61SBarry Smith     ierr     = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
79801a81e61SBarry Smith     for (j=0; j<nz; j++) {
79901a81e61SBarry Smith       col = cols[j];
80001a81e61SBarry Smith       len = rows->len[row_indx]++;
80101a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
80201a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
80301a81e61SBarry Smith     }
80401a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
80501a81e61SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
80601a81e61SBarry Smith   }
80701a81e61SBarry Smith 
80801a81e61SBarry Smith 
80901a81e61SBarry Smith   /************************************************************/
81001a81e61SBarry Smith   /************** Set up the column structure *****************/
81101a81e61SBarry Smith   /*********************************************************/
81201a81e61SBarry Smith 
81301a81e61SBarry Smith   if (AT) {
81401a81e61SBarry Smith 
81501a81e61SBarry Smith     /* count number of nonzeros in every column */
81601a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
81701a81e61SBarry Smith       ierr = MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
81801a81e61SBarry Smith       ierr = MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
81901a81e61SBarry Smith     }
82001a81e61SBarry Smith 
82101a81e61SBarry Smith     /* allocate buffers */
82201a81e61SBarry Smith     len = 0;
82301a81e61SBarry Smith     for (i=0; i<mnl; i++) {
82401a81e61SBarry Smith       if (len < num_ptr[i]) len = num_ptr[i];
82501a81e61SBarry Smith     }
82601a81e61SBarry Smith 
82701a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
82801a81e61SBarry Smith       row_indx = i-rstart;
82901a81e61SBarry Smith       len      = num_ptr[row_indx];
83001a81e61SBarry Smith       rows->rptrs[row_indx] = (int*)malloc(len*sizeof(int));
83101a81e61SBarry Smith     }
83201a81e61SBarry Smith 
83301a81e61SBarry Smith     /* copy the matrix (i.e., the structure) */
83401a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
83501a81e61SBarry Smith       row_indx = i - rstart;
83601a81e61SBarry Smith       ierr     = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
83701a81e61SBarry Smith       for (j=0; j<nz; j++) {
83801a81e61SBarry Smith 	col = cols[j];
83901a81e61SBarry Smith 	len = rows->rlen[row_indx]++;
84001a81e61SBarry Smith 	rows->rptrs[row_indx][len] = mapping[col];
84101a81e61SBarry Smith       }
84201a81e61SBarry Smith       ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
84301a81e61SBarry Smith     }
84401a81e61SBarry Smith   }
84501a81e61SBarry Smith 
84601a81e61SBarry Smith   ierr = PetscFree(num_ptr);CHKERRQ(ierr);
84701a81e61SBarry Smith   ierr = PetscFree(mapping);CHKERRQ(ierr);
84801a81e61SBarry Smith 
84901a81e61SBarry Smith   order_pointers(M);
85001a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
85101a81e61SBarry Smith 
85201a81e61SBarry Smith   *B = M;
85301a81e61SBarry Smith 
85401a81e61SBarry Smith   PetscFunctionReturn(0);
85501a81e61SBarry Smith }
85601a81e61SBarry Smith 
85701a81e61SBarry Smith /**********************************************************************/
85801a81e61SBarry Smith 
85901a81e61SBarry Smith /*
86001a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
86101a81e61SBarry Smith    This assumes that the the SPAI matrix B is stored in
86201a81e61SBarry Smith    COMPRESSED-ROW format.
86301a81e61SBarry Smith */
86401a81e61SBarry Smith #undef __FUNCT__
86501a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat"
86601a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB)
86701a81e61SBarry Smith {
8684b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
86901a81e61SBarry Smith   PetscErrorCode ierr;
87001a81e61SBarry Smith   int            m,n,M,N;
87101a81e61SBarry Smith   int            d_nz,o_nz;
87201a81e61SBarry Smith   int            *d_nnz,*o_nnz;
87301a81e61SBarry Smith   int            i,k,global_row,global_col,first_diag_col,last_diag_col;
87401a81e61SBarry Smith   PetscScalar    val;
87501a81e61SBarry Smith 
87601a81e61SBarry Smith   PetscFunctionBegin;
87701a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
87801a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
87901a81e61SBarry Smith 
88001a81e61SBarry Smith   m = n = B->mnls[rank];
88101a81e61SBarry Smith   d_nz = o_nz = 0;
88201a81e61SBarry Smith 
88301a81e61SBarry Smith   /* Determine preallocation for MatCreateMPIAIJ */
884*2a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&d_nnz);CHKERRQ(ierr);
885*2a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&o_nnz);CHKERRQ(ierr);
88601a81e61SBarry Smith   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
88701a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
88801a81e61SBarry Smith   last_diag_col = first_diag_col + B->mnls[rank];
88901a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
89001a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
89101a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
89201a81e61SBarry Smith       if ((global_col >= first_diag_col) && (global_col <= last_diag_col))
89301a81e61SBarry Smith 	d_nnz[i]++;
89401a81e61SBarry Smith       else
89501a81e61SBarry Smith 	o_nnz[i]++;
89601a81e61SBarry Smith     }
89701a81e61SBarry Smith   }
89801a81e61SBarry Smith 
89901a81e61SBarry Smith   M = N = B->n;
90001a81e61SBarry Smith   /* Here we only know how to create AIJ format */
901f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,PB);CHKERRQ(ierr);
902f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr);
90301a81e61SBarry Smith   ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr);
90401a81e61SBarry Smith   ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr);
90501a81e61SBarry Smith   ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr);
90601a81e61SBarry Smith 
90701a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
90801a81e61SBarry Smith     global_row = B->start_indices[rank]+i;
90901a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
91001a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
91101a81e61SBarry Smith       val = B->lines->A[i][k];
91201a81e61SBarry Smith       ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr);
91301a81e61SBarry Smith     }
91401a81e61SBarry Smith   }
91501a81e61SBarry Smith 
91601a81e61SBarry Smith   ierr = PetscFree(d_nnz);CHKERRQ(ierr);
91701a81e61SBarry Smith   ierr = PetscFree(o_nnz);CHKERRQ(ierr);
91801a81e61SBarry Smith 
91901a81e61SBarry Smith   ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92001a81e61SBarry Smith   ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92101a81e61SBarry Smith 
92201a81e61SBarry Smith   PetscFunctionReturn(0);
92301a81e61SBarry Smith }
92401a81e61SBarry Smith 
92501a81e61SBarry Smith /**********************************************************************/
92601a81e61SBarry Smith 
92701a81e61SBarry Smith /*
92801a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
92901a81e61SBarry Smith */
93001a81e61SBarry Smith #undef __FUNCT__
93101a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec"
93201a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv)
93301a81e61SBarry Smith {
93401a81e61SBarry Smith   PetscErrorCode ierr;
9354b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
9364b6b5fe1SBarry Smith   int            m,M,i,*mnls,*start_indices,*global_indices;
93701a81e61SBarry Smith 
93801a81e61SBarry Smith   PetscFunctionBegin;
93901a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
94001a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
94101a81e61SBarry Smith 
94201a81e61SBarry Smith   m = v->mnl;
94301a81e61SBarry Smith   M = v->n;
94401a81e61SBarry Smith 
94501a81e61SBarry Smith 
94601a81e61SBarry Smith   ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr);
94701a81e61SBarry Smith 
94801a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&mnls);CHKERRQ(ierr);
949*2a4967abSBarry Smith   ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRQ(ierr);
95001a81e61SBarry Smith 
95101a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&start_indices);CHKERRQ(ierr);
95201a81e61SBarry Smith   start_indices[0] = 0;
95301a81e61SBarry Smith   for (i=1; i<size; i++)
95401a81e61SBarry Smith     start_indices[i] = start_indices[i-1] +mnls[i-1];
95501a81e61SBarry Smith 
95601a81e61SBarry Smith   ierr = PetscMalloc(v->mnl*sizeof(int),&global_indices);CHKERRQ(ierr);
95701a81e61SBarry Smith   for (i=0; i<v->mnl; i++)
95801a81e61SBarry Smith     global_indices[i] = start_indices[rank] + i;
95901a81e61SBarry Smith 
96001a81e61SBarry Smith   ierr = PetscFree(mnls);CHKERRQ(ierr);
96101a81e61SBarry Smith   ierr = PetscFree(start_indices);CHKERRQ(ierr);
96201a81e61SBarry Smith 
96301a81e61SBarry Smith   ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr);
96401a81e61SBarry Smith   ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr);
96501a81e61SBarry Smith   ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr);
96601a81e61SBarry Smith 
9674b6b5fe1SBarry Smith   ierr = PetscFree(global_indices);CHKERRQ(ierr);
96801a81e61SBarry Smith   PetscFunctionReturn(0);
96901a81e61SBarry Smith }
97001a81e61SBarry Smith 
97101a81e61SBarry Smith 
97201a81e61SBarry Smith 
97301a81e61SBarry Smith 
974