xref: /petsc/src/ksp/pc/impls/spai/ispai.c (revision 5fd668637986a8d8518383a9159eebc368e1d5b4)
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 */
9201a81e61SBarry Smith   /* int cache_size */     /* one of (1,2,3,4,5,6) indicting size of cache */
9301a81e61SBarry Smith                            /* cache_size == 0 indicates no caching */
9401a81e61SBarry Smith   /* int    verbose    */  /* verbose == 0 specifies that SPAI is silent
9501a81e61SBarry Smith                               verbose == 1 prints timing and matrix statistics */
9601a81e61SBarry Smith 
9701a81e61SBarry Smith   ierr = bspai(ispai->B,&ispai->M,
9801a81e61SBarry Smith                stdout,
9901a81e61SBarry Smith                ispai->epsilon,
10001a81e61SBarry Smith                ispai->nbsteps,
10101a81e61SBarry Smith                ispai->max,
10201a81e61SBarry Smith                ispai->maxnew,
10301a81e61SBarry Smith                ispai->block_size,
10401a81e61SBarry Smith                ispai->cache_size,
10501a81e61SBarry Smith                ispai->verbose);CHKERRQ(ierr);
10601a81e61SBarry Smith 
1077adad957SLisandro Dalcin   ierr = ConvertMatrixToMat(((PetscObject)pc)->comm,ispai->M,&ispai->PM);CHKERRQ(ierr);
10801a81e61SBarry Smith 
10901a81e61SBarry Smith   /* free the SPAI matrices */
11001a81e61SBarry Smith   sp_free_matrix(ispai->B);
11101a81e61SBarry Smith   sp_free_matrix(ispai->M);
11201a81e61SBarry Smith 
11301a81e61SBarry Smith   PetscFunctionReturn(0);
11401a81e61SBarry Smith }
11501a81e61SBarry Smith 
11601a81e61SBarry Smith /**********************************************************************/
11701a81e61SBarry Smith 
11801a81e61SBarry Smith #undef __FUNCT__
11901a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI"
12001a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y)
12101a81e61SBarry Smith {
12201a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
12301a81e61SBarry Smith   PetscErrorCode ierr;
12401a81e61SBarry Smith 
12501a81e61SBarry Smith   PetscFunctionBegin;
12601a81e61SBarry Smith   /* Now using PETSc's multiply */
12701a81e61SBarry Smith   ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr);
12801a81e61SBarry Smith   PetscFunctionReturn(0);
12901a81e61SBarry Smith }
13001a81e61SBarry Smith 
13101a81e61SBarry Smith /**********************************************************************/
13201a81e61SBarry Smith 
13301a81e61SBarry Smith #undef __FUNCT__
13401a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI"
13501a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc)
13601a81e61SBarry Smith {
13701a81e61SBarry Smith   PetscErrorCode ierr;
13801a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
13901a81e61SBarry Smith 
14001a81e61SBarry Smith   PetscFunctionBegin;
14169a3b875SJed Brown   ierr = MatDestroy(&ispai->PM);CHKERRQ(ierr);
14201a81e61SBarry Smith   ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr);
143c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
14401a81e61SBarry Smith   PetscFunctionReturn(0);
14501a81e61SBarry Smith }
14601a81e61SBarry Smith 
14701a81e61SBarry Smith /**********************************************************************/
14801a81e61SBarry Smith 
14901a81e61SBarry Smith #undef __FUNCT__
15001a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI"
15101a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer)
15201a81e61SBarry Smith {
15301a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
15401a81e61SBarry Smith   PetscErrorCode ierr;
155ace3abfcSBarry Smith   PetscBool      iascii;
15601a81e61SBarry Smith 
15701a81e61SBarry Smith   PetscFunctionBegin;
158251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
15901a81e61SBarry Smith   if (iascii) {
16001a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    SPAI preconditioner\n");CHKERRQ(ierr);
161a83599f4SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    epsilon %G\n",   ispai->epsilon);CHKERRQ(ierr);
16201a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    nbsteps %d\n",   ispai->nbsteps);CHKERRQ(ierr);
16301a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    max %d\n",       ispai->max);CHKERRQ(ierr);
16401a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    maxnew %d\n",    ispai->maxnew);CHKERRQ(ierr);
16501a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    block_size %d\n",ispai->block_size);CHKERRQ(ierr);
16601a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    cache_size %d\n",ispai->cache_size);CHKERRQ(ierr);
16701a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    verbose %d\n",   ispai->verbose);CHKERRQ(ierr);
16801a81e61SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"    sp %d\n",        ispai->sp);CHKERRQ(ierr);
16901a81e61SBarry Smith   }
17001a81e61SBarry Smith   PetscFunctionReturn(0);
17101a81e61SBarry Smith }
17201a81e61SBarry Smith 
17301a81e61SBarry Smith EXTERN_C_BEGIN
17401a81e61SBarry Smith #undef __FUNCT__
17501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI"
1767087cfbeSBarry Smith PetscErrorCode  PCSPAISetEpsilon_SPAI(PC pc,double epsilon1)
17701a81e61SBarry Smith {
17801a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
179*5fd66863SKarl Rupp 
18001a81e61SBarry Smith   PetscFunctionBegin;
18101a81e61SBarry Smith   ispai->epsilon = epsilon1;
18201a81e61SBarry Smith   PetscFunctionReturn(0);
18301a81e61SBarry Smith }
18401a81e61SBarry Smith EXTERN_C_END
18501a81e61SBarry Smith 
18601a81e61SBarry Smith /**********************************************************************/
18701a81e61SBarry Smith 
18801a81e61SBarry Smith EXTERN_C_BEGIN
18901a81e61SBarry Smith #undef __FUNCT__
19001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI"
1917087cfbeSBarry Smith PetscErrorCode  PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1)
19201a81e61SBarry Smith {
19301a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
194*5fd66863SKarl Rupp 
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"
2077087cfbeSBarry Smith PetscErrorCode  PCSPAISetMax_SPAI(PC pc,int max1)
20801a81e61SBarry Smith {
20901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
210*5fd66863SKarl Rupp 
21101a81e61SBarry Smith   PetscFunctionBegin;
21201a81e61SBarry Smith   ispai->max = max1;
21301a81e61SBarry Smith   PetscFunctionReturn(0);
21401a81e61SBarry Smith }
21501a81e61SBarry Smith EXTERN_C_END
21601a81e61SBarry Smith 
21701a81e61SBarry Smith /**********************************************************************/
21801a81e61SBarry Smith 
21901a81e61SBarry Smith EXTERN_C_BEGIN
22001a81e61SBarry Smith #undef __FUNCT__
22101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI"
2227087cfbeSBarry Smith PetscErrorCode  PCSPAISetMaxNew_SPAI(PC pc,int maxnew1)
22301a81e61SBarry Smith {
22401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
225*5fd66863SKarl Rupp 
22601a81e61SBarry Smith   PetscFunctionBegin;
22701a81e61SBarry Smith   ispai->maxnew = maxnew1;
22801a81e61SBarry Smith   PetscFunctionReturn(0);
22901a81e61SBarry Smith }
23001a81e61SBarry Smith EXTERN_C_END
23101a81e61SBarry Smith 
23201a81e61SBarry Smith /**********************************************************************/
23301a81e61SBarry Smith 
23401a81e61SBarry Smith EXTERN_C_BEGIN
23501a81e61SBarry Smith #undef __FUNCT__
23601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize_SPAI"
2377087cfbeSBarry Smith PetscErrorCode  PCSPAISetBlockSize_SPAI(PC pc,int block_size1)
23801a81e61SBarry Smith {
23901a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
240*5fd66863SKarl Rupp 
24101a81e61SBarry Smith   PetscFunctionBegin;
24201a81e61SBarry Smith   ispai->block_size = block_size1;
24301a81e61SBarry Smith   PetscFunctionReturn(0);
24401a81e61SBarry Smith }
24501a81e61SBarry Smith EXTERN_C_END
24601a81e61SBarry Smith 
24701a81e61SBarry Smith /**********************************************************************/
24801a81e61SBarry Smith 
24901a81e61SBarry Smith EXTERN_C_BEGIN
25001a81e61SBarry Smith #undef __FUNCT__
25101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI"
2527087cfbeSBarry Smith PetscErrorCode  PCSPAISetCacheSize_SPAI(PC pc,int cache_size)
25301a81e61SBarry Smith {
25401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
255*5fd66863SKarl Rupp 
25601a81e61SBarry Smith   PetscFunctionBegin;
25701a81e61SBarry Smith   ispai->cache_size = cache_size;
25801a81e61SBarry Smith   PetscFunctionReturn(0);
25901a81e61SBarry Smith }
26001a81e61SBarry Smith EXTERN_C_END
26101a81e61SBarry Smith 
26201a81e61SBarry Smith /**********************************************************************/
26301a81e61SBarry Smith 
26401a81e61SBarry Smith EXTERN_C_BEGIN
26501a81e61SBarry Smith #undef __FUNCT__
26601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI"
2677087cfbeSBarry Smith PetscErrorCode  PCSPAISetVerbose_SPAI(PC pc,int verbose)
26801a81e61SBarry Smith {
26901a81e61SBarry Smith   PC_SPAI    *ispai = (PC_SPAI*)pc->data;
270*5fd66863SKarl Rupp 
27101a81e61SBarry Smith   PetscFunctionBegin;
27201a81e61SBarry Smith   ispai->verbose = verbose;
27301a81e61SBarry Smith   PetscFunctionReturn(0);
27401a81e61SBarry Smith }
27501a81e61SBarry Smith EXTERN_C_END
27601a81e61SBarry Smith 
27701a81e61SBarry Smith /**********************************************************************/
27801a81e61SBarry Smith 
27901a81e61SBarry Smith EXTERN_C_BEGIN
28001a81e61SBarry Smith #undef __FUNCT__
28101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI"
2827087cfbeSBarry Smith PetscErrorCode  PCSPAISetSp_SPAI(PC pc,int sp)
28301a81e61SBarry Smith {
28401a81e61SBarry Smith   PC_SPAI *ispai = (PC_SPAI*)pc->data;
285*5fd66863SKarl Rupp 
28601a81e61SBarry Smith   PetscFunctionBegin;
28701a81e61SBarry Smith   ispai->sp = sp;
28801a81e61SBarry Smith   PetscFunctionReturn(0);
28901a81e61SBarry Smith }
29001a81e61SBarry Smith EXTERN_C_END
29101a81e61SBarry Smith 
29201a81e61SBarry Smith /* -------------------------------------------------------------------*/
29301a81e61SBarry Smith 
29401a81e61SBarry Smith #undef __FUNCT__
29501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon"
29601a81e61SBarry Smith /*@
29701a81e61SBarry Smith   PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner
29801a81e61SBarry Smith 
29901a81e61SBarry Smith   Input Parameters:
30001a81e61SBarry Smith + pc - the preconditioner
30101a81e61SBarry Smith - eps - epsilon (default .4)
30201a81e61SBarry Smith 
30301a81e61SBarry Smith   Notes:  Espilon must be between 0 and 1. It controls the
30401a81e61SBarry Smith                  quality of the approximation of M to the inverse of
30501a81e61SBarry Smith                  A. Higher values of epsilon lead to more work, more
30601a81e61SBarry Smith                  fill, and usually better preconditioners. In many
30701a81e61SBarry Smith                  cases the best choice of epsilon is the one that
30801a81e61SBarry Smith                  divides the total solution time equally between the
30901a81e61SBarry Smith                  preconditioner and the solver.
31001a81e61SBarry Smith 
31101a81e61SBarry Smith   Level: intermediate
31201a81e61SBarry Smith 
31301a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
31401a81e61SBarry Smith   @*/
3157087cfbeSBarry Smith PetscErrorCode  PCSPAISetEpsilon(PC pc,double epsilon1)
31601a81e61SBarry Smith {
3174ac538c5SBarry Smith   PetscErrorCode ierr;
318*5fd66863SKarl Rupp 
31901a81e61SBarry Smith   PetscFunctionBegin;
3204ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));CHKERRQ(ierr);
32101a81e61SBarry Smith   PetscFunctionReturn(0);
32201a81e61SBarry Smith }
32301a81e61SBarry Smith 
32401a81e61SBarry Smith /**********************************************************************/
32501a81e61SBarry Smith 
32601a81e61SBarry Smith #undef __FUNCT__
32701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps"
32801a81e61SBarry Smith /*@
32901a81e61SBarry Smith   PCSPAISetNBSteps - set maximum number of improvement steps per row in
33001a81e61SBarry Smith         the SPAI preconditioner
33101a81e61SBarry Smith 
33201a81e61SBarry Smith   Input Parameters:
33301a81e61SBarry Smith + pc - the preconditioner
33401a81e61SBarry Smith - n - number of steps (default 5)
33501a81e61SBarry Smith 
33601a81e61SBarry Smith   Notes:  SPAI constructs to approximation to every column of
33701a81e61SBarry Smith                  the exact inverse of A in a series of improvement
33801a81e61SBarry Smith                  steps. The quality of the approximation is determined
33901a81e61SBarry Smith                  by epsilon. If an approximation achieving an accuracy
34001a81e61SBarry Smith                  of epsilon is not obtained after ns steps, SPAI simply
34101a81e61SBarry Smith                  uses the best approximation constructed so far.
34201a81e61SBarry Smith 
34301a81e61SBarry Smith   Level: intermediate
34401a81e61SBarry Smith 
34501a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew()
34601a81e61SBarry Smith @*/
3477087cfbeSBarry Smith PetscErrorCode  PCSPAISetNBSteps(PC pc,int nbsteps1)
34801a81e61SBarry Smith {
3494ac538c5SBarry Smith   PetscErrorCode ierr;
350*5fd66863SKarl Rupp 
35101a81e61SBarry Smith   PetscFunctionBegin;
3524ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));CHKERRQ(ierr);
35301a81e61SBarry Smith   PetscFunctionReturn(0);
35401a81e61SBarry Smith }
35501a81e61SBarry Smith 
35601a81e61SBarry Smith /**********************************************************************/
35701a81e61SBarry Smith 
35801a81e61SBarry Smith /* added 1/7/99 g.h. */
35901a81e61SBarry Smith #undef __FUNCT__
36001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax"
36101a81e61SBarry Smith /*@
36201a81e61SBarry Smith   PCSPAISetMax - set the size of various working buffers in
36301a81e61SBarry Smith         the SPAI preconditioner
36401a81e61SBarry Smith 
36501a81e61SBarry Smith   Input Parameters:
36601a81e61SBarry Smith + pc - the preconditioner
36701a81e61SBarry Smith - n - size (default is 5000)
36801a81e61SBarry Smith 
36901a81e61SBarry Smith   Level: intermediate
37001a81e61SBarry Smith 
37101a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
37201a81e61SBarry Smith @*/
3737087cfbeSBarry Smith PetscErrorCode  PCSPAISetMax(PC pc,int max1)
37401a81e61SBarry Smith {
3754ac538c5SBarry Smith   PetscErrorCode ierr;
376*5fd66863SKarl Rupp 
37701a81e61SBarry Smith   PetscFunctionBegin;
3784ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));CHKERRQ(ierr);
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 @*/
3987087cfbeSBarry Smith PetscErrorCode  PCSPAISetMaxNew(PC pc,int maxnew1)
39901a81e61SBarry Smith {
4004ac538c5SBarry Smith   PetscErrorCode ierr;
401*5fd66863SKarl Rupp 
40201a81e61SBarry Smith   PetscFunctionBegin;
4034ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));CHKERRQ(ierr);
40401a81e61SBarry Smith   PetscFunctionReturn(0);
40501a81e61SBarry Smith }
40601a81e61SBarry Smith 
40701a81e61SBarry Smith /**********************************************************************/
40801a81e61SBarry Smith 
40901a81e61SBarry Smith #undef __FUNCT__
41001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize"
41101a81e61SBarry Smith /*@
41201a81e61SBarry Smith   PCSPAISetBlockSize - set the block size for the SPAI preconditioner
41301a81e61SBarry Smith 
41401a81e61SBarry Smith   Input Parameters:
41501a81e61SBarry Smith + pc - the preconditioner
41601a81e61SBarry Smith - n - block size (default 1)
41701a81e61SBarry Smith 
41801a81e61SBarry Smith   Notes: A block
41901a81e61SBarry Smith                  size of 1 treats A as a matrix of scalar elements. A
42001a81e61SBarry Smith                  block size of s > 1 treats A as a matrix of sxs
42101a81e61SBarry Smith                  blocks. A block size of 0 treats A as a matrix with
42201a81e61SBarry Smith                  variable sized blocks, which are determined by
42301a81e61SBarry Smith                  searching for dense square diagonal blocks in A.
42401a81e61SBarry Smith                  This can be very effective for finite-element
42501a81e61SBarry Smith                  matrices.
42601a81e61SBarry Smith 
42701a81e61SBarry Smith                  SPAI will convert A to block form, use a block
42801a81e61SBarry Smith                  version of the preconditioner algorithm, and then
42901a81e61SBarry Smith                  convert the result back to scalar form.
43001a81e61SBarry Smith 
43101a81e61SBarry Smith                  In many cases the a block-size parameter other than 1
43201a81e61SBarry Smith                  can lead to very significant improvement in
43301a81e61SBarry Smith                  performance.
43401a81e61SBarry Smith 
43501a81e61SBarry Smith 
43601a81e61SBarry Smith   Level: intermediate
43701a81e61SBarry Smith 
43801a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
43901a81e61SBarry Smith @*/
4407087cfbeSBarry Smith PetscErrorCode  PCSPAISetBlockSize(PC pc,int block_size1)
44101a81e61SBarry Smith {
4424ac538c5SBarry Smith   PetscErrorCode ierr;
443*5fd66863SKarl Rupp 
44401a81e61SBarry Smith   PetscFunctionBegin;
4454ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));CHKERRQ(ierr);
44601a81e61SBarry Smith   PetscFunctionReturn(0);
44701a81e61SBarry Smith }
44801a81e61SBarry Smith 
44901a81e61SBarry Smith /**********************************************************************/
45001a81e61SBarry Smith 
45101a81e61SBarry Smith #undef __FUNCT__
45201a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize"
45301a81e61SBarry Smith /*@
45401a81e61SBarry Smith   PCSPAISetCacheSize - specify cache size in the SPAI preconditioner
45501a81e61SBarry Smith 
45601a81e61SBarry Smith   Input Parameters:
45701a81e61SBarry Smith + pc - the preconditioner
45801a81e61SBarry Smith - n -  cache size {0,1,2,3,4,5} (default 5)
45901a81e61SBarry Smith 
46001a81e61SBarry Smith   Notes:    SPAI uses a hash table to cache messages and avoid
46101a81e61SBarry Smith                  redundant communication. If suggest always using
46201a81e61SBarry Smith                  5. This parameter is irrelevant in the serial
46301a81e61SBarry Smith                  version.
46401a81e61SBarry Smith 
46501a81e61SBarry Smith   Level: intermediate
46601a81e61SBarry Smith 
46701a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
46801a81e61SBarry Smith @*/
4697087cfbeSBarry Smith PetscErrorCode  PCSPAISetCacheSize(PC pc,int cache_size)
47001a81e61SBarry Smith {
4714ac538c5SBarry Smith   PetscErrorCode ierr;
472*5fd66863SKarl Rupp 
47301a81e61SBarry Smith   PetscFunctionBegin;
4744ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));CHKERRQ(ierr);
47501a81e61SBarry Smith   PetscFunctionReturn(0);
47601a81e61SBarry Smith }
47701a81e61SBarry Smith 
47801a81e61SBarry Smith /**********************************************************************/
47901a81e61SBarry Smith 
48001a81e61SBarry Smith #undef __FUNCT__
48101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose"
48201a81e61SBarry Smith /*@
48301a81e61SBarry Smith   PCSPAISetVerbose - verbosity level for the SPAI preconditioner
48401a81e61SBarry Smith 
48501a81e61SBarry Smith   Input Parameters:
48601a81e61SBarry Smith + pc - the preconditioner
48701a81e61SBarry Smith - n - level (default 1)
48801a81e61SBarry Smith 
48901a81e61SBarry Smith   Notes: print parameters, timings and matrix statistics
49001a81e61SBarry Smith 
49101a81e61SBarry Smith   Level: intermediate
49201a81e61SBarry Smith 
49301a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
49401a81e61SBarry Smith @*/
4957087cfbeSBarry Smith PetscErrorCode  PCSPAISetVerbose(PC pc,int verbose)
49601a81e61SBarry Smith {
4974ac538c5SBarry Smith   PetscErrorCode ierr;
498*5fd66863SKarl Rupp 
49901a81e61SBarry Smith   PetscFunctionBegin;
5004ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));CHKERRQ(ierr);
50101a81e61SBarry Smith   PetscFunctionReturn(0);
50201a81e61SBarry Smith }
50301a81e61SBarry Smith 
50401a81e61SBarry Smith /**********************************************************************/
50501a81e61SBarry Smith 
50601a81e61SBarry Smith #undef __FUNCT__
50701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp"
50801a81e61SBarry Smith /*@
50901a81e61SBarry Smith   PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner
51001a81e61SBarry Smith 
51101a81e61SBarry Smith   Input Parameters:
51201a81e61SBarry Smith + pc - the preconditioner
51301a81e61SBarry Smith - n - 0 or 1
51401a81e61SBarry Smith 
51501a81e61SBarry Smith   Notes: If A has a symmetric nonzero pattern use -sp 1 to
51601a81e61SBarry Smith                  improve performance by eliminating some communication
51701a81e61SBarry Smith                  in the parallel version. Even if A does not have a
51801a81e61SBarry Smith                  symmetric nonzero pattern -sp 1 may well lead to good
51901a81e61SBarry Smith                  results, but the code will not follow the published
52001a81e61SBarry Smith                  SPAI algorithm exactly.
52101a81e61SBarry Smith 
52201a81e61SBarry Smith 
52301a81e61SBarry Smith   Level: intermediate
52401a81e61SBarry Smith 
52501a81e61SBarry Smith .seealso: PCSPAI, PCSetType()
52601a81e61SBarry Smith @*/
5277087cfbeSBarry Smith PetscErrorCode  PCSPAISetSp(PC pc,int sp)
52801a81e61SBarry Smith {
5294ac538c5SBarry Smith   PetscErrorCode ierr;
530*5fd66863SKarl Rupp 
53101a81e61SBarry Smith   PetscFunctionBegin;
5324ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));CHKERRQ(ierr);
53301a81e61SBarry Smith   PetscFunctionReturn(0);
53401a81e61SBarry Smith }
53501a81e61SBarry Smith 
53601a81e61SBarry Smith /**********************************************************************/
53701a81e61SBarry Smith 
53801a81e61SBarry Smith /**********************************************************************/
53901a81e61SBarry Smith 
54001a81e61SBarry Smith #undef __FUNCT__
54101a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI"
54201a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc)
54301a81e61SBarry Smith {
54401a81e61SBarry Smith   PC_SPAI        *ispai = (PC_SPAI*)pc->data;
54501a81e61SBarry Smith   PetscErrorCode ierr;
54601a81e61SBarry Smith   int            nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp;
54701a81e61SBarry Smith   double         epsilon1;
548ace3abfcSBarry Smith   PetscBool      flg;
54901a81e61SBarry Smith 
55001a81e61SBarry Smith   PetscFunctionBegin;
55101a81e61SBarry Smith   ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr);
55201a81e61SBarry Smith     ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr);
55301a81e61SBarry Smith     if (flg) {
55401a81e61SBarry Smith       ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr);
55501a81e61SBarry Smith     }
55601a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr);
55701a81e61SBarry Smith     if (flg) {
55801a81e61SBarry Smith       ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr);
55901a81e61SBarry Smith     }
56001a81e61SBarry Smith     /* added 1/7/99 g.h. */
56101a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr);
56201a81e61SBarry Smith     if (flg) {
56301a81e61SBarry Smith       ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr);
56401a81e61SBarry Smith     }
56501a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr);
56601a81e61SBarry Smith     if (flg) {
56701a81e61SBarry Smith       ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr);
56801a81e61SBarry Smith     }
56901a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr);
57001a81e61SBarry Smith     if (flg) {
57101a81e61SBarry Smith       ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr);
57201a81e61SBarry Smith     }
57301a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr);
57401a81e61SBarry Smith     if (flg) {
57501a81e61SBarry Smith       ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr);
57601a81e61SBarry Smith     }
57701a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr);
57801a81e61SBarry Smith     if (flg) {
57901a81e61SBarry Smith       ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr);
58001a81e61SBarry Smith     }
58101a81e61SBarry Smith     ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr);
58201a81e61SBarry Smith     if (flg) {
58301a81e61SBarry Smith       ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr);
58401a81e61SBarry Smith     }
58501a81e61SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
58601a81e61SBarry Smith   PetscFunctionReturn(0);
58701a81e61SBarry Smith }
58801a81e61SBarry Smith 
58901a81e61SBarry Smith /**********************************************************************/
59001a81e61SBarry Smith 
59101a81e61SBarry Smith /*MC
59201a81e61SBarry Smith    PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard
59301a81e61SBarry Smith      as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3)
59401a81e61SBarry Smith 
59501a81e61SBarry Smith    Options Database Keys:
596c5ae99e3SSatish Balay +  -pc_spai_epsilon <eps> - set tolerance
59701a81e61SBarry Smith .  -pc_spai_nbstep <n> - set nbsteps
59801a81e61SBarry Smith .  -pc_spai_max <m> - set max
59901a81e61SBarry Smith .  -pc_spai_max_new <m> - set maxnew
60001a81e61SBarry Smith .  -pc_spai_block_size <n> - set block size
60101a81e61SBarry Smith .  -pc_spai_cache_size <n> - set cache size
60201a81e61SBarry Smith .  -pc_spai_sp <m> - set sp
60301a81e61SBarry Smith -  -pc_spai_set_verbose <true,false> - verbose output
60401a81e61SBarry Smith 
60501a81e61SBarry Smith    Notes: This only works with AIJ matrices.
60601a81e61SBarry Smith 
60701a81e61SBarry Smith    Level: beginner
60801a81e61SBarry Smith 
60901a81e61SBarry Smith    Concepts: approximate inverse
61001a81e61SBarry Smith 
61101a81e61SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
61201a81e61SBarry Smith     PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(),
61301a81e61SBarry Smith     PCSPAISetVerbose(), PCSPAISetSp()
61401a81e61SBarry Smith M*/
61501a81e61SBarry Smith 
61601a81e61SBarry Smith EXTERN_C_BEGIN
61701a81e61SBarry Smith #undef __FUNCT__
61801a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI"
6197087cfbeSBarry Smith PetscErrorCode  PCCreate_SPAI(PC pc)
62001a81e61SBarry Smith {
62101a81e61SBarry Smith   PC_SPAI        *ispai;
62201a81e61SBarry Smith   PetscErrorCode ierr;
62301a81e61SBarry Smith 
62401a81e61SBarry Smith   PetscFunctionBegin;
62538f2d2fdSLisandro Dalcin   ierr               = PetscNewLog(pc,PC_SPAI,&ispai);CHKERRQ(ierr);
6262a4967abSBarry Smith   pc->data           = ispai;
62701a81e61SBarry Smith 
62801a81e61SBarry Smith   pc->ops->destroy         = PCDestroy_SPAI;
62901a81e61SBarry Smith   pc->ops->apply           = PCApply_SPAI;
63001a81e61SBarry Smith   pc->ops->applyrichardson = 0;
63101a81e61SBarry Smith   pc->ops->setup           = PCSetUp_SPAI;
63201a81e61SBarry Smith   pc->ops->view            = PCView_SPAI;
63301a81e61SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_SPAI;
63401a81e61SBarry Smith 
63501a81e61SBarry Smith   ispai->epsilon    = .4;
63601a81e61SBarry Smith   ispai->nbsteps    = 5;
63701a81e61SBarry Smith   ispai->max        = 5000;
63801a81e61SBarry Smith   ispai->maxnew     = 5;
63901a81e61SBarry Smith   ispai->block_size = 1;
64001a81e61SBarry Smith   ispai->cache_size = 5;
64101a81e61SBarry Smith   ispai->verbose    = 0;
64201a81e61SBarry Smith 
64301a81e61SBarry Smith   ispai->sp         = 1;
6447adad957SLisandro Dalcin   ierr = MPI_Comm_dup(((PetscObject)pc)->comm,&(ispai->comm_spai));CHKERRQ(ierr);
64501a81e61SBarry Smith 
64601a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C",
64701a81e61SBarry Smith                     "PCSPAISetEpsilon_SPAI",
64801a81e61SBarry Smith                      PCSPAISetEpsilon_SPAI);CHKERRQ(ierr);
64901a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C",
65001a81e61SBarry Smith                     "PCSPAISetNBSteps_SPAI",
65101a81e61SBarry Smith                      PCSPAISetNBSteps_SPAI);CHKERRQ(ierr);
65201a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C",
65301a81e61SBarry Smith                     "PCSPAISetMax_SPAI",
65401a81e61SBarry Smith                      PCSPAISetMax_SPAI);CHKERRQ(ierr);
65501a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC",
65601a81e61SBarry Smith                     "PCSPAISetMaxNew_SPAI",
65701a81e61SBarry Smith                      PCSPAISetMaxNew_SPAI);CHKERRQ(ierr);
65801a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C",
65901a81e61SBarry Smith                     "PCSPAISetBlockSize_SPAI",
66001a81e61SBarry Smith                      PCSPAISetBlockSize_SPAI);CHKERRQ(ierr);
66101a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C",
66201a81e61SBarry Smith                     "PCSPAISetCacheSize_SPAI",
66301a81e61SBarry Smith                      PCSPAISetCacheSize_SPAI);CHKERRQ(ierr);
66401a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C",
66501a81e61SBarry Smith                     "PCSPAISetVerbose_SPAI",
66601a81e61SBarry Smith                      PCSPAISetVerbose_SPAI);CHKERRQ(ierr);
66701a81e61SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C",
66801a81e61SBarry Smith                     "PCSPAISetSp_SPAI",
66901a81e61SBarry Smith                      PCSPAISetSp_SPAI);CHKERRQ(ierr);
67001a81e61SBarry Smith 
67101a81e61SBarry Smith   PetscFunctionReturn(0);
67201a81e61SBarry Smith }
67301a81e61SBarry Smith EXTERN_C_END
67401a81e61SBarry Smith 
67501a81e61SBarry Smith /**********************************************************************/
67601a81e61SBarry Smith 
67701a81e61SBarry Smith /*
67801a81e61SBarry Smith    Converts from a PETSc matrix to an SPAI matrix
67901a81e61SBarry Smith */
68001a81e61SBarry Smith #undef __FUNCT__
68101a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix"
68201a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B)
68301a81e61SBarry Smith {
68401a81e61SBarry Smith   matrix                  *M;
68501a81e61SBarry Smith   int                     i,j,col;
68601a81e61SBarry Smith   int                     row_indx;
68701a81e61SBarry Smith   int                     len,pe,local_indx,start_indx;
68801a81e61SBarry Smith   int                     *mapping;
68901a81e61SBarry Smith   PetscErrorCode          ierr;
69001a81e61SBarry Smith   const int               *cols;
69101a81e61SBarry Smith   const double            *vals;
6922a4967abSBarry Smith   int                     *num_ptr,n,mnl,nnl,nz,rstart,rend;
6932a4967abSBarry Smith   PetscMPIInt             size,rank;
69401a81e61SBarry Smith   struct compressed_lines *rows;
69501a81e61SBarry Smith 
69601a81e61SBarry Smith   PetscFunctionBegin;
69701a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
69801a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
69901a81e61SBarry Smith   ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr);
70001a81e61SBarry Smith   ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr);
70101a81e61SBarry Smith 
70201a81e61SBarry Smith   /*
70301a81e61SBarry Smith     not sure why a barrier is required. commenting out
70401a81e61SBarry Smith   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
70501a81e61SBarry Smith   */
70601a81e61SBarry Smith 
70729b2f704SSatish Balay   M = new_matrix((SPAI_Comm)comm);
70801a81e61SBarry Smith 
70901a81e61SBarry Smith   M->n = n;
71001a81e61SBarry Smith   M->bs = 1;
71101a81e61SBarry Smith   M->max_block_size = 1;
71201a81e61SBarry Smith 
71301a81e61SBarry Smith   M->mnls          = (int*)malloc(sizeof(int)*size);
71401a81e61SBarry Smith   M->start_indices = (int*)malloc(sizeof(int)*size);
71501a81e61SBarry Smith   M->pe            = (int*)malloc(sizeof(int)*n);
71601a81e61SBarry Smith   M->block_sizes   = (int*)malloc(sizeof(int)*n);
71701a81e61SBarry Smith   for (i=0; i<n; i++) M->block_sizes[i] = 1;
71801a81e61SBarry Smith 
71901a81e61SBarry Smith   ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr);
72001a81e61SBarry Smith 
72101a81e61SBarry Smith   M->start_indices[0] = 0;
72201a81e61SBarry Smith   for (i=1; i<size; i++) {
72301a81e61SBarry Smith     M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1];
72401a81e61SBarry Smith   }
72501a81e61SBarry Smith 
72601a81e61SBarry Smith   M->mnl = M->mnls[M->myid];
72701a81e61SBarry Smith   M->my_start_index = M->start_indices[M->myid];
72801a81e61SBarry Smith 
72901a81e61SBarry Smith   for (i=0; i<size; i++) {
73001a81e61SBarry Smith     start_indx = M->start_indices[i];
73101a81e61SBarry Smith     for (j=0; j<M->mnls[i]; j++)
73201a81e61SBarry Smith       M->pe[start_indx+j] = i;
73301a81e61SBarry Smith   }
73401a81e61SBarry Smith 
73501a81e61SBarry Smith   if (AT) {
73601a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr);
73701a81e61SBarry Smith   } else {
73801a81e61SBarry Smith     M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr);
73901a81e61SBarry Smith   }
74001a81e61SBarry Smith 
74101a81e61SBarry Smith   rows = M->lines;
74201a81e61SBarry Smith 
74301a81e61SBarry Smith   /* Determine the mapping from global indices to pointers */
74401a81e61SBarry Smith   ierr       = PetscMalloc(M->n*sizeof(int),&mapping);CHKERRQ(ierr);
74501a81e61SBarry Smith   pe         = 0;
74601a81e61SBarry Smith   local_indx = 0;
74701a81e61SBarry Smith   for (i=0; i<M->n; i++) {
74801a81e61SBarry Smith     if (local_indx >= M->mnls[pe]) {
74901a81e61SBarry Smith       pe++;
75001a81e61SBarry Smith       local_indx = 0;
75101a81e61SBarry Smith     }
75201a81e61SBarry Smith     mapping[i] = local_indx + M->start_indices[pe];
75301a81e61SBarry Smith     local_indx++;
75401a81e61SBarry Smith   }
75501a81e61SBarry Smith 
75601a81e61SBarry Smith 
75701a81e61SBarry Smith   ierr = PetscMalloc(mnl*sizeof(int),&num_ptr);CHKERRQ(ierr);
75801a81e61SBarry Smith 
75901a81e61SBarry Smith   /*********************************************************/
76001a81e61SBarry Smith   /************** Set up the row structure *****************/
76101a81e61SBarry Smith   /*********************************************************/
76201a81e61SBarry Smith 
76301a81e61SBarry Smith   /* count number of nonzeros in every row */
76401a81e61SBarry Smith   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
76501a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
76601a81e61SBarry Smith     ierr = MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
76701a81e61SBarry Smith     ierr = MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
76801a81e61SBarry Smith   }
76901a81e61SBarry Smith 
77001a81e61SBarry Smith   /* allocate buffers */
77101a81e61SBarry Smith   len = 0;
77201a81e61SBarry Smith   for (i=0; i<mnl; i++) {
77301a81e61SBarry Smith     if (len < num_ptr[i]) len = num_ptr[i];
77401a81e61SBarry Smith   }
77501a81e61SBarry Smith 
77601a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
77701a81e61SBarry Smith     row_indx             = i-rstart;
77801a81e61SBarry Smith     len                  = num_ptr[row_indx];
77901a81e61SBarry Smith     rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int));
78001a81e61SBarry Smith     rows->A[row_indx]    = (double*)malloc(len*sizeof(double));
78101a81e61SBarry Smith   }
78201a81e61SBarry Smith 
78301a81e61SBarry Smith   /* copy the matrix */
78401a81e61SBarry Smith   for (i=rstart; i<rend; i++) {
78501a81e61SBarry Smith     row_indx = i - rstart;
78601a81e61SBarry Smith     ierr     = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
78701a81e61SBarry Smith     for (j=0; j<nz; j++) {
78801a81e61SBarry Smith       col = cols[j];
78901a81e61SBarry Smith       len = rows->len[row_indx]++;
79001a81e61SBarry Smith       rows->ptrs[row_indx][len] = mapping[col];
79101a81e61SBarry Smith       rows->A[row_indx][len]    = vals[j];
79201a81e61SBarry Smith     }
79301a81e61SBarry Smith     rows->slen[row_indx] = rows->len[row_indx];
79401a81e61SBarry Smith     ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr);
79501a81e61SBarry Smith   }
79601a81e61SBarry Smith 
79701a81e61SBarry Smith 
79801a81e61SBarry Smith   /************************************************************/
79901a81e61SBarry Smith   /************** Set up the column structure *****************/
80001a81e61SBarry Smith   /*********************************************************/
80101a81e61SBarry Smith 
80201a81e61SBarry Smith   if (AT) {
80301a81e61SBarry Smith 
80401a81e61SBarry Smith     /* count number of nonzeros in every column */
80501a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
80601a81e61SBarry Smith       ierr = MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
80701a81e61SBarry Smith       ierr = MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
80801a81e61SBarry Smith     }
80901a81e61SBarry Smith 
81001a81e61SBarry Smith     /* allocate buffers */
81101a81e61SBarry Smith     len = 0;
81201a81e61SBarry Smith     for (i=0; i<mnl; i++) {
81301a81e61SBarry Smith       if (len < num_ptr[i]) len = num_ptr[i];
81401a81e61SBarry Smith     }
81501a81e61SBarry Smith 
81601a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
81701a81e61SBarry Smith       row_indx = i-rstart;
81801a81e61SBarry Smith       len      = num_ptr[row_indx];
81901a81e61SBarry Smith       rows->rptrs[row_indx] = (int*)malloc(len*sizeof(int));
82001a81e61SBarry Smith     }
82101a81e61SBarry Smith 
82201a81e61SBarry Smith     /* copy the matrix (i.e., the structure) */
82301a81e61SBarry Smith     for (i=rstart; i<rend; i++) {
82401a81e61SBarry Smith       row_indx = i - rstart;
82501a81e61SBarry Smith       ierr     = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
82601a81e61SBarry Smith       for (j=0; j<nz; j++) {
82701a81e61SBarry Smith         col = cols[j];
82801a81e61SBarry Smith         len = rows->rlen[row_indx]++;
82901a81e61SBarry Smith         rows->rptrs[row_indx][len] = mapping[col];
83001a81e61SBarry Smith       }
83101a81e61SBarry Smith       ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr);
83201a81e61SBarry Smith     }
83301a81e61SBarry Smith   }
83401a81e61SBarry Smith 
83501a81e61SBarry Smith   ierr = PetscFree(num_ptr);CHKERRQ(ierr);
83601a81e61SBarry Smith   ierr = PetscFree(mapping);CHKERRQ(ierr);
83701a81e61SBarry Smith 
83801a81e61SBarry Smith   order_pointers(M);
83901a81e61SBarry Smith   M->maxnz = calc_maxnz(M);
84001a81e61SBarry Smith 
84101a81e61SBarry Smith   *B = M;
84201a81e61SBarry Smith 
84301a81e61SBarry Smith   PetscFunctionReturn(0);
84401a81e61SBarry Smith }
84501a81e61SBarry Smith 
84601a81e61SBarry Smith /**********************************************************************/
84701a81e61SBarry Smith 
84801a81e61SBarry Smith /*
84901a81e61SBarry Smith    Converts from an SPAI matrix B  to a PETSc matrix PB.
85001a81e61SBarry Smith    This assumes that the the SPAI matrix B is stored in
85101a81e61SBarry Smith    COMPRESSED-ROW format.
85201a81e61SBarry Smith */
85301a81e61SBarry Smith #undef __FUNCT__
85401a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat"
85501a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB)
85601a81e61SBarry Smith {
8574b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
85801a81e61SBarry Smith   PetscErrorCode ierr;
85901a81e61SBarry Smith   int            m,n,M,N;
86001a81e61SBarry Smith   int            d_nz,o_nz;
86101a81e61SBarry Smith   int            *d_nnz,*o_nnz;
86201a81e61SBarry Smith   int            i,k,global_row,global_col,first_diag_col,last_diag_col;
86301a81e61SBarry Smith   PetscScalar    val;
86401a81e61SBarry Smith 
86501a81e61SBarry Smith   PetscFunctionBegin;
86601a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
86701a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
86801a81e61SBarry Smith 
86901a81e61SBarry Smith   m = n = B->mnls[rank];
87001a81e61SBarry Smith   d_nz = o_nz = 0;
87101a81e61SBarry Smith 
87201a81e61SBarry Smith   /* Determine preallocation for MatCreateMPIAIJ */
8732a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&d_nnz);CHKERRQ(ierr);
8742a4967abSBarry Smith   ierr = PetscMalloc(m*sizeof(PetscInt),&o_nnz);CHKERRQ(ierr);
87501a81e61SBarry Smith   for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0;
87601a81e61SBarry Smith   first_diag_col = B->start_indices[rank];
87701a81e61SBarry Smith   last_diag_col = first_diag_col + B->mnls[rank];
87801a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
87901a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
88001a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
881db4deed7SKarl Rupp       if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++;
882db4deed7SKarl Rupp       else o_nnz[i]++;
88301a81e61SBarry Smith     }
88401a81e61SBarry Smith   }
88501a81e61SBarry Smith 
88601a81e61SBarry Smith   M = N = B->n;
88701a81e61SBarry Smith   /* Here we only know how to create AIJ format */
888f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,PB);CHKERRQ(ierr);
889f69a0ea3SMatthew Knepley   ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr);
89001a81e61SBarry Smith   ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr);
89101a81e61SBarry Smith   ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr);
89201a81e61SBarry Smith   ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr);
89301a81e61SBarry Smith 
89401a81e61SBarry Smith   for (i=0; i<B->mnls[rank]; i++) {
89501a81e61SBarry Smith     global_row = B->start_indices[rank]+i;
89601a81e61SBarry Smith     for (k=0; k<B->lines->len[i]; k++) {
89701a81e61SBarry Smith       global_col = B->lines->ptrs[i][k];
89801a81e61SBarry Smith       val = B->lines->A[i][k];
89901a81e61SBarry Smith       ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr);
90001a81e61SBarry Smith     }
90101a81e61SBarry Smith   }
90201a81e61SBarry Smith 
90301a81e61SBarry Smith   ierr = PetscFree(d_nnz);CHKERRQ(ierr);
90401a81e61SBarry Smith   ierr = PetscFree(o_nnz);CHKERRQ(ierr);
90501a81e61SBarry Smith 
90601a81e61SBarry Smith   ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90701a81e61SBarry Smith   ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90801a81e61SBarry Smith 
90901a81e61SBarry Smith   PetscFunctionReturn(0);
91001a81e61SBarry Smith }
91101a81e61SBarry Smith 
91201a81e61SBarry Smith /**********************************************************************/
91301a81e61SBarry Smith 
91401a81e61SBarry Smith /*
91501a81e61SBarry Smith    Converts from an SPAI vector v  to a PETSc vec Pv.
91601a81e61SBarry Smith */
91701a81e61SBarry Smith #undef __FUNCT__
91801a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec"
91901a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv)
92001a81e61SBarry Smith {
92101a81e61SBarry Smith   PetscErrorCode ierr;
9224b6b5fe1SBarry Smith   PetscMPIInt    size,rank;
9234b6b5fe1SBarry Smith   int            m,M,i,*mnls,*start_indices,*global_indices;
92401a81e61SBarry Smith 
92501a81e61SBarry Smith   PetscFunctionBegin;
92601a81e61SBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
92701a81e61SBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
92801a81e61SBarry Smith 
92901a81e61SBarry Smith   m = v->mnl;
93001a81e61SBarry Smith   M = v->n;
93101a81e61SBarry Smith 
93201a81e61SBarry Smith 
93301a81e61SBarry Smith   ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr);
93401a81e61SBarry Smith 
93501a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&mnls);CHKERRQ(ierr);
9362a4967abSBarry Smith   ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRQ(ierr);
93701a81e61SBarry Smith 
93801a81e61SBarry Smith   ierr = PetscMalloc(size*sizeof(int),&start_indices);CHKERRQ(ierr);
93901a81e61SBarry Smith   start_indices[0] = 0;
94001a81e61SBarry Smith   for (i=1; i<size; i++)
94101a81e61SBarry Smith     start_indices[i] = start_indices[i-1] +mnls[i-1];
94201a81e61SBarry Smith 
94301a81e61SBarry Smith   ierr = PetscMalloc(v->mnl*sizeof(int),&global_indices);CHKERRQ(ierr);
94401a81e61SBarry Smith   for (i=0; i<v->mnl; i++)
94501a81e61SBarry Smith     global_indices[i] = start_indices[rank] + i;
94601a81e61SBarry Smith 
94701a81e61SBarry Smith   ierr = PetscFree(mnls);CHKERRQ(ierr);
94801a81e61SBarry Smith   ierr = PetscFree(start_indices);CHKERRQ(ierr);
94901a81e61SBarry Smith 
95001a81e61SBarry Smith   ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr);
95101a81e61SBarry Smith   ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr);
95201a81e61SBarry Smith   ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr);
95301a81e61SBarry Smith 
9544b6b5fe1SBarry Smith   ierr = PetscFree(global_indices);CHKERRQ(ierr);
95501a81e61SBarry Smith   PetscFunctionReturn(0);
95601a81e61SBarry Smith }
95701a81e61SBarry Smith 
95801a81e61SBarry Smith 
95901a81e61SBarry Smith 
96001a81e61SBarry Smith 
961