101a81e61SBarry Smith 201a81e61SBarry Smith /* 301a81e61SBarry Smith 3/99 Modified by Stephen Barnard to support SPAI version 3.0 401a81e61SBarry Smith */ 501a81e61SBarry Smith 601a81e61SBarry Smith /* 701a81e61SBarry Smith Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner 801a81e61SBarry Smith Code written by Stephen Barnard. 901a81e61SBarry Smith 1001a81e61SBarry Smith Note: there is some BAD memory bleeding below! 1101a81e61SBarry Smith 1201a81e61SBarry Smith This code needs work 1301a81e61SBarry Smith 1401a81e61SBarry Smith 1) get rid of all memory bleeding 1501a81e61SBarry Smith 2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix 1601a81e61SBarry Smith rather than having the sp flag for PC_SPAI 172a4967abSBarry Smith 3) fix to set the block size based on the matrix block size 1801a81e61SBarry Smith 1901a81e61SBarry Smith */ 20762360ffSBarry Smith #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */ 2101a81e61SBarry Smith 22af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 231c7d2463SJed Brown #include <../src/ksp/pc/impls/spai/petscspai.h> 2401a81e61SBarry Smith 2501a81e61SBarry Smith /* 2601a81e61SBarry Smith These are the SPAI include files 2701a81e61SBarry Smith */ 2801a81e61SBarry Smith EXTERN_C_BEGIN 29b1bda85cSSatish Balay #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */ 30c6db04a5SJed Brown #include <spai.h> 31c6db04a5SJed Brown #include <matrix.h> 3201a81e61SBarry Smith EXTERN_C_END 3301a81e61SBarry Smith 3409573ac7SBarry Smith extern PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**); 3509573ac7SBarry Smith extern PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix*,Mat*); 3609573ac7SBarry Smith extern PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv); 3709573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char*,char*,char*); 3801a81e61SBarry Smith 3901a81e61SBarry Smith typedef struct { 4001a81e61SBarry Smith 4101a81e61SBarry Smith matrix *B; /* matrix in SPAI format */ 4201a81e61SBarry Smith matrix *BT; /* transpose of matrix in SPAI format */ 4301a81e61SBarry Smith matrix *M; /* the approximate inverse in SPAI format */ 4401a81e61SBarry Smith 4501a81e61SBarry Smith Mat PM; /* the approximate inverse PETSc format */ 4601a81e61SBarry Smith 4701a81e61SBarry Smith double epsilon; /* tolerance */ 4801a81e61SBarry Smith int nbsteps; /* max number of "improvement" steps per line */ 4901a81e61SBarry Smith int max; /* max dimensions of is_I, q, etc. */ 5001a81e61SBarry Smith int maxnew; /* max number of new entries per step */ 5101a81e61SBarry Smith int block_size; /* constant block size */ 5201a81e61SBarry Smith int cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */ 5301a81e61SBarry Smith int verbose; /* SPAI prints timing and statistics */ 5401a81e61SBarry Smith 5501a81e61SBarry Smith int sp; /* symmetric nonzero pattern */ 5601a81e61SBarry Smith MPI_Comm comm_spai; /* communicator to be used with spai */ 5701a81e61SBarry Smith } PC_SPAI; 5801a81e61SBarry Smith 5901a81e61SBarry Smith /**********************************************************************/ 6001a81e61SBarry Smith 6101a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc) 6201a81e61SBarry Smith { 6301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 6401a81e61SBarry Smith PetscErrorCode ierr; 6501a81e61SBarry Smith Mat AT; 6601a81e61SBarry Smith 6701a81e61SBarry Smith PetscFunctionBegin; 6801a81e61SBarry Smith init_SPAI(); 6901a81e61SBarry Smith 7001a81e61SBarry Smith if (ispai->sp) { 7101a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr); 7201a81e61SBarry Smith } else { 7301a81e61SBarry Smith /* Use the transpose to get the column nonzero structure. */ 74fc4dec0aSBarry Smith ierr = MatTranspose(pc->pmat,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 7501a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr); 7669a3b875SJed Brown ierr = MatDestroy(&AT);CHKERRQ(ierr); 7701a81e61SBarry Smith } 7801a81e61SBarry Smith 7901a81e61SBarry Smith /* Destroy the transpose */ 8001a81e61SBarry Smith /* Don't know how to do it. PETSc developers? */ 8101a81e61SBarry Smith 8201a81e61SBarry Smith /* construct SPAI preconditioner */ 8301a81e61SBarry Smith /* FILE *messages */ /* file for warning messages */ 8401a81e61SBarry Smith /* double epsilon */ /* tolerance */ 8501a81e61SBarry Smith /* int nbsteps */ /* max number of "improvement" steps per line */ 8601a81e61SBarry Smith /* int max */ /* max dimensions of is_I, q, etc. */ 8701a81e61SBarry Smith /* int maxnew */ /* max number of new entries per step */ 8801a81e61SBarry Smith /* int block_size */ /* block_size == 1 specifies scalar elments 8901a81e61SBarry Smith block_size == n specifies nxn constant-block elements 9001a81e61SBarry Smith block_size == 0 specifies variable-block elements */ 912fa5cd67SKarl Rupp /* int cache_size */ /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */ 9201a81e61SBarry Smith /* int verbose */ /* verbose == 0 specifies that SPAI is silent 9301a81e61SBarry Smith verbose == 1 prints timing and matrix statistics */ 9401a81e61SBarry Smith 9501a81e61SBarry Smith ierr = bspai(ispai->B,&ispai->M, 9601a81e61SBarry Smith stdout, 9701a81e61SBarry Smith ispai->epsilon, 9801a81e61SBarry Smith ispai->nbsteps, 9901a81e61SBarry Smith ispai->max, 10001a81e61SBarry Smith ispai->maxnew, 10101a81e61SBarry Smith ispai->block_size, 10201a81e61SBarry Smith ispai->cache_size, 10301a81e61SBarry Smith ispai->verbose);CHKERRQ(ierr); 10401a81e61SBarry Smith 105ce94432eSBarry Smith ierr = ConvertMatrixToMat(PetscObjectComm((PetscObject)pc),ispai->M,&ispai->PM);CHKERRQ(ierr); 10601a81e61SBarry Smith 10701a81e61SBarry Smith /* free the SPAI matrices */ 10801a81e61SBarry Smith sp_free_matrix(ispai->B); 10901a81e61SBarry Smith sp_free_matrix(ispai->M); 11001a81e61SBarry Smith PetscFunctionReturn(0); 11101a81e61SBarry Smith } 11201a81e61SBarry Smith 11301a81e61SBarry Smith /**********************************************************************/ 11401a81e61SBarry Smith 11501a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y) 11601a81e61SBarry Smith { 11701a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 11801a81e61SBarry Smith PetscErrorCode ierr; 11901a81e61SBarry Smith 12001a81e61SBarry Smith PetscFunctionBegin; 12101a81e61SBarry Smith /* Now using PETSc's multiply */ 12201a81e61SBarry Smith ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr); 12301a81e61SBarry Smith PetscFunctionReturn(0); 12401a81e61SBarry Smith } 12501a81e61SBarry Smith 12648e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_SPAI(PC pc,Mat X,Mat Y) 12748e38a8aSPierre Jolivet { 12848e38a8aSPierre Jolivet PC_SPAI *ispai = (PC_SPAI*)pc->data; 12948e38a8aSPierre Jolivet PetscErrorCode ierr; 13048e38a8aSPierre Jolivet 13148e38a8aSPierre Jolivet PetscFunctionBegin; 13248e38a8aSPierre Jolivet /* Now using PETSc's multiply */ 13348e38a8aSPierre Jolivet ierr = MatMatMult(ispai->PM,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&Y);CHKERRQ(ierr); 13448e38a8aSPierre Jolivet PetscFunctionReturn(0); 13548e38a8aSPierre Jolivet } 13648e38a8aSPierre Jolivet 13701a81e61SBarry Smith /**********************************************************************/ 13801a81e61SBarry Smith 13901a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc) 14001a81e61SBarry Smith { 14101a81e61SBarry Smith PetscErrorCode ierr; 14201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 14301a81e61SBarry Smith 14401a81e61SBarry Smith PetscFunctionBegin; 14569a3b875SJed Brown ierr = MatDestroy(&ispai->PM);CHKERRQ(ierr); 146ffc4695bSBarry Smith ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRMPI(ierr); 147c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 14801a81e61SBarry Smith PetscFunctionReturn(0); 14901a81e61SBarry Smith } 15001a81e61SBarry Smith 15101a81e61SBarry Smith /**********************************************************************/ 15201a81e61SBarry Smith 15301a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer) 15401a81e61SBarry Smith { 15501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 15601a81e61SBarry Smith PetscErrorCode ierr; 157ace3abfcSBarry Smith PetscBool iascii; 15801a81e61SBarry Smith 15901a81e61SBarry Smith PetscFunctionBegin; 160251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 16101a81e61SBarry Smith if (iascii) { 16257622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," epsilon %g\n", (double)ispai->epsilon);CHKERRQ(ierr); 16301a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," nbsteps %d\n", ispai->nbsteps);CHKERRQ(ierr); 16401a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," max %d\n", ispai->max);CHKERRQ(ierr); 16501a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxnew %d\n", ispai->maxnew);CHKERRQ(ierr); 16601a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block_size %d\n",ispai->block_size);CHKERRQ(ierr); 16701a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," cache_size %d\n",ispai->cache_size);CHKERRQ(ierr); 16801a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," verbose %d\n", ispai->verbose);CHKERRQ(ierr); 16901a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," sp %d\n", ispai->sp);CHKERRQ(ierr); 17001a81e61SBarry Smith } 17101a81e61SBarry Smith PetscFunctionReturn(0); 17201a81e61SBarry Smith } 17301a81e61SBarry Smith 174f7a08781SBarry Smith static 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 18301a81e61SBarry Smith /**********************************************************************/ 18401a81e61SBarry Smith 185f7a08781SBarry Smith static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1) 18601a81e61SBarry Smith { 18701a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1885fd66863SKarl Rupp 18901a81e61SBarry Smith PetscFunctionBegin; 19001a81e61SBarry Smith ispai->nbsteps = nbsteps1; 19101a81e61SBarry Smith PetscFunctionReturn(0); 19201a81e61SBarry Smith } 19301a81e61SBarry Smith 19401a81e61SBarry Smith /**********************************************************************/ 19501a81e61SBarry Smith 19601a81e61SBarry Smith /* added 1/7/99 g.h. */ 197f7a08781SBarry Smith static PetscErrorCode PCSPAISetMax_SPAI(PC pc,int max1) 19801a81e61SBarry Smith { 19901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2005fd66863SKarl Rupp 20101a81e61SBarry Smith PetscFunctionBegin; 20201a81e61SBarry Smith ispai->max = max1; 20301a81e61SBarry Smith PetscFunctionReturn(0); 20401a81e61SBarry Smith } 20501a81e61SBarry Smith 20601a81e61SBarry Smith /**********************************************************************/ 20701a81e61SBarry Smith 208f7a08781SBarry Smith static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc,int maxnew1) 20901a81e61SBarry Smith { 21001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2115fd66863SKarl Rupp 21201a81e61SBarry Smith PetscFunctionBegin; 21301a81e61SBarry Smith ispai->maxnew = maxnew1; 21401a81e61SBarry Smith PetscFunctionReturn(0); 21501a81e61SBarry Smith } 21601a81e61SBarry Smith 21701a81e61SBarry Smith /**********************************************************************/ 21801a81e61SBarry Smith 219f7a08781SBarry Smith static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc,int block_size1) 22001a81e61SBarry Smith { 22101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2225fd66863SKarl Rupp 22301a81e61SBarry Smith PetscFunctionBegin; 22401a81e61SBarry Smith ispai->block_size = block_size1; 22501a81e61SBarry Smith PetscFunctionReturn(0); 22601a81e61SBarry Smith } 22701a81e61SBarry Smith 22801a81e61SBarry Smith /**********************************************************************/ 22901a81e61SBarry Smith 230f7a08781SBarry Smith static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc,int cache_size) 23101a81e61SBarry Smith { 23201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2335fd66863SKarl Rupp 23401a81e61SBarry Smith PetscFunctionBegin; 23501a81e61SBarry Smith ispai->cache_size = cache_size; 23601a81e61SBarry Smith PetscFunctionReturn(0); 23701a81e61SBarry Smith } 23801a81e61SBarry Smith 23901a81e61SBarry Smith /**********************************************************************/ 24001a81e61SBarry Smith 241f7a08781SBarry Smith static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc,int verbose) 24201a81e61SBarry Smith { 24301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2445fd66863SKarl Rupp 24501a81e61SBarry Smith PetscFunctionBegin; 24601a81e61SBarry Smith ispai->verbose = verbose; 24701a81e61SBarry Smith PetscFunctionReturn(0); 24801a81e61SBarry Smith } 24901a81e61SBarry Smith 25001a81e61SBarry Smith /**********************************************************************/ 25101a81e61SBarry Smith 252f7a08781SBarry Smith static PetscErrorCode PCSPAISetSp_SPAI(PC pc,int sp) 25301a81e61SBarry Smith { 25401a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2555fd66863SKarl Rupp 25601a81e61SBarry Smith PetscFunctionBegin; 25701a81e61SBarry Smith ispai->sp = sp; 25801a81e61SBarry Smith PetscFunctionReturn(0); 25901a81e61SBarry Smith } 26001a81e61SBarry Smith 26101a81e61SBarry Smith /* -------------------------------------------------------------------*/ 26201a81e61SBarry Smith 26301a81e61SBarry Smith /*@ 26401a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 26501a81e61SBarry Smith 26601a81e61SBarry Smith Input Parameters: 26701a81e61SBarry Smith + pc - the preconditioner 26801a81e61SBarry Smith - eps - epsilon (default .4) 26901a81e61SBarry Smith 27095452b02SPatrick Sanan Notes: 27195452b02SPatrick Sanan Espilon must be between 0 and 1. It controls the 27201a81e61SBarry Smith quality of the approximation of M to the inverse of 27301a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 27401a81e61SBarry Smith fill, and usually better preconditioners. In many 27501a81e61SBarry Smith cases the best choice of epsilon is the one that 27601a81e61SBarry Smith divides the total solution time equally between the 27701a81e61SBarry Smith preconditioner and the solver. 27801a81e61SBarry Smith 27901a81e61SBarry Smith Level: intermediate 28001a81e61SBarry Smith 28101a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 28201a81e61SBarry Smith @*/ 2837087cfbeSBarry Smith PetscErrorCode PCSPAISetEpsilon(PC pc,double epsilon1) 28401a81e61SBarry Smith { 2854ac538c5SBarry Smith PetscErrorCode ierr; 2865fd66863SKarl Rupp 28701a81e61SBarry Smith PetscFunctionBegin; 2884ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));CHKERRQ(ierr); 28901a81e61SBarry Smith PetscFunctionReturn(0); 29001a81e61SBarry Smith } 29101a81e61SBarry Smith 29201a81e61SBarry Smith /**********************************************************************/ 29301a81e61SBarry Smith 29401a81e61SBarry Smith /*@ 29501a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 29601a81e61SBarry Smith the SPAI preconditioner 29701a81e61SBarry Smith 29801a81e61SBarry Smith Input Parameters: 29901a81e61SBarry Smith + pc - the preconditioner 30001a81e61SBarry Smith - n - number of steps (default 5) 30101a81e61SBarry Smith 30295452b02SPatrick Sanan Notes: 30395452b02SPatrick Sanan SPAI constructs to approximation to every column of 30401a81e61SBarry Smith the exact inverse of A in a series of improvement 30501a81e61SBarry Smith steps. The quality of the approximation is determined 30601a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 30701a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 30801a81e61SBarry Smith uses the best approximation constructed so far. 30901a81e61SBarry Smith 31001a81e61SBarry Smith Level: intermediate 31101a81e61SBarry Smith 31201a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew() 31301a81e61SBarry Smith @*/ 3147087cfbeSBarry Smith PetscErrorCode PCSPAISetNBSteps(PC pc,int nbsteps1) 31501a81e61SBarry Smith { 3164ac538c5SBarry Smith PetscErrorCode ierr; 3175fd66863SKarl Rupp 31801a81e61SBarry Smith PetscFunctionBegin; 3194ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));CHKERRQ(ierr); 32001a81e61SBarry Smith PetscFunctionReturn(0); 32101a81e61SBarry Smith } 32201a81e61SBarry Smith 32301a81e61SBarry Smith /**********************************************************************/ 32401a81e61SBarry Smith 32501a81e61SBarry Smith /* added 1/7/99 g.h. */ 32601a81e61SBarry Smith /*@ 32701a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 32801a81e61SBarry Smith the SPAI preconditioner 32901a81e61SBarry Smith 33001a81e61SBarry Smith Input Parameters: 33101a81e61SBarry Smith + pc - the preconditioner 33201a81e61SBarry Smith - n - size (default is 5000) 33301a81e61SBarry Smith 33401a81e61SBarry Smith Level: intermediate 33501a81e61SBarry Smith 33601a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 33701a81e61SBarry Smith @*/ 3387087cfbeSBarry Smith PetscErrorCode PCSPAISetMax(PC pc,int max1) 33901a81e61SBarry Smith { 3404ac538c5SBarry Smith PetscErrorCode ierr; 3415fd66863SKarl Rupp 34201a81e61SBarry Smith PetscFunctionBegin; 3434ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));CHKERRQ(ierr); 34401a81e61SBarry Smith PetscFunctionReturn(0); 34501a81e61SBarry Smith } 34601a81e61SBarry Smith 34701a81e61SBarry Smith /**********************************************************************/ 34801a81e61SBarry Smith 34901a81e61SBarry Smith /*@ 35001a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 35101a81e61SBarry Smith in SPAI preconditioner 35201a81e61SBarry Smith 35301a81e61SBarry Smith Input Parameters: 35401a81e61SBarry Smith + pc - the preconditioner 35501a81e61SBarry Smith - n - maximum number (default 5) 35601a81e61SBarry Smith 35701a81e61SBarry Smith Level: intermediate 35801a81e61SBarry Smith 35901a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps() 36001a81e61SBarry Smith @*/ 3617087cfbeSBarry Smith PetscErrorCode PCSPAISetMaxNew(PC pc,int maxnew1) 36201a81e61SBarry Smith { 3634ac538c5SBarry Smith PetscErrorCode ierr; 3645fd66863SKarl Rupp 36501a81e61SBarry Smith PetscFunctionBegin; 3664ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));CHKERRQ(ierr); 36701a81e61SBarry Smith PetscFunctionReturn(0); 36801a81e61SBarry Smith } 36901a81e61SBarry Smith 37001a81e61SBarry Smith /**********************************************************************/ 37101a81e61SBarry Smith 37201a81e61SBarry Smith /*@ 37301a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 37401a81e61SBarry Smith 37501a81e61SBarry Smith Input Parameters: 37601a81e61SBarry Smith + pc - the preconditioner 37701a81e61SBarry Smith - n - block size (default 1) 37801a81e61SBarry Smith 37995452b02SPatrick Sanan Notes: 38095452b02SPatrick Sanan A block 38101a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 38201a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 38301a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 38401a81e61SBarry Smith variable sized blocks, which are determined by 38501a81e61SBarry Smith searching for dense square diagonal blocks in A. 38601a81e61SBarry Smith This can be very effective for finite-element 38701a81e61SBarry Smith matrices. 38801a81e61SBarry Smith 38901a81e61SBarry Smith SPAI will convert A to block form, use a block 39001a81e61SBarry Smith version of the preconditioner algorithm, and then 39101a81e61SBarry Smith convert the result back to scalar form. 39201a81e61SBarry Smith 39301a81e61SBarry Smith In many cases the a block-size parameter other than 1 39401a81e61SBarry Smith can lead to very significant improvement in 39501a81e61SBarry Smith performance. 39601a81e61SBarry Smith 39701a81e61SBarry Smith 39801a81e61SBarry Smith Level: intermediate 39901a81e61SBarry Smith 40001a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 40101a81e61SBarry Smith @*/ 4027087cfbeSBarry Smith PetscErrorCode PCSPAISetBlockSize(PC pc,int block_size1) 40301a81e61SBarry Smith { 4044ac538c5SBarry Smith PetscErrorCode ierr; 4055fd66863SKarl Rupp 40601a81e61SBarry Smith PetscFunctionBegin; 4074ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));CHKERRQ(ierr); 40801a81e61SBarry Smith PetscFunctionReturn(0); 40901a81e61SBarry Smith } 41001a81e61SBarry Smith 41101a81e61SBarry Smith /**********************************************************************/ 41201a81e61SBarry Smith 41301a81e61SBarry Smith /*@ 41401a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 41501a81e61SBarry Smith 41601a81e61SBarry Smith Input Parameters: 41701a81e61SBarry Smith + pc - the preconditioner 41801a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 41901a81e61SBarry Smith 42095452b02SPatrick Sanan Notes: 42195452b02SPatrick Sanan SPAI uses a hash table to cache messages and avoid 42201a81e61SBarry Smith redundant communication. If suggest always using 42301a81e61SBarry Smith 5. This parameter is irrelevant in the serial 42401a81e61SBarry Smith version. 42501a81e61SBarry Smith 42601a81e61SBarry Smith Level: intermediate 42701a81e61SBarry Smith 42801a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 42901a81e61SBarry Smith @*/ 4307087cfbeSBarry Smith PetscErrorCode PCSPAISetCacheSize(PC pc,int cache_size) 43101a81e61SBarry Smith { 4324ac538c5SBarry Smith PetscErrorCode ierr; 4335fd66863SKarl Rupp 43401a81e61SBarry Smith PetscFunctionBegin; 4354ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));CHKERRQ(ierr); 43601a81e61SBarry Smith PetscFunctionReturn(0); 43701a81e61SBarry Smith } 43801a81e61SBarry Smith 43901a81e61SBarry Smith /**********************************************************************/ 44001a81e61SBarry Smith 44101a81e61SBarry Smith /*@ 44201a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 44301a81e61SBarry Smith 44401a81e61SBarry Smith Input Parameters: 44501a81e61SBarry Smith + pc - the preconditioner 44601a81e61SBarry Smith - n - level (default 1) 44701a81e61SBarry Smith 44895452b02SPatrick Sanan Notes: 44995452b02SPatrick Sanan print parameters, timings and matrix statistics 45001a81e61SBarry Smith 45101a81e61SBarry Smith Level: intermediate 45201a81e61SBarry Smith 45301a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 45401a81e61SBarry Smith @*/ 4557087cfbeSBarry Smith PetscErrorCode PCSPAISetVerbose(PC pc,int verbose) 45601a81e61SBarry Smith { 4574ac538c5SBarry Smith PetscErrorCode ierr; 4585fd66863SKarl Rupp 45901a81e61SBarry Smith PetscFunctionBegin; 4604ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));CHKERRQ(ierr); 46101a81e61SBarry Smith PetscFunctionReturn(0); 46201a81e61SBarry Smith } 46301a81e61SBarry Smith 46401a81e61SBarry Smith /**********************************************************************/ 46501a81e61SBarry Smith 46601a81e61SBarry Smith /*@ 46701a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 46801a81e61SBarry Smith 46901a81e61SBarry Smith Input Parameters: 47001a81e61SBarry Smith + pc - the preconditioner 47101a81e61SBarry Smith - n - 0 or 1 47201a81e61SBarry Smith 47395452b02SPatrick Sanan Notes: 47495452b02SPatrick Sanan If A has a symmetric nonzero pattern use -sp 1 to 47501a81e61SBarry Smith improve performance by eliminating some communication 47601a81e61SBarry Smith in the parallel version. Even if A does not have a 47701a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 47801a81e61SBarry Smith results, but the code will not follow the published 47901a81e61SBarry Smith SPAI algorithm exactly. 48001a81e61SBarry Smith 48101a81e61SBarry Smith 48201a81e61SBarry Smith Level: intermediate 48301a81e61SBarry Smith 48401a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 48501a81e61SBarry Smith @*/ 4867087cfbeSBarry Smith PetscErrorCode PCSPAISetSp(PC pc,int sp) 48701a81e61SBarry Smith { 4884ac538c5SBarry Smith PetscErrorCode ierr; 4895fd66863SKarl Rupp 49001a81e61SBarry Smith PetscFunctionBegin; 4914ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));CHKERRQ(ierr); 49201a81e61SBarry Smith PetscFunctionReturn(0); 49301a81e61SBarry Smith } 49401a81e61SBarry Smith 49501a81e61SBarry Smith /**********************************************************************/ 49601a81e61SBarry Smith 49701a81e61SBarry Smith /**********************************************************************/ 49801a81e61SBarry Smith 4994416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PetscOptionItems *PetscOptionsObject,PC pc) 50001a81e61SBarry Smith { 50101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 50201a81e61SBarry Smith PetscErrorCode ierr; 50301a81e61SBarry Smith int nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp; 50401a81e61SBarry Smith double epsilon1; 505ace3abfcSBarry Smith PetscBool flg; 50601a81e61SBarry Smith 50701a81e61SBarry Smith PetscFunctionBegin; 508e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SPAI options");CHKERRQ(ierr); 50901a81e61SBarry Smith ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr); 51001a81e61SBarry Smith if (flg) { 51101a81e61SBarry Smith ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr); 51201a81e61SBarry Smith } 51301a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr); 51401a81e61SBarry Smith if (flg) { 51501a81e61SBarry Smith ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr); 51601a81e61SBarry Smith } 51701a81e61SBarry Smith /* added 1/7/99 g.h. */ 51801a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr); 51901a81e61SBarry Smith if (flg) { 52001a81e61SBarry Smith ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr); 52101a81e61SBarry Smith } 52201a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr); 52301a81e61SBarry Smith if (flg) { 52401a81e61SBarry Smith ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr); 52501a81e61SBarry Smith } 52601a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr); 52701a81e61SBarry Smith if (flg) { 52801a81e61SBarry Smith ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr); 52901a81e61SBarry Smith } 53001a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr); 53101a81e61SBarry Smith if (flg) { 53201a81e61SBarry Smith ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr); 53301a81e61SBarry Smith } 53401a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr); 53501a81e61SBarry Smith if (flg) { 53601a81e61SBarry Smith ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr); 53701a81e61SBarry Smith } 53801a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr); 53901a81e61SBarry Smith if (flg) { 54001a81e61SBarry Smith ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr); 54101a81e61SBarry Smith } 54201a81e61SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 54301a81e61SBarry Smith PetscFunctionReturn(0); 54401a81e61SBarry Smith } 54501a81e61SBarry Smith 54601a81e61SBarry Smith /**********************************************************************/ 54701a81e61SBarry Smith 54801a81e61SBarry Smith /*MC 54901a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 55001a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 55101a81e61SBarry Smith 55201a81e61SBarry Smith Options Database Keys: 553c5ae99e3SSatish Balay + -pc_spai_epsilon <eps> - set tolerance 55401a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 55501a81e61SBarry Smith . -pc_spai_max <m> - set max 55601a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 55701a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 55801a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 55901a81e61SBarry Smith . -pc_spai_sp <m> - set sp 56001a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 56101a81e61SBarry Smith 56295452b02SPatrick Sanan Notes: 56395452b02SPatrick Sanan This only works with AIJ matrices. 56401a81e61SBarry Smith 56501a81e61SBarry Smith Level: beginner 56601a81e61SBarry Smith 56701a81e61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 56801a81e61SBarry Smith PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(), 56901a81e61SBarry Smith PCSPAISetVerbose(), PCSPAISetSp() 57001a81e61SBarry Smith M*/ 57101a81e61SBarry Smith 5728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc) 57301a81e61SBarry Smith { 57401a81e61SBarry Smith PC_SPAI *ispai; 57501a81e61SBarry Smith PetscErrorCode ierr; 57601a81e61SBarry Smith 57701a81e61SBarry Smith PetscFunctionBegin; 578b00a9115SJed Brown ierr = PetscNewLog(pc,&ispai);CHKERRQ(ierr); 5792a4967abSBarry Smith pc->data = ispai; 58001a81e61SBarry Smith 58101a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 58201a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 58348e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_SPAI; 58401a81e61SBarry Smith pc->ops->applyrichardson = 0; 58501a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 58601a81e61SBarry Smith pc->ops->view = PCView_SPAI; 58701a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 58801a81e61SBarry Smith 58901a81e61SBarry Smith ispai->epsilon = .4; 59001a81e61SBarry Smith ispai->nbsteps = 5; 59101a81e61SBarry Smith ispai->max = 5000; 59201a81e61SBarry Smith ispai->maxnew = 5; 59301a81e61SBarry Smith ispai->block_size = 1; 59401a81e61SBarry Smith ispai->cache_size = 5; 59501a81e61SBarry Smith ispai->verbose = 0; 59601a81e61SBarry Smith 59701a81e61SBarry Smith ispai->sp = 1; 59855b25c41SPierre Jolivet ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai));CHKERRMPI(ierr); 59901a81e61SBarry Smith 600bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",PCSPAISetEpsilon_SPAI);CHKERRQ(ierr); 601bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",PCSPAISetNBSteps_SPAI);CHKERRQ(ierr); 602bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",PCSPAISetMax_SPAI);CHKERRQ(ierr); 603bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",PCSPAISetMaxNew_SPAI);CHKERRQ(ierr); 604bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",PCSPAISetBlockSize_SPAI);CHKERRQ(ierr); 605bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",PCSPAISetCacheSize_SPAI);CHKERRQ(ierr); 606bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",PCSPAISetVerbose_SPAI);CHKERRQ(ierr); 607bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",PCSPAISetSp_SPAI);CHKERRQ(ierr); 60801a81e61SBarry Smith PetscFunctionReturn(0); 60901a81e61SBarry Smith } 61001a81e61SBarry Smith 61101a81e61SBarry Smith /**********************************************************************/ 61201a81e61SBarry Smith 61301a81e61SBarry Smith /* 61401a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 61501a81e61SBarry Smith */ 61601a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B) 61701a81e61SBarry Smith { 61801a81e61SBarry Smith matrix *M; 61901a81e61SBarry Smith int i,j,col; 62001a81e61SBarry Smith int row_indx; 62101a81e61SBarry Smith int len,pe,local_indx,start_indx; 62201a81e61SBarry Smith int *mapping; 62301a81e61SBarry Smith PetscErrorCode ierr; 62401a81e61SBarry Smith const int *cols; 62501a81e61SBarry Smith const double *vals; 6262122902bSSatish Balay int n,mnl,nnl,nz,rstart,rend; 6272a4967abSBarry Smith PetscMPIInt size,rank; 62801a81e61SBarry Smith struct compressed_lines *rows; 62901a81e61SBarry Smith 63001a81e61SBarry Smith PetscFunctionBegin; 631ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 632ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 63301a81e61SBarry Smith ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr); 63401a81e61SBarry Smith ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr); 63501a81e61SBarry Smith 63601a81e61SBarry Smith /* 63701a81e61SBarry Smith not sure why a barrier is required. commenting out 638ffc4695bSBarry Smith ierr = MPI_Barrier(comm);CHKERRMPI(ierr); 63901a81e61SBarry Smith */ 64001a81e61SBarry Smith 64129b2f704SSatish Balay M = new_matrix((SPAI_Comm)comm); 64201a81e61SBarry Smith 64301a81e61SBarry Smith M->n = n; 64401a81e61SBarry Smith M->bs = 1; 64501a81e61SBarry Smith M->max_block_size = 1; 64601a81e61SBarry Smith 64701a81e61SBarry Smith M->mnls = (int*)malloc(sizeof(int)*size); 64801a81e61SBarry Smith M->start_indices = (int*)malloc(sizeof(int)*size); 64901a81e61SBarry Smith M->pe = (int*)malloc(sizeof(int)*n); 65001a81e61SBarry Smith M->block_sizes = (int*)malloc(sizeof(int)*n); 65101a81e61SBarry Smith for (i=0; i<n; i++) M->block_sizes[i] = 1; 65201a81e61SBarry Smith 653ffc4695bSBarry Smith ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRMPI(ierr); 65401a81e61SBarry Smith 65501a81e61SBarry Smith M->start_indices[0] = 0; 6562fa5cd67SKarl Rupp for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1]; 65701a81e61SBarry Smith 65801a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 65901a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 66001a81e61SBarry Smith 66101a81e61SBarry Smith for (i=0; i<size; i++) { 66201a81e61SBarry Smith start_indx = M->start_indices[i]; 6632fa5cd67SKarl Rupp for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i; 66401a81e61SBarry Smith } 66501a81e61SBarry Smith 66601a81e61SBarry Smith if (AT) { 66701a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr); 66801a81e61SBarry Smith } else { 66901a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr); 67001a81e61SBarry Smith } 67101a81e61SBarry Smith 67201a81e61SBarry Smith rows = M->lines; 67301a81e61SBarry Smith 67401a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 675785e854fSJed Brown ierr = PetscMalloc1(M->n,&mapping);CHKERRQ(ierr); 67601a81e61SBarry Smith pe = 0; 67701a81e61SBarry Smith local_indx = 0; 67801a81e61SBarry Smith for (i=0; i<M->n; i++) { 67901a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 68001a81e61SBarry Smith pe++; 68101a81e61SBarry Smith local_indx = 0; 68201a81e61SBarry Smith } 68301a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 68401a81e61SBarry Smith local_indx++; 68501a81e61SBarry Smith } 68601a81e61SBarry Smith 68701a81e61SBarry Smith /*********************************************************/ 68801a81e61SBarry Smith /************** Set up the row structure *****************/ 68901a81e61SBarry Smith /*********************************************************/ 69001a81e61SBarry Smith 69101a81e61SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 69201a81e61SBarry Smith for (i=rstart; i<rend; i++) { 69301a81e61SBarry Smith row_indx = i - rstart; 69401a81e61SBarry Smith ierr = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 6952122902bSSatish Balay /* allocate buffers */ 6962122902bSSatish Balay rows->ptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 6972122902bSSatish Balay rows->A[row_indx] = (double*)malloc(nz*sizeof(double)); 6982122902bSSatish Balay /* copy the matrix */ 69901a81e61SBarry Smith for (j=0; j<nz; j++) { 70001a81e61SBarry Smith col = cols[j]; 70101a81e61SBarry Smith len = rows->len[row_indx]++; 7022fa5cd67SKarl Rupp 70301a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 70401a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 70501a81e61SBarry Smith } 70601a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 7072fa5cd67SKarl Rupp 70801a81e61SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 70901a81e61SBarry Smith } 71001a81e61SBarry Smith 71101a81e61SBarry Smith 71201a81e61SBarry Smith /************************************************************/ 71301a81e61SBarry Smith /************** Set up the column structure *****************/ 71401a81e61SBarry Smith /*********************************************************/ 71501a81e61SBarry Smith 71601a81e61SBarry Smith if (AT) { 71701a81e61SBarry Smith 71801a81e61SBarry Smith for (i=rstart; i<rend; i++) { 71901a81e61SBarry Smith row_indx = i - rstart; 72001a81e61SBarry Smith ierr = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 7212122902bSSatish Balay /* allocate buffers */ 7222122902bSSatish Balay rows->rptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 7232122902bSSatish Balay /* copy the matrix (i.e., the structure) */ 72401a81e61SBarry Smith for (j=0; j<nz; j++) { 72501a81e61SBarry Smith col = cols[j]; 72601a81e61SBarry Smith len = rows->rlen[row_indx]++; 7272fa5cd67SKarl Rupp 72801a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 72901a81e61SBarry Smith } 73001a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 73101a81e61SBarry Smith } 73201a81e61SBarry Smith } 73301a81e61SBarry Smith 73401a81e61SBarry Smith ierr = PetscFree(mapping);CHKERRQ(ierr); 73501a81e61SBarry Smith 73601a81e61SBarry Smith order_pointers(M); 73701a81e61SBarry Smith M->maxnz = calc_maxnz(M); 73801a81e61SBarry Smith *B = M; 73901a81e61SBarry Smith PetscFunctionReturn(0); 74001a81e61SBarry Smith } 74101a81e61SBarry Smith 74201a81e61SBarry Smith /**********************************************************************/ 74301a81e61SBarry Smith 74401a81e61SBarry Smith /* 74501a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 746df3898eeSBarry Smith This assumes that the SPAI matrix B is stored in 74701a81e61SBarry Smith COMPRESSED-ROW format. 74801a81e61SBarry Smith */ 74901a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB) 75001a81e61SBarry Smith { 7514b6b5fe1SBarry Smith PetscMPIInt size,rank; 75201a81e61SBarry Smith PetscErrorCode ierr; 75301a81e61SBarry Smith int m,n,M,N; 75401a81e61SBarry Smith int d_nz,o_nz; 75501a81e61SBarry Smith int *d_nnz,*o_nnz; 75601a81e61SBarry Smith int i,k,global_row,global_col,first_diag_col,last_diag_col; 75701a81e61SBarry Smith PetscScalar val; 75801a81e61SBarry Smith 75901a81e61SBarry Smith PetscFunctionBegin; 760ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 761ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 76201a81e61SBarry Smith 76301a81e61SBarry Smith m = n = B->mnls[rank]; 76401a81e61SBarry Smith d_nz = o_nz = 0; 76501a81e61SBarry Smith 766*06946f3aSJose E. Roman /* Determine preallocation for MatCreateAIJ */ 767785e854fSJed Brown ierr = PetscMalloc1(m,&d_nnz);CHKERRQ(ierr); 768785e854fSJed Brown ierr = PetscMalloc1(m,&o_nnz);CHKERRQ(ierr); 76901a81e61SBarry Smith for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0; 77001a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 77101a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 77201a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 77301a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 77401a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 775db4deed7SKarl Rupp if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++; 776db4deed7SKarl Rupp else o_nnz[i]++; 77701a81e61SBarry Smith } 77801a81e61SBarry Smith } 77901a81e61SBarry Smith 78001a81e61SBarry Smith M = N = B->n; 78101a81e61SBarry Smith /* Here we only know how to create AIJ format */ 782f69a0ea3SMatthew Knepley ierr = MatCreate(comm,PB);CHKERRQ(ierr); 783f69a0ea3SMatthew Knepley ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr); 78401a81e61SBarry Smith ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr); 78501a81e61SBarry Smith ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr); 78601a81e61SBarry Smith ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr); 78701a81e61SBarry Smith 78801a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 78901a81e61SBarry Smith global_row = B->start_indices[rank]+i; 79001a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 79101a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 7922fa5cd67SKarl Rupp 79301a81e61SBarry Smith val = B->lines->A[i][k]; 79401a81e61SBarry Smith ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr); 79501a81e61SBarry Smith } 79601a81e61SBarry Smith } 79701a81e61SBarry Smith 79801a81e61SBarry Smith ierr = PetscFree(d_nnz);CHKERRQ(ierr); 79901a81e61SBarry Smith ierr = PetscFree(o_nnz);CHKERRQ(ierr); 80001a81e61SBarry Smith 80101a81e61SBarry Smith ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 80201a81e61SBarry Smith ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 80301a81e61SBarry Smith PetscFunctionReturn(0); 80401a81e61SBarry Smith } 80501a81e61SBarry Smith 80601a81e61SBarry Smith /**********************************************************************/ 80701a81e61SBarry Smith 80801a81e61SBarry Smith /* 80901a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 81001a81e61SBarry Smith */ 81101a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv) 81201a81e61SBarry Smith { 81301a81e61SBarry Smith PetscErrorCode ierr; 8144b6b5fe1SBarry Smith PetscMPIInt size,rank; 8154b6b5fe1SBarry Smith int m,M,i,*mnls,*start_indices,*global_indices; 81601a81e61SBarry Smith 81701a81e61SBarry Smith PetscFunctionBegin; 818ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 819ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 82001a81e61SBarry Smith 82101a81e61SBarry Smith m = v->mnl; 82201a81e61SBarry Smith M = v->n; 82301a81e61SBarry Smith 82401a81e61SBarry Smith 82501a81e61SBarry Smith ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr); 82601a81e61SBarry Smith 827785e854fSJed Brown ierr = PetscMalloc1(size,&mnls);CHKERRQ(ierr); 828ffc4695bSBarry Smith ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRMPI(ierr); 82901a81e61SBarry Smith 830785e854fSJed Brown ierr = PetscMalloc1(size,&start_indices);CHKERRQ(ierr); 8312fa5cd67SKarl Rupp 83201a81e61SBarry Smith start_indices[0] = 0; 8332fa5cd67SKarl Rupp for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1]; 83401a81e61SBarry Smith 835785e854fSJed Brown ierr = PetscMalloc1(v->mnl,&global_indices);CHKERRQ(ierr); 8362fa5cd67SKarl Rupp for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i; 83701a81e61SBarry Smith 83801a81e61SBarry Smith ierr = PetscFree(mnls);CHKERRQ(ierr); 83901a81e61SBarry Smith ierr = PetscFree(start_indices);CHKERRQ(ierr); 84001a81e61SBarry Smith 84101a81e61SBarry Smith ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr); 84201a81e61SBarry Smith ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr); 84301a81e61SBarry Smith ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr); 84401a81e61SBarry Smith 8454b6b5fe1SBarry Smith ierr = PetscFree(global_indices);CHKERRQ(ierr); 84601a81e61SBarry Smith PetscFunctionReturn(0); 84701a81e61SBarry Smith } 848