xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision 00de8ff0695ff394d09a2c60082aeaab5870b6e2)
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 */
2001a81e61SBarry Smith 
21b45d2f2cSJed Brown #include <petsc-private/pcimpl.h>        /*I "petscpc.h" I*/
2275031707SSatish Balay #include "petscspai.h"
2301a81e61SBarry Smith 
2401a81e61SBarry Smith /*
2501a81e61SBarry Smith     These are the SPAI include files
2601a81e61SBarry Smith */
2701a81e61SBarry Smith EXTERN_C_BEGIN
28b7f9c950SSatish Balay #define MPI /* required for setting SPAI_Comm correctly in basics.h */
29c6db04a5SJed Brown #include <spai.h>
30c6db04a5SJed Brown #include <matrix.h>
3101a81e61SBarry Smith EXTERN_C_END
3201a81e61SBarry Smith 
3309573ac7SBarry Smith extern PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**);
3409573ac7SBarry Smith extern PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix*,Mat*);
3509573ac7SBarry Smith extern PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv);
3609573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char*,char*,char*);
3701a81e61SBarry Smith 
3801a81e61SBarry Smith typedef struct {
3901a81e61SBarry Smith 
4001a81e61SBarry Smith   matrix *B;                /* matrix in SPAI format */
4101a81e61SBarry Smith   matrix *BT;               /* transpose of matrix in SPAI format */
4201a81e61SBarry Smith   matrix *M;                /* the approximate inverse in SPAI format */
4301a81e61SBarry Smith 
4401a81e61SBarry Smith   Mat PM;                   /* the approximate inverse PETSc format */
4501a81e61SBarry Smith 
4601a81e61SBarry Smith   double epsilon;           /* tolerance */
4701a81e61SBarry Smith   int    nbsteps;           /* max number of "improvement" steps per line */
4801a81e61SBarry Smith   int    max;               /* max dimensions of is_I, q, etc. */
4901a81e61SBarry Smith   int    maxnew;            /* max number of new entries per step */
5001a81e61SBarry Smith   int    block_size;        /* constant block size */
5101a81e61SBarry Smith   int    cache_size;        /* one of (1,2,3,4,5,6) indicting size of cache */
5201a81e61SBarry Smith   int    verbose;           /* SPAI prints timing and statistics */
5301a81e61SBarry Smith 
5401a81e61SBarry Smith   int      sp;              /* symmetric nonzero pattern */
5501a81e61SBarry Smith   MPI_Comm comm_spai;     /* communicator to be used with spai */
5601a81e61SBarry Smith } PC_SPAI;
5701a81e61SBarry Smith 
5801a81e61SBarry Smith /**********************************************************************/
5901a81e61SBarry Smith 
6001a81e61SBarry Smith #undef __FUNCT__
6101a81e61SBarry Smith #define __FUNCT__ "PCSetUp_SPAI"
6201a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc)
6301a81e61SBarry Smith {
6401a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
6501a81e61SBarry Smith   PetscErrorCode ierr;
6601a81e61SBarry Smith   Mat            AT;
6701a81e61SBarry Smith 
6801a81e61SBarry Smith   PetscFunctionBegin;
6901a81e61SBarry Smith   init_SPAI();
7001a81e61SBarry Smith 
7101a81e61SBarry Smith   if (ispai->sp) {
7201a81e61SBarry Smith     ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr);
7301a81e61SBarry Smith   } else {
7401a81e61SBarry Smith     /* Use the transpose to get the column nonzero structure. */
75fc4dec0aSBarry Smith     ierr = MatTranspose(pc->pmat,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr);
7601a81e61SBarry Smith     ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr);
7769a3b875SJed Brown     ierr = MatDestroy(&AT);CHKERRQ(ierr);
7801a81e61SBarry Smith   }
7901a81e61SBarry Smith 
8001a81e61SBarry Smith   /* Destroy the transpose */
8101a81e61SBarry Smith   /* Don't know how to do it. PETSc developers? */
8201a81e61SBarry Smith 
8301a81e61SBarry Smith   /* construct SPAI preconditioner */
8401a81e61SBarry Smith   /* FILE *messages */     /* file for warning messages */
8501a81e61SBarry Smith   /* double epsilon */     /* tolerance */
8601a81e61SBarry Smith   /* int nbsteps */        /* max number of "improvement" steps per line */
8701a81e61SBarry Smith   /* int max */            /* max dimensions of is_I, q, etc. */
8801a81e61SBarry Smith   /* int maxnew */         /* max number of new entries per step */
8901a81e61SBarry Smith   /* int block_size */     /* block_size == 1 specifies scalar elments
9001a81e61SBarry Smith                               block_size == n specifies nxn constant-block elements
9101a81e61SBarry Smith                               block_size == 0 specifies variable-block elements */
922fa5cd67SKarl Rupp   /* int cache_size */     /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */
9301a81e61SBarry Smith   /* int    verbose    */  /* verbose == 0 specifies that SPAI is silent
9401a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9501a81e61SBarry Smith 
9601a81e61SBarry Smith   ierr = bspai(ispai->B,&ispai->M,
9701a81e61SBarry Smith                stdout,
9801a81e61SBarry Smith                ispai->epsilon,
9901a81e61SBarry Smith                ispai->nbsteps,
10001a81e61SBarry Smith                ispai->max,
10101a81e61SBarry Smith                ispai->maxnew,
10201a81e61SBarry Smith                ispai->block_size,
10301a81e61SBarry Smith                ispai->cache_size,
10401a81e61SBarry Smith                ispai->verbose);CHKERRQ(ierr);
10501a81e61SBarry Smith 
106ce94432eSBarry Smith   ierr = ConvertMatrixToMat(PetscObjectComm((PetscObject)pc),ispai->M,&ispai->PM);CHKERRQ(ierr);
10701a81e61SBarry Smith 
10801a81e61SBarry Smith   /* free the SPAI matrices */
10901a81e61SBarry Smith   sp_free_matrix(ispai->B);
11001a81e61SBarry Smith   sp_free_matrix(ispai->M);
11101a81e61SBarry Smith   PetscFunctionReturn(0);
11201a81e61SBarry Smith }
11301a81e61SBarry Smith 
11401a81e61SBarry Smith /**********************************************************************/
11501a81e61SBarry Smith 
11601a81e61SBarry Smith #undef __FUNCT__
11701a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI"
11801a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y)
11901a81e61SBarry Smith {
12001a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
12101a81e61SBarry Smith   PetscErrorCode ierr;
12201a81e61SBarry Smith 
12301a81e61SBarry Smith   PetscFunctionBegin;
12401a81e61SBarry Smith   /* Now using PETSc's multiply */
12501a81e61SBarry Smith   ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr);
12601a81e61SBarry Smith   PetscFunctionReturn(0);
12701a81e61SBarry Smith }
12801a81e61SBarry Smith 
12901a81e61SBarry Smith /**********************************************************************/
13001a81e61SBarry Smith 
13101a81e61SBarry Smith #undef __FUNCT__
13201a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI"
13301a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc)
13401a81e61SBarry Smith {
13501a81e61SBarry Smith   PetscErrorCode ierr;
13601a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
13701a81e61SBarry Smith 
13801a81e61SBarry Smith   PetscFunctionBegin;
13969a3b875SJed Brown   ierr = MatDestroy(&ispai->PM);CHKERRQ(ierr);
14001a81e61SBarry Smith   ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr);
141c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
14201a81e61SBarry Smith   PetscFunctionReturn(0);
14301a81e61SBarry Smith }
14401a81e61SBarry Smith 
14501a81e61SBarry Smith /**********************************************************************/
14601a81e61SBarry Smith 
14701a81e61SBarry Smith #undef __FUNCT__
14801a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI"
14901a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer)
15001a81e61SBarry Smith {
15101a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
15201a81e61SBarry Smith   PetscErrorCode ierr;
153ace3abfcSBarry Smith   PetscBool      iascii;
15401a81e61SBarry Smith 
15501a81e61SBarry Smith   PetscFunctionBegin;
156251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
15701a81e61SBarry Smith   if (iascii) {
15801a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    SPAI preconditioner\n");CHKERRQ(ierr);
159a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    epsilon %G\n",   ispai->epsilon);CHKERRQ(ierr);
16001a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps);CHKERRQ(ierr);
16101a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max);CHKERRQ(ierr);
16201a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew);CHKERRQ(ierr);
16301a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size);CHKERRQ(ierr);
16401a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size);CHKERRQ(ierr);
16501a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose);CHKERRQ(ierr);
16601a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp);CHKERRQ(ierr);
16701a81e61SBarry Smith   }
16801a81e61SBarry Smith   PetscFunctionReturn(0);
16901a81e61SBarry Smith }
17001a81e61SBarry Smith 
17101a81e61SBarry Smith EXTERN_C_BEGIN
17201a81e61SBarry Smith #undef __FUNCT__
17301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI"
1747087cfbeSBarry Smith PetscErrorCode  PCSPAISetEpsilon_SPAI(PC pc,double epsilon1)
17501a81e61SBarry Smith {
17601a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
1775fd66863SKarl Rupp 
17801a81e61SBarry Smith   PetscFunctionBegin;
17901a81e61SBarry Smith   ispai->epsilon = epsilon1;
18001a81e61SBarry Smith   PetscFunctionReturn(0);
18101a81e61SBarry Smith }
18201a81e61SBarry Smith EXTERN_C_END
18301a81e61SBarry Smith 
18401a81e61SBarry Smith /**********************************************************************/
18501a81e61SBarry Smith 
18601a81e61SBarry Smith EXTERN_C_BEGIN
18701a81e61SBarry Smith #undef __FUNCT__
18801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI"
1897087cfbeSBarry Smith PetscErrorCode  PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1)
19001a81e61SBarry Smith {
19101a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
1925fd66863SKarl Rupp 
19301a81e61SBarry Smith   PetscFunctionBegin;
19401a81e61SBarry Smith   ispai->nbsteps = nbsteps1;
19501a81e61SBarry Smith   PetscFunctionReturn(0);
19601a81e61SBarry Smith }
19701a81e61SBarry Smith EXTERN_C_END
19801a81e61SBarry Smith 
19901a81e61SBarry Smith /**********************************************************************/
20001a81e61SBarry Smith 
20101a81e61SBarry Smith /* added 1/7/99 g.h. */
20201a81e61SBarry Smith EXTERN_C_BEGIN
20301a81e61SBarry Smith #undef __FUNCT__
20401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax_SPAI"
2057087cfbeSBarry Smith PetscErrorCode  PCSPAISetMax_SPAI(PC pc,int max1)
20601a81e61SBarry Smith {
20701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2085fd66863SKarl Rupp 
20901a81e61SBarry Smith   PetscFunctionBegin;
21001a81e61SBarry Smith   ispai->max = max1;
21101a81e61SBarry Smith   PetscFunctionReturn(0);
21201a81e61SBarry Smith }
21301a81e61SBarry Smith EXTERN_C_END
21401a81e61SBarry Smith 
21501a81e61SBarry Smith /**********************************************************************/
21601a81e61SBarry Smith 
21701a81e61SBarry Smith EXTERN_C_BEGIN
21801a81e61SBarry Smith #undef __FUNCT__
21901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI"
2207087cfbeSBarry Smith PetscErrorCode  PCSPAISetMaxNew_SPAI(PC pc,int maxnew1)
22101a81e61SBarry Smith {
22201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2235fd66863SKarl Rupp 
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"
2357087cfbeSBarry Smith PetscErrorCode  PCSPAISetBlockSize_SPAI(PC pc,int block_size1)
23601a81e61SBarry Smith {
23701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2385fd66863SKarl Rupp 
23901a81e61SBarry Smith   PetscFunctionBegin;
24001a81e61SBarry Smith   ispai->block_size = block_size1;
24101a81e61SBarry Smith   PetscFunctionReturn(0);
24201a81e61SBarry Smith }
24301a81e61SBarry Smith EXTERN_C_END
24401a81e61SBarry Smith 
24501a81e61SBarry Smith /**********************************************************************/
24601a81e61SBarry Smith 
24701a81e61SBarry Smith EXTERN_C_BEGIN
24801a81e61SBarry Smith #undef __FUNCT__
24901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI"
2507087cfbeSBarry Smith PetscErrorCode  PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
25101a81e61SBarry Smith {
25201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2535fd66863SKarl Rupp 
25401a81e61SBarry Smith   PetscFunctionBegin;
25501a81e61SBarry Smith   ispai->cache_size = cache_size;
25601a81e61SBarry Smith   PetscFunctionReturn(0);
25701a81e61SBarry Smith }
25801a81e61SBarry Smith EXTERN_C_END
25901a81e61SBarry Smith 
26001a81e61SBarry Smith /**********************************************************************/
26101a81e61SBarry Smith 
26201a81e61SBarry Smith EXTERN_C_BEGIN
26301a81e61SBarry Smith #undef __FUNCT__
26401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI"
2657087cfbeSBarry Smith PetscErrorCode  PCSPAISetVerbose_SPAI(PC pc,int verbose)
26601a81e61SBarry Smith {
26701a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2685fd66863SKarl Rupp 
26901a81e61SBarry Smith   PetscFunctionBegin;
27001a81e61SBarry Smith   ispai->verbose = verbose;
27101a81e61SBarry Smith   PetscFunctionReturn(0);
27201a81e61SBarry Smith }
27301a81e61SBarry Smith EXTERN_C_END
27401a81e61SBarry Smith 
27501a81e61SBarry Smith /**********************************************************************/
27601a81e61SBarry Smith 
27701a81e61SBarry Smith EXTERN_C_BEGIN
27801a81e61SBarry Smith #undef __FUNCT__
27901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI"
2807087cfbeSBarry Smith PetscErrorCode  PCSPAISetSp_SPAI(PC pc,int sp)
28101a81e61SBarry Smith {
28201a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
2835fd66863SKarl Rupp 
28401a81e61SBarry Smith   PetscFunctionBegin;
28501a81e61SBarry Smith   ispai->sp = sp;
28601a81e61SBarry Smith   PetscFunctionReturn(0);
28701a81e61SBarry Smith }
28801a81e61SBarry Smith EXTERN_C_END
28901a81e61SBarry Smith 
29001a81e61SBarry Smith /* -------------------------------------------------------------------*/
29101a81e61SBarry Smith 
29201a81e61SBarry Smith #undef __FUNCT__
29301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon"
29401a81e61SBarry Smith /*@
29501a81e61SBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner
29601a81e61SBarry Smith 
29701a81e61SBarry Smith   Input Parameters:
29801a81e61SBarry Smith + pc - the preconditioner
29901a81e61SBarry Smith - eps - epsilon (default .4)
30001a81e61SBarry Smith 
30101a81e61SBarry Smith   Notes:  Espilon must be between 0 and 1. It controls the
30201a81e61SBarry Smith                  quality of the approximation of M to the inverse of
30301a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
30401a81e61SBarry Smith                  fill, and usually better preconditioners. In many
30501a81e61SBarry Smith                  cases the best choice of epsilon is the one that
30601a81e61SBarry Smith                  divides the total solution time equally between the
30701a81e61SBarry Smith                  preconditioner and the solver.
30801a81e61SBarry Smith 
30901a81e61SBarry Smith   Level: intermediate
31001a81e61SBarry Smith 
31101a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
31201a81e61SBarry Smith   @*/
3137087cfbeSBarry Smith PetscErrorCode  PCSPAISetEpsilon(PC pc,double epsilon1)
31401a81e61SBarry Smith {
3154ac538c5SBarry Smith   PetscErrorCode ierr;
3165fd66863SKarl Rupp 
31701a81e61SBarry Smith   PetscFunctionBegin;
3184ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));CHKERRQ(ierr);
31901a81e61SBarry Smith   PetscFunctionReturn(0);
32001a81e61SBarry Smith }
32101a81e61SBarry Smith 
32201a81e61SBarry Smith /**********************************************************************/
32301a81e61SBarry Smith 
32401a81e61SBarry Smith #undef __FUNCT__
32501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps"
32601a81e61SBarry Smith /*@
32701a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
32801a81e61SBarry Smith         the SPAI preconditioner
32901a81e61SBarry Smith 
33001a81e61SBarry Smith   Input Parameters:
33101a81e61SBarry Smith + pc - the preconditioner
33201a81e61SBarry Smith - n - number of steps (default 5)
33301a81e61SBarry Smith 
33401a81e61SBarry Smith   Notes:  SPAI constructs to approximation to every column of
33501a81e61SBarry Smith                  the exact inverse of A in a series of improvement
33601a81e61SBarry Smith                  steps. The quality of the approximation is determined
33701a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
33801a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
33901a81e61SBarry Smith                  uses the best approximation constructed so far.
34001a81e61SBarry Smith 
34101a81e61SBarry Smith   Level: intermediate
34201a81e61SBarry Smith 
34301a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew()
34401a81e61SBarry Smith @*/
3457087cfbeSBarry Smith PetscErrorCode  PCSPAISetNBSteps(PC pc,int nbsteps1)
34601a81e61SBarry Smith {
3474ac538c5SBarry Smith   PetscErrorCode ierr;
3485fd66863SKarl Rupp 
34901a81e61SBarry Smith   PetscFunctionBegin;
3504ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));CHKERRQ(ierr);
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 @*/
3717087cfbeSBarry Smith PetscErrorCode  PCSPAISetMax(PC pc,int max1)
37201a81e61SBarry Smith {
3734ac538c5SBarry Smith   PetscErrorCode ierr;
3745fd66863SKarl Rupp 
37501a81e61SBarry Smith   PetscFunctionBegin;
3764ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));CHKERRQ(ierr);
37701a81e61SBarry Smith   PetscFunctionReturn(0);
37801a81e61SBarry Smith }
37901a81e61SBarry Smith 
38001a81e61SBarry Smith /**********************************************************************/
38101a81e61SBarry Smith 
38201a81e61SBarry Smith #undef __FUNCT__
38301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew"
38401a81e61SBarry Smith /*@
38501a81e61SBarry Smith   PCSPAISetMaxNew - set maximum number of new nonzero candidates per step
38601a81e61SBarry Smith    in SPAI preconditioner
38701a81e61SBarry Smith 
38801a81e61SBarry Smith   Input Parameters:
38901a81e61SBarry Smith + pc - the preconditioner
39001a81e61SBarry Smith - n - maximum number (default 5)
39101a81e61SBarry Smith 
39201a81e61SBarry Smith   Level: intermediate
39301a81e61SBarry Smith 
39401a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps()
39501a81e61SBarry Smith @*/
3967087cfbeSBarry Smith PetscErrorCode  PCSPAISetMaxNew(PC pc,int maxnew1)
39701a81e61SBarry Smith {
3984ac538c5SBarry Smith   PetscErrorCode ierr;
3995fd66863SKarl Rupp 
40001a81e61SBarry Smith   PetscFunctionBegin;
4014ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));CHKERRQ(ierr);
40201a81e61SBarry Smith   PetscFunctionReturn(0);
40301a81e61SBarry Smith }
40401a81e61SBarry Smith 
40501a81e61SBarry Smith /**********************************************************************/
40601a81e61SBarry Smith 
40701a81e61SBarry Smith #undef __FUNCT__
40801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize"
40901a81e61SBarry Smith /*@
41001a81e61SBarry Smith   PCSPAISetBlockSize - set the block size for the SPAI preconditioner
41101a81e61SBarry Smith 
41201a81e61SBarry Smith   Input Parameters:
41301a81e61SBarry Smith + pc - the preconditioner
41401a81e61SBarry Smith - n - block size (default 1)
41501a81e61SBarry Smith 
41601a81e61SBarry Smith   Notes: A block
41701a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
41801a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
41901a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
42001a81e61SBarry Smith                  variable sized blocks, which are determined by
42101a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
42201a81e61SBarry Smith                  This can be very effective for finite-element
42301a81e61SBarry Smith                  matrices.
42401a81e61SBarry Smith 
42501a81e61SBarry Smith                  SPAI will convert A to block form, use a block
42601a81e61SBarry Smith                  version of the preconditioner algorithm, and then
42701a81e61SBarry Smith                  convert the result back to scalar form.
42801a81e61SBarry Smith 
42901a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
43001a81e61SBarry Smith                  can lead to very significant improvement in
43101a81e61SBarry Smith                  performance.
43201a81e61SBarry Smith 
43301a81e61SBarry Smith 
43401a81e61SBarry Smith   Level: intermediate
43501a81e61SBarry Smith 
43601a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
43701a81e61SBarry Smith @*/
4387087cfbeSBarry Smith PetscErrorCode  PCSPAISetBlockSize(PC pc,int block_size1)
43901a81e61SBarry Smith {
4404ac538c5SBarry Smith   PetscErrorCode ierr;
4415fd66863SKarl Rupp 
44201a81e61SBarry Smith   PetscFunctionBegin;
4434ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));CHKERRQ(ierr);
44401a81e61SBarry Smith   PetscFunctionReturn(0);
44501a81e61SBarry Smith }
44601a81e61SBarry Smith 
44701a81e61SBarry Smith /**********************************************************************/
44801a81e61SBarry Smith 
44901a81e61SBarry Smith #undef __FUNCT__
45001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize"
45101a81e61SBarry Smith /*@
45201a81e61SBarry Smith   PCSPAISetCacheSize - specify cache size in the SPAI preconditioner
45301a81e61SBarry Smith 
45401a81e61SBarry Smith   Input Parameters:
45501a81e61SBarry Smith + pc - the preconditioner
45601a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
45701a81e61SBarry Smith 
45801a81e61SBarry Smith   Notes:    SPAI uses a hash table to cache messages and avoid
45901a81e61SBarry Smith                  redundant communication. If suggest always using
46001a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
46101a81e61SBarry Smith                  version.
46201a81e61SBarry Smith 
46301a81e61SBarry Smith   Level: intermediate
46401a81e61SBarry Smith 
46501a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
46601a81e61SBarry Smith @*/
4677087cfbeSBarry Smith PetscErrorCode  PCSPAISetCacheSize(PC pc,int cache_size)
46801a81e61SBarry Smith {
4694ac538c5SBarry Smith   PetscErrorCode ierr;
4705fd66863SKarl Rupp 
47101a81e61SBarry Smith   PetscFunctionBegin;
4724ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));CHKERRQ(ierr);
47301a81e61SBarry Smith   PetscFunctionReturn(0);
47401a81e61SBarry Smith }
47501a81e61SBarry Smith 
47601a81e61SBarry Smith /**********************************************************************/
47701a81e61SBarry Smith 
47801a81e61SBarry Smith #undef __FUNCT__
47901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose"
48001a81e61SBarry Smith /*@
48101a81e61SBarry Smith   PCSPAISetVerbose - verbosity level for the SPAI preconditioner
48201a81e61SBarry Smith 
48301a81e61SBarry Smith   Input Parameters:
48401a81e61SBarry Smith + pc - the preconditioner
48501a81e61SBarry Smith - n - level (default 1)
48601a81e61SBarry Smith 
48701a81e61SBarry Smith   Notes: print parameters, timings and matrix statistics
48801a81e61SBarry Smith 
48901a81e61SBarry Smith   Level: intermediate
49001a81e61SBarry Smith 
49101a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
49201a81e61SBarry Smith @*/
4937087cfbeSBarry Smith PetscErrorCode  PCSPAISetVerbose(PC pc,int verbose)
49401a81e61SBarry Smith {
4954ac538c5SBarry Smith   PetscErrorCode ierr;
4965fd66863SKarl Rupp 
49701a81e61SBarry Smith   PetscFunctionBegin;
4984ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));CHKERRQ(ierr);
49901a81e61SBarry Smith   PetscFunctionReturn(0);
50001a81e61SBarry Smith }
50101a81e61SBarry Smith 
50201a81e61SBarry Smith /**********************************************************************/
50301a81e61SBarry Smith 
50401a81e61SBarry Smith #undef __FUNCT__
50501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp"
50601a81e61SBarry Smith /*@
50701a81e61SBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner
50801a81e61SBarry Smith 
50901a81e61SBarry Smith   Input Parameters:
51001a81e61SBarry Smith + pc - the preconditioner
51101a81e61SBarry Smith - n - 0 or 1
51201a81e61SBarry Smith 
51301a81e61SBarry Smith   Notes: If A has a symmetric nonzero pattern use -sp 1 to
51401a81e61SBarry Smith                  improve performance by eliminating some communication
51501a81e61SBarry Smith                  in the parallel version. Even if A does not have a
51601a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
51701a81e61SBarry Smith                  results, but the code will not follow the published
51801a81e61SBarry Smith                  SPAI algorithm exactly.
51901a81e61SBarry Smith 
52001a81e61SBarry Smith 
52101a81e61SBarry Smith   Level: intermediate
52201a81e61SBarry Smith 
52301a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
52401a81e61SBarry Smith @*/
5257087cfbeSBarry Smith PetscErrorCode  PCSPAISetSp(PC pc,int sp)
52601a81e61SBarry Smith {
5274ac538c5SBarry Smith   PetscErrorCode ierr;
5285fd66863SKarl Rupp 
52901a81e61SBarry Smith   PetscFunctionBegin;
5304ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));CHKERRQ(ierr);
53101a81e61SBarry Smith   PetscFunctionReturn(0);
53201a81e61SBarry Smith }
53301a81e61SBarry Smith 
53401a81e61SBarry Smith /**********************************************************************/
53501a81e61SBarry Smith 
53601a81e61SBarry Smith /**********************************************************************/
53701a81e61SBarry Smith 
53801a81e61SBarry Smith #undef __FUNCT__
53901a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI"
54001a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc)
54101a81e61SBarry Smith {
54201a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
54301a81e61SBarry Smith   PetscErrorCode ierr;
54401a81e61SBarry Smith   int            nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
54501a81e61SBarry Smith   double         epsilon1;
546ace3abfcSBarry Smith   PetscBool      flg;
54701a81e61SBarry Smith 
54801a81e61SBarry Smith   PetscFunctionBegin;
54901a81e61SBarry Smith   ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr);
55001a81e61SBarry Smith   ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr);
55101a81e61SBarry Smith   if (flg) {
55201a81e61SBarry Smith     ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr);
55301a81e61SBarry Smith   }
55401a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr);
55501a81e61SBarry Smith   if (flg) {
55601a81e61SBarry Smith     ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr);
55701a81e61SBarry Smith   }
55801a81e61SBarry Smith   /* added 1/7/99 g.h. */
55901a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr);
56001a81e61SBarry Smith   if (flg) {
56101a81e61SBarry Smith     ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr);
56201a81e61SBarry Smith   }
56301a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr);
56401a81e61SBarry Smith   if (flg) {
56501a81e61SBarry Smith     ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr);
56601a81e61SBarry Smith   }
56701a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr);
56801a81e61SBarry Smith   if (flg) {
56901a81e61SBarry Smith     ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr);
57001a81e61SBarry Smith   }
57101a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr);
57201a81e61SBarry Smith   if (flg) {
57301a81e61SBarry Smith     ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr);
57401a81e61SBarry Smith   }
57501a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr);
57601a81e61SBarry Smith   if (flg) {
57701a81e61SBarry Smith     ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr);
57801a81e61SBarry Smith   }
57901a81e61SBarry Smith   ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr);
58001a81e61SBarry Smith   if (flg) {
58101a81e61SBarry Smith     ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr);
58201a81e61SBarry Smith   }
58301a81e61SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
58401a81e61SBarry Smith   PetscFunctionReturn(0);
58501a81e61SBarry Smith }
58601a81e61SBarry Smith 
58701a81e61SBarry Smith /**********************************************************************/
58801a81e61SBarry Smith 
58901a81e61SBarry Smith /*MC
59001a81e61SBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard
59101a81e61SBarry Smith      as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)
59201a81e61SBarry Smith 
59301a81e61SBarry Smith    Options Database Keys:
594c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
59501a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
59601a81e61SBarry Smith .  -pc_spai_max <m> - set max
59701a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
59801a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
59901a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
60001a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
60101a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
60201a81e61SBarry Smith 
60301a81e61SBarry Smith    Notes: This only works with AIJ matrices.
60401a81e61SBarry Smith 
60501a81e61SBarry Smith    Level: beginner
60601a81e61SBarry Smith 
60701a81e61SBarry Smith    Concepts: approximate inverse
60801a81e61SBarry Smith 
60901a81e61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
61001a81e61SBarry Smith     PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(),
61101a81e61SBarry Smith     PCSPAISetVerbose(), PCSPAISetSp()
61201a81e61SBarry Smith M*/
61301a81e61SBarry Smith 
61401a81e61SBarry Smith EXTERN_C_BEGIN
61501a81e61SBarry Smith #undef __FUNCT__
61601a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI"
6177087cfbeSBarry Smith PetscErrorCode  PCCreate_SPAI(PC pc)
61801a81e61SBarry Smith {
61901a81e61SBarry Smith   PC_SPAI        *ispai;
62001a81e61SBarry Smith   PetscErrorCode ierr;
62101a81e61SBarry Smith 
62201a81e61SBarry Smith   PetscFunctionBegin;
62338f2d2fdSLisandro Dalcin   ierr     = PetscNewLog(pc,PC_SPAI,&ispai);CHKERRQ(ierr);
6242a4967abSBarry Smith   pc->data = ispai;
62501a81e61SBarry Smith 
62601a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
62701a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
62801a81e61SBarry Smith   pc->ops->applyrichardson = 0;
62901a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
63001a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
63101a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
63201a81e61SBarry Smith 
63301a81e61SBarry Smith   ispai->epsilon    = .4;
63401a81e61SBarry Smith   ispai->nbsteps    = 5;
63501a81e61SBarry Smith   ispai->max        = 5000;
63601a81e61SBarry Smith   ispai->maxnew     = 5;
63701a81e61SBarry Smith   ispai->block_size = 1;
63801a81e61SBarry Smith   ispai->cache_size = 5;
63901a81e61SBarry Smith   ispai->verbose    = 0;
64001a81e61SBarry Smith 
64101a81e61SBarry Smith   ispai->sp = 1;
642ce94432eSBarry Smith   ierr      = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai));CHKERRQ(ierr);
64301a81e61SBarry Smith 
644*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C","PCSPAISetEpsilon_SPAI",PCSPAISetEpsilon_SPAI);CHKERRQ(ierr);
645*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C","PCSPAISetNBSteps_SPAI",PCSPAISetNBSteps_SPAI);CHKERRQ(ierr);
646*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C","PCSPAISetMax_SPAI",PCSPAISetMax_SPAI);CHKERRQ(ierr);
647*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C","PCSPAISetMaxNew_SPAI",PCSPAISetMaxNew_SPAI);CHKERRQ(ierr);
648*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C","PCSPAISetBlockSize_SPAI",PCSPAISetBlockSize_SPAI);CHKERRQ(ierr);
649*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C","PCSPAISetCacheSize_SPAI",PCSPAISetCacheSize_SPAI);CHKERRQ(ierr);
650*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C","PCSPAISetVerbose_SPAI",PCSPAISetVerbose_SPAI);CHKERRQ(ierr);
651*00de8ff0SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C","PCSPAISetSp_SPAI",PCSPAISetSp_SPAI);CHKERRQ(ierr);
65201a81e61SBarry Smith   PetscFunctionReturn(0);
65301a81e61SBarry Smith }
65401a81e61SBarry Smith EXTERN_C_END
65501a81e61SBarry Smith 
65601a81e61SBarry Smith /**********************************************************************/
65701a81e61SBarry Smith 
65801a81e61SBarry Smith /*
65901a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
66001a81e61SBarry Smith */
66101a81e61SBarry Smith #undef __FUNCT__
66201a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix"
66301a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B)
66401a81e61SBarry Smith {
66501a81e61SBarry Smith   matrix                  *M;
66601a81e61SBarry Smith   int                     i,j,col;
66701a81e61SBarry Smith   int                     row_indx;
66801a81e61SBarry Smith   int                     len,pe,local_indx,start_indx;
66901a81e61SBarry Smith   int                     *mapping;
67001a81e61SBarry Smith   PetscErrorCode          ierr;
67101a81e61SBarry Smith   const int               *cols;
67201a81e61SBarry Smith   const double            *vals;
6732a4967abSBarry Smith   int                     *num_ptr,n,mnl,nnl,nz,rstart,rend;
6742a4967abSBarry Smith   PetscMPIInt             size,rank;
67501a81e61SBarry Smith   struct compressed_lines *rows;
67601a81e61SBarry Smith 
67701a81e61SBarry Smith   PetscFunctionBegin;
67801a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
67901a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
68001a81e61SBarry Smith   ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr);
68101a81e61SBarry Smith   ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr);
68201a81e61SBarry Smith 
68301a81e61SBarry Smith   /*
68401a81e61SBarry Smith     not sure why a barrier is required. commenting out
68501a81e61SBarry Smith   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
68601a81e61SBarry Smith   */
68701a81e61SBarry Smith 
68829b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
68901a81e61SBarry Smith 
69001a81e61SBarry Smith   M->n              = n;
69101a81e61SBarry Smith   M->bs             = 1;
69201a81e61SBarry Smith   M->max_block_size = 1;
69301a81e61SBarry Smith 
69401a81e61SBarry Smith   M->mnls          = (int*)malloc(sizeof(int)*size);
69501a81e61SBarry Smith   M->start_indices = (int*)malloc(sizeof(int)*size);
69601a81e61SBarry Smith   M->pe            = (int*)malloc(sizeof(int)*n);
69701a81e61SBarry Smith   M->block_sizes   = (int*)malloc(sizeof(int)*n);
69801a81e61SBarry Smith   for (i=0; i<n; i++) M->block_sizes[i] = 1;
69901a81e61SBarry Smith 
70001a81e61SBarry Smith   ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr);
70101a81e61SBarry Smith 
70201a81e61SBarry Smith   M->start_indices[0] = 0;
7032fa5cd67SKarl Rupp   for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
70401a81e61SBarry Smith 
70501a81e61SBarry Smith   M->mnl            = M->mnls[M->myid];
70601a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
70701a81e61SBarry Smith 
70801a81e61SBarry Smith   for (i=0; i<size; i++) {
70901a81e61SBarry Smith     start_indx = M->start_indices[i];
7102fa5cd67SKarl Rupp     for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i;
71101a81e61SBarry Smith   }
71201a81e61SBarry Smith 
71301a81e61SBarry Smith   if (AT) {
71401a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr);
71501a81e61SBarry Smith   } else {
71601a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr);
71701a81e61SBarry Smith   }
71801a81e61SBarry Smith 
71901a81e61SBarry Smith   rows = M->lines;
72001a81e61SBarry Smith 
72101a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
72201a81e61SBarry Smith   ierr       = PetscMalloc(M->n*sizeof(int),&mapping);CHKERRQ(ierr);
72301a81e61SBarry Smith   pe         = 0;
72401a81e61SBarry Smith   local_indx = 0;
72501a81e61SBarry Smith   for (i=0; i<M->n; i++) {
72601a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
72701a81e61SBarry Smith       pe++;
72801a81e61SBarry Smith       local_indx = 0;
72901a81e61SBarry Smith     }
73001a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
73101a81e61SBarry Smith     local_indx++;
73201a81e61SBarry Smith   }
73301a81e61SBarry Smith 
73401a81e61SBarry Smith 
73501a81e61SBarry Smith   ierr = PetscMalloc(mnl*sizeof(int),&num_ptr);CHKERRQ(ierr);
73601a81e61SBarry Smith 
73701a81e61SBarry Smith   /*********************************************************/
73801a81e61SBarry Smith   /************** Set up the row structure *****************/
73901a81e61SBarry Smith   /*********************************************************/
74001a81e61SBarry Smith 
74101a81e61SBarry Smith   /* count number of nonzeros in every row */
74201a81e61SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
74301a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
7440298fd71SBarry Smith     ierr = MatGetRow(A,i,&num_ptr[i-rstart],NULL,NULL);CHKERRQ(ierr);
7450298fd71SBarry Smith     ierr = MatRestoreRow(A,i,&num_ptr[i-rstart],NULL,NULL);CHKERRQ(ierr);
74601a81e61SBarry Smith   }
74701a81e61SBarry Smith 
74801a81e61SBarry Smith   /* allocate buffers */
74901a81e61SBarry Smith   len = 0;
75001a81e61SBarry Smith   for (i=0; i<mnl; i++) {
75101a81e61SBarry Smith     if (len < num_ptr[i]) len = num_ptr[i];
75201a81e61SBarry Smith   }
75301a81e61SBarry Smith 
75401a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
75501a81e61SBarry Smith     row_indx             = i-rstart;
75601a81e61SBarry Smith     len                  = num_ptr[row_indx];
75701a81e61SBarry Smith     rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int));
75801a81e61SBarry Smith     rows->A[row_indx]    = (double*)malloc(len*sizeof(double));
75901a81e61SBarry Smith   }
76001a81e61SBarry Smith 
76101a81e61SBarry Smith   /* copy the matrix */
76201a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
76301a81e61SBarry Smith     row_indx = i - rstart;
76401a81e61SBarry Smith     ierr     = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
76501a81e61SBarry Smith     for (j=0; j<nz; j++) {
76601a81e61SBarry Smith       col = cols[j];
76701a81e61SBarry Smith       len = rows->len[row_indx]++;
7682fa5cd67SKarl Rupp 
76901a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
77001a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
77101a81e61SBarry Smith     }
77201a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
7732fa5cd67SKarl Rupp 
77401a81e61SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
77501a81e61SBarry Smith   }
77601a81e61SBarry Smith 
77701a81e61SBarry Smith 
77801a81e61SBarry Smith   /************************************************************/
77901a81e61SBarry Smith   /************** Set up the column structure *****************/
78001a81e61SBarry Smith   /*********************************************************/
78101a81e61SBarry Smith 
78201a81e61SBarry Smith   if (AT) {
78301a81e61SBarry Smith 
78401a81e61SBarry Smith     /* count number of nonzeros in every column */
78501a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
7860298fd71SBarry Smith       ierr = MatGetRow(AT,i,&num_ptr[i-rstart],NULL,NULL);CHKERRQ(ierr);
7870298fd71SBarry Smith       ierr = MatRestoreRow(AT,i,&num_ptr[i-rstart],NULL,NULL);CHKERRQ(ierr);
78801a81e61SBarry Smith     }
78901a81e61SBarry Smith 
79001a81e61SBarry Smith     /* allocate buffers */
79101a81e61SBarry Smith     len = 0;
79201a81e61SBarry Smith     for (i=0; i<mnl; i++) {
79301a81e61SBarry Smith       if (len < num_ptr[i]) len = num_ptr[i];
79401a81e61SBarry Smith     }
79501a81e61SBarry Smith 
79601a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
79701a81e61SBarry Smith       row_indx = i-rstart;
79801a81e61SBarry Smith       len      = num_ptr[row_indx];
7992fa5cd67SKarl Rupp 
80001a81e61SBarry Smith       rows->rptrs[row_indx] = (int*)malloc(len*sizeof(int));
80101a81e61SBarry Smith     }
80201a81e61SBarry Smith 
80301a81e61SBarry Smith     /* copy the matrix (i.e., the structure) */
80401a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
80501a81e61SBarry Smith       row_indx = i - rstart;
80601a81e61SBarry Smith       ierr     = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
80701a81e61SBarry Smith       for (j=0; j<nz; j++) {
80801a81e61SBarry Smith         col = cols[j];
80901a81e61SBarry Smith         len = rows->rlen[row_indx]++;
8102fa5cd67SKarl Rupp 
81101a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
81201a81e61SBarry Smith       }
81301a81e61SBarry Smith       ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
81401a81e61SBarry Smith     }
81501a81e61SBarry Smith   }
81601a81e61SBarry Smith 
81701a81e61SBarry Smith   ierr = PetscFree(num_ptr);CHKERRQ(ierr);
81801a81e61SBarry Smith   ierr = PetscFree(mapping);CHKERRQ(ierr);
81901a81e61SBarry Smith 
82001a81e61SBarry Smith   order_pointers(M);
82101a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
82201a81e61SBarry Smith   *B       = M;
82301a81e61SBarry Smith   PetscFunctionReturn(0);
82401a81e61SBarry Smith }
82501a81e61SBarry Smith 
82601a81e61SBarry Smith /**********************************************************************/
82701a81e61SBarry Smith 
82801a81e61SBarry Smith /*
82901a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
83001a81e61SBarry Smith    This assumes that the the SPAI matrix B is stored in
83101a81e61SBarry Smith    COMPRESSED-ROW format.
83201a81e61SBarry Smith */
83301a81e61SBarry Smith #undef __FUNCT__
83401a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat"
83501a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB)
83601a81e61SBarry Smith {
8374b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
83801a81e61SBarry Smith   PetscErrorCode ierr;
83901a81e61SBarry Smith   int            m,n,M,N;
84001a81e61SBarry Smith   int            d_nz,o_nz;
84101a81e61SBarry Smith   int            *d_nnz,*o_nnz;
84201a81e61SBarry Smith   int            i,k,global_row,global_col,first_diag_col,last_diag_col;
84301a81e61SBarry Smith   PetscScalar    val;
84401a81e61SBarry Smith 
84501a81e61SBarry Smith   PetscFunctionBegin;
84601a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
84701a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
84801a81e61SBarry Smith 
84901a81e61SBarry Smith   m    = n = B->mnls[rank];
85001a81e61SBarry Smith   d_nz = o_nz = 0;
85101a81e61SBarry Smith 
85201a81e61SBarry Smith   /* Determine preallocation for MatCreateMPIAIJ */
8532a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&d_nnz);CHKERRQ(ierr);
8542a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&o_nnz);CHKERRQ(ierr);
85501a81e61SBarry Smith   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
85601a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
85701a81e61SBarry Smith   last_diag_col  = first_diag_col + B->mnls[rank];
85801a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
85901a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
86001a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
861db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
862db4deed7SKarl Rupp       else o_nnz[i]++;
86301a81e61SBarry Smith     }
86401a81e61SBarry Smith   }
86501a81e61SBarry Smith 
86601a81e61SBarry Smith   M = N = B->n;
86701a81e61SBarry Smith   /* Here we only know how to create AIJ format */
868f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,PB);CHKERRQ(ierr);
869f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr);
87001a81e61SBarry Smith   ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr);
87101a81e61SBarry Smith   ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr);
87201a81e61SBarry Smith   ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr);
87301a81e61SBarry Smith 
87401a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
87501a81e61SBarry Smith     global_row = B->start_indices[rank]+i;
87601a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
87701a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
8782fa5cd67SKarl Rupp 
87901a81e61SBarry Smith       val  = B->lines->A[i][k];
88001a81e61SBarry Smith       ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr);
88101a81e61SBarry Smith     }
88201a81e61SBarry Smith   }
88301a81e61SBarry Smith 
88401a81e61SBarry Smith   ierr = PetscFree(d_nnz);CHKERRQ(ierr);
88501a81e61SBarry Smith   ierr = PetscFree(o_nnz);CHKERRQ(ierr);
88601a81e61SBarry Smith 
88701a81e61SBarry Smith   ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88801a81e61SBarry Smith   ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88901a81e61SBarry Smith   PetscFunctionReturn(0);
89001a81e61SBarry Smith }
89101a81e61SBarry Smith 
89201a81e61SBarry Smith /**********************************************************************/
89301a81e61SBarry Smith 
89401a81e61SBarry Smith /*
89501a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
89601a81e61SBarry Smith */
89701a81e61SBarry Smith #undef __FUNCT__
89801a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec"
89901a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv)
90001a81e61SBarry Smith {
90101a81e61SBarry Smith   PetscErrorCode ierr;
9024b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
9034b6b5fe1SBarry Smith   int            m,M,i,*mnls,*start_indices,*global_indices;
90401a81e61SBarry Smith 
90501a81e61SBarry Smith   PetscFunctionBegin;
90601a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
90701a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
90801a81e61SBarry Smith 
90901a81e61SBarry Smith   m = v->mnl;
91001a81e61SBarry Smith   M = v->n;
91101a81e61SBarry Smith 
91201a81e61SBarry Smith 
91301a81e61SBarry Smith   ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr);
91401a81e61SBarry Smith 
91501a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&mnls);CHKERRQ(ierr);
9162a4967abSBarry Smith   ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRQ(ierr);
91701a81e61SBarry Smith 
91801a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&start_indices);CHKERRQ(ierr);
9192fa5cd67SKarl Rupp 
92001a81e61SBarry Smith   start_indices[0] = 0;
9212fa5cd67SKarl Rupp   for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1];
92201a81e61SBarry Smith 
92301a81e61SBarry Smith   ierr = PetscMalloc(v->mnl*sizeof(int),&global_indices);CHKERRQ(ierr);
9242fa5cd67SKarl Rupp   for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i;
92501a81e61SBarry Smith 
92601a81e61SBarry Smith   ierr = PetscFree(mnls);CHKERRQ(ierr);
92701a81e61SBarry Smith   ierr = PetscFree(start_indices);CHKERRQ(ierr);
92801a81e61SBarry Smith 
92901a81e61SBarry Smith   ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr);
93001a81e61SBarry Smith   ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr);
93101a81e61SBarry Smith   ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr);
93201a81e61SBarry Smith 
9334b6b5fe1SBarry Smith   ierr = PetscFree(global_indices);CHKERRQ(ierr);
93401a81e61SBarry Smith   PetscFunctionReturn(0);
93501a81e61SBarry Smith }
93601a81e61SBarry Smith 
93701a81e61SBarry Smith 
93801a81e61SBarry Smith 
93901a81e61SBarry Smith 
940