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 12601a81e61SBarry Smith /**********************************************************************/ 12701a81e61SBarry Smith 12801a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc) 12901a81e61SBarry Smith { 13001a81e61SBarry Smith PetscErrorCode ierr; 13101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 13201a81e61SBarry Smith 13301a81e61SBarry Smith PetscFunctionBegin; 13469a3b875SJed Brown ierr = MatDestroy(&ispai->PM);CHKERRQ(ierr); 13501a81e61SBarry Smith ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr); 136c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 13701a81e61SBarry Smith PetscFunctionReturn(0); 13801a81e61SBarry Smith } 13901a81e61SBarry Smith 14001a81e61SBarry Smith /**********************************************************************/ 14101a81e61SBarry Smith 14201a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer) 14301a81e61SBarry Smith { 14401a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 14501a81e61SBarry Smith PetscErrorCode ierr; 146ace3abfcSBarry Smith PetscBool iascii; 14701a81e61SBarry Smith 14801a81e61SBarry Smith PetscFunctionBegin; 149251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 15001a81e61SBarry Smith if (iascii) { 15157622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," epsilon %g\n", (double)ispai->epsilon);CHKERRQ(ierr); 15201a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," nbsteps %d\n", ispai->nbsteps);CHKERRQ(ierr); 15301a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," max %d\n", ispai->max);CHKERRQ(ierr); 15401a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxnew %d\n", ispai->maxnew);CHKERRQ(ierr); 15501a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block_size %d\n",ispai->block_size);CHKERRQ(ierr); 15601a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," cache_size %d\n",ispai->cache_size);CHKERRQ(ierr); 15701a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," verbose %d\n", ispai->verbose);CHKERRQ(ierr); 15801a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," sp %d\n", ispai->sp);CHKERRQ(ierr); 15901a81e61SBarry Smith } 16001a81e61SBarry Smith PetscFunctionReturn(0); 16101a81e61SBarry Smith } 16201a81e61SBarry Smith 163f7a08781SBarry Smith static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc,double epsilon1) 16401a81e61SBarry Smith { 16501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1665fd66863SKarl Rupp 16701a81e61SBarry Smith PetscFunctionBegin; 16801a81e61SBarry Smith ispai->epsilon = epsilon1; 16901a81e61SBarry Smith PetscFunctionReturn(0); 17001a81e61SBarry Smith } 17101a81e61SBarry Smith 17201a81e61SBarry Smith /**********************************************************************/ 17301a81e61SBarry Smith 174f7a08781SBarry Smith static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1) 17501a81e61SBarry Smith { 17601a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1775fd66863SKarl Rupp 17801a81e61SBarry Smith PetscFunctionBegin; 17901a81e61SBarry Smith ispai->nbsteps = nbsteps1; 18001a81e61SBarry Smith PetscFunctionReturn(0); 18101a81e61SBarry Smith } 18201a81e61SBarry Smith 18301a81e61SBarry Smith /**********************************************************************/ 18401a81e61SBarry Smith 18501a81e61SBarry Smith /* added 1/7/99 g.h. */ 186f7a08781SBarry Smith static PetscErrorCode PCSPAISetMax_SPAI(PC pc,int max1) 18701a81e61SBarry Smith { 18801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1895fd66863SKarl Rupp 19001a81e61SBarry Smith PetscFunctionBegin; 19101a81e61SBarry Smith ispai->max = max1; 19201a81e61SBarry Smith PetscFunctionReturn(0); 19301a81e61SBarry Smith } 19401a81e61SBarry Smith 19501a81e61SBarry Smith /**********************************************************************/ 19601a81e61SBarry Smith 197f7a08781SBarry Smith static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc,int maxnew1) 19801a81e61SBarry Smith { 19901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2005fd66863SKarl Rupp 20101a81e61SBarry Smith PetscFunctionBegin; 20201a81e61SBarry Smith ispai->maxnew = maxnew1; 20301a81e61SBarry Smith PetscFunctionReturn(0); 20401a81e61SBarry Smith } 20501a81e61SBarry Smith 20601a81e61SBarry Smith /**********************************************************************/ 20701a81e61SBarry Smith 208f7a08781SBarry Smith static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc,int block_size1) 20901a81e61SBarry Smith { 21001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2115fd66863SKarl Rupp 21201a81e61SBarry Smith PetscFunctionBegin; 21301a81e61SBarry Smith ispai->block_size = block_size1; 21401a81e61SBarry Smith PetscFunctionReturn(0); 21501a81e61SBarry Smith } 21601a81e61SBarry Smith 21701a81e61SBarry Smith /**********************************************************************/ 21801a81e61SBarry Smith 219f7a08781SBarry Smith static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc,int cache_size) 22001a81e61SBarry Smith { 22101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2225fd66863SKarl Rupp 22301a81e61SBarry Smith PetscFunctionBegin; 22401a81e61SBarry Smith ispai->cache_size = cache_size; 22501a81e61SBarry Smith PetscFunctionReturn(0); 22601a81e61SBarry Smith } 22701a81e61SBarry Smith 22801a81e61SBarry Smith /**********************************************************************/ 22901a81e61SBarry Smith 230f7a08781SBarry Smith static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc,int verbose) 23101a81e61SBarry Smith { 23201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2335fd66863SKarl Rupp 23401a81e61SBarry Smith PetscFunctionBegin; 23501a81e61SBarry Smith ispai->verbose = verbose; 23601a81e61SBarry Smith PetscFunctionReturn(0); 23701a81e61SBarry Smith } 23801a81e61SBarry Smith 23901a81e61SBarry Smith /**********************************************************************/ 24001a81e61SBarry Smith 241f7a08781SBarry Smith static PetscErrorCode PCSPAISetSp_SPAI(PC pc,int sp) 24201a81e61SBarry Smith { 24301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2445fd66863SKarl Rupp 24501a81e61SBarry Smith PetscFunctionBegin; 24601a81e61SBarry Smith ispai->sp = sp; 24701a81e61SBarry Smith PetscFunctionReturn(0); 24801a81e61SBarry Smith } 24901a81e61SBarry Smith 25001a81e61SBarry Smith /* -------------------------------------------------------------------*/ 25101a81e61SBarry Smith 25201a81e61SBarry Smith /*@ 25301a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 25401a81e61SBarry Smith 25501a81e61SBarry Smith Input Parameters: 25601a81e61SBarry Smith + pc - the preconditioner 25701a81e61SBarry Smith - eps - epsilon (default .4) 25801a81e61SBarry Smith 259*95452b02SPatrick Sanan Notes: 260*95452b02SPatrick Sanan Espilon must be between 0 and 1. It controls the 26101a81e61SBarry Smith quality of the approximation of M to the inverse of 26201a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 26301a81e61SBarry Smith fill, and usually better preconditioners. In many 26401a81e61SBarry Smith cases the best choice of epsilon is the one that 26501a81e61SBarry Smith divides the total solution time equally between the 26601a81e61SBarry Smith preconditioner and the solver. 26701a81e61SBarry Smith 26801a81e61SBarry Smith Level: intermediate 26901a81e61SBarry Smith 27001a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 27101a81e61SBarry Smith @*/ 2727087cfbeSBarry Smith PetscErrorCode PCSPAISetEpsilon(PC pc,double epsilon1) 27301a81e61SBarry Smith { 2744ac538c5SBarry Smith PetscErrorCode ierr; 2755fd66863SKarl Rupp 27601a81e61SBarry Smith PetscFunctionBegin; 2774ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));CHKERRQ(ierr); 27801a81e61SBarry Smith PetscFunctionReturn(0); 27901a81e61SBarry Smith } 28001a81e61SBarry Smith 28101a81e61SBarry Smith /**********************************************************************/ 28201a81e61SBarry Smith 28301a81e61SBarry Smith /*@ 28401a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 28501a81e61SBarry Smith the SPAI preconditioner 28601a81e61SBarry Smith 28701a81e61SBarry Smith Input Parameters: 28801a81e61SBarry Smith + pc - the preconditioner 28901a81e61SBarry Smith - n - number of steps (default 5) 29001a81e61SBarry Smith 291*95452b02SPatrick Sanan Notes: 292*95452b02SPatrick Sanan SPAI constructs to approximation to every column of 29301a81e61SBarry Smith the exact inverse of A in a series of improvement 29401a81e61SBarry Smith steps. The quality of the approximation is determined 29501a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 29601a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 29701a81e61SBarry Smith uses the best approximation constructed so far. 29801a81e61SBarry Smith 29901a81e61SBarry Smith Level: intermediate 30001a81e61SBarry Smith 30101a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew() 30201a81e61SBarry Smith @*/ 3037087cfbeSBarry Smith PetscErrorCode PCSPAISetNBSteps(PC pc,int nbsteps1) 30401a81e61SBarry Smith { 3054ac538c5SBarry Smith PetscErrorCode ierr; 3065fd66863SKarl Rupp 30701a81e61SBarry Smith PetscFunctionBegin; 3084ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));CHKERRQ(ierr); 30901a81e61SBarry Smith PetscFunctionReturn(0); 31001a81e61SBarry Smith } 31101a81e61SBarry Smith 31201a81e61SBarry Smith /**********************************************************************/ 31301a81e61SBarry Smith 31401a81e61SBarry Smith /* added 1/7/99 g.h. */ 31501a81e61SBarry Smith /*@ 31601a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 31701a81e61SBarry Smith the SPAI preconditioner 31801a81e61SBarry Smith 31901a81e61SBarry Smith Input Parameters: 32001a81e61SBarry Smith + pc - the preconditioner 32101a81e61SBarry Smith - n - size (default is 5000) 32201a81e61SBarry Smith 32301a81e61SBarry Smith Level: intermediate 32401a81e61SBarry Smith 32501a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 32601a81e61SBarry Smith @*/ 3277087cfbeSBarry Smith PetscErrorCode PCSPAISetMax(PC pc,int max1) 32801a81e61SBarry Smith { 3294ac538c5SBarry Smith PetscErrorCode ierr; 3305fd66863SKarl Rupp 33101a81e61SBarry Smith PetscFunctionBegin; 3324ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));CHKERRQ(ierr); 33301a81e61SBarry Smith PetscFunctionReturn(0); 33401a81e61SBarry Smith } 33501a81e61SBarry Smith 33601a81e61SBarry Smith /**********************************************************************/ 33701a81e61SBarry Smith 33801a81e61SBarry Smith /*@ 33901a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 34001a81e61SBarry Smith in SPAI preconditioner 34101a81e61SBarry Smith 34201a81e61SBarry Smith Input Parameters: 34301a81e61SBarry Smith + pc - the preconditioner 34401a81e61SBarry Smith - n - maximum number (default 5) 34501a81e61SBarry Smith 34601a81e61SBarry Smith Level: intermediate 34701a81e61SBarry Smith 34801a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps() 34901a81e61SBarry Smith @*/ 3507087cfbeSBarry Smith PetscErrorCode PCSPAISetMaxNew(PC pc,int maxnew1) 35101a81e61SBarry Smith { 3524ac538c5SBarry Smith PetscErrorCode ierr; 3535fd66863SKarl Rupp 35401a81e61SBarry Smith PetscFunctionBegin; 3554ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));CHKERRQ(ierr); 35601a81e61SBarry Smith PetscFunctionReturn(0); 35701a81e61SBarry Smith } 35801a81e61SBarry Smith 35901a81e61SBarry Smith /**********************************************************************/ 36001a81e61SBarry Smith 36101a81e61SBarry Smith /*@ 36201a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 36301a81e61SBarry Smith 36401a81e61SBarry Smith Input Parameters: 36501a81e61SBarry Smith + pc - the preconditioner 36601a81e61SBarry Smith - n - block size (default 1) 36701a81e61SBarry Smith 368*95452b02SPatrick Sanan Notes: 369*95452b02SPatrick Sanan A block 37001a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 37101a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 37201a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 37301a81e61SBarry Smith variable sized blocks, which are determined by 37401a81e61SBarry Smith searching for dense square diagonal blocks in A. 37501a81e61SBarry Smith This can be very effective for finite-element 37601a81e61SBarry Smith matrices. 37701a81e61SBarry Smith 37801a81e61SBarry Smith SPAI will convert A to block form, use a block 37901a81e61SBarry Smith version of the preconditioner algorithm, and then 38001a81e61SBarry Smith convert the result back to scalar form. 38101a81e61SBarry Smith 38201a81e61SBarry Smith In many cases the a block-size parameter other than 1 38301a81e61SBarry Smith can lead to very significant improvement in 38401a81e61SBarry Smith performance. 38501a81e61SBarry Smith 38601a81e61SBarry Smith 38701a81e61SBarry Smith Level: intermediate 38801a81e61SBarry Smith 38901a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 39001a81e61SBarry Smith @*/ 3917087cfbeSBarry Smith PetscErrorCode PCSPAISetBlockSize(PC pc,int block_size1) 39201a81e61SBarry Smith { 3934ac538c5SBarry Smith PetscErrorCode ierr; 3945fd66863SKarl Rupp 39501a81e61SBarry Smith PetscFunctionBegin; 3964ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));CHKERRQ(ierr); 39701a81e61SBarry Smith PetscFunctionReturn(0); 39801a81e61SBarry Smith } 39901a81e61SBarry Smith 40001a81e61SBarry Smith /**********************************************************************/ 40101a81e61SBarry Smith 40201a81e61SBarry Smith /*@ 40301a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 40401a81e61SBarry Smith 40501a81e61SBarry Smith Input Parameters: 40601a81e61SBarry Smith + pc - the preconditioner 40701a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 40801a81e61SBarry Smith 409*95452b02SPatrick Sanan Notes: 410*95452b02SPatrick Sanan SPAI uses a hash table to cache messages and avoid 41101a81e61SBarry Smith redundant communication. If suggest always using 41201a81e61SBarry Smith 5. This parameter is irrelevant in the serial 41301a81e61SBarry Smith version. 41401a81e61SBarry Smith 41501a81e61SBarry Smith Level: intermediate 41601a81e61SBarry Smith 41701a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 41801a81e61SBarry Smith @*/ 4197087cfbeSBarry Smith PetscErrorCode PCSPAISetCacheSize(PC pc,int cache_size) 42001a81e61SBarry Smith { 4214ac538c5SBarry Smith PetscErrorCode ierr; 4225fd66863SKarl Rupp 42301a81e61SBarry Smith PetscFunctionBegin; 4244ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));CHKERRQ(ierr); 42501a81e61SBarry Smith PetscFunctionReturn(0); 42601a81e61SBarry Smith } 42701a81e61SBarry Smith 42801a81e61SBarry Smith /**********************************************************************/ 42901a81e61SBarry Smith 43001a81e61SBarry Smith /*@ 43101a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 43201a81e61SBarry Smith 43301a81e61SBarry Smith Input Parameters: 43401a81e61SBarry Smith + pc - the preconditioner 43501a81e61SBarry Smith - n - level (default 1) 43601a81e61SBarry Smith 437*95452b02SPatrick Sanan Notes: 438*95452b02SPatrick Sanan print parameters, timings and matrix statistics 43901a81e61SBarry Smith 44001a81e61SBarry Smith Level: intermediate 44101a81e61SBarry Smith 44201a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 44301a81e61SBarry Smith @*/ 4447087cfbeSBarry Smith PetscErrorCode PCSPAISetVerbose(PC pc,int verbose) 44501a81e61SBarry Smith { 4464ac538c5SBarry Smith PetscErrorCode ierr; 4475fd66863SKarl Rupp 44801a81e61SBarry Smith PetscFunctionBegin; 4494ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));CHKERRQ(ierr); 45001a81e61SBarry Smith PetscFunctionReturn(0); 45101a81e61SBarry Smith } 45201a81e61SBarry Smith 45301a81e61SBarry Smith /**********************************************************************/ 45401a81e61SBarry Smith 45501a81e61SBarry Smith /*@ 45601a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 45701a81e61SBarry Smith 45801a81e61SBarry Smith Input Parameters: 45901a81e61SBarry Smith + pc - the preconditioner 46001a81e61SBarry Smith - n - 0 or 1 46101a81e61SBarry Smith 462*95452b02SPatrick Sanan Notes: 463*95452b02SPatrick Sanan If A has a symmetric nonzero pattern use -sp 1 to 46401a81e61SBarry Smith improve performance by eliminating some communication 46501a81e61SBarry Smith in the parallel version. Even if A does not have a 46601a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 46701a81e61SBarry Smith results, but the code will not follow the published 46801a81e61SBarry Smith SPAI algorithm exactly. 46901a81e61SBarry Smith 47001a81e61SBarry Smith 47101a81e61SBarry Smith Level: intermediate 47201a81e61SBarry Smith 47301a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 47401a81e61SBarry Smith @*/ 4757087cfbeSBarry Smith PetscErrorCode PCSPAISetSp(PC pc,int sp) 47601a81e61SBarry Smith { 4774ac538c5SBarry Smith PetscErrorCode ierr; 4785fd66863SKarl Rupp 47901a81e61SBarry Smith PetscFunctionBegin; 4804ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));CHKERRQ(ierr); 48101a81e61SBarry Smith PetscFunctionReturn(0); 48201a81e61SBarry Smith } 48301a81e61SBarry Smith 48401a81e61SBarry Smith /**********************************************************************/ 48501a81e61SBarry Smith 48601a81e61SBarry Smith /**********************************************************************/ 48701a81e61SBarry Smith 4884416b707SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PetscOptionItems *PetscOptionsObject,PC pc) 48901a81e61SBarry Smith { 49001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 49101a81e61SBarry Smith PetscErrorCode ierr; 49201a81e61SBarry Smith int nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp; 49301a81e61SBarry Smith double epsilon1; 494ace3abfcSBarry Smith PetscBool flg; 49501a81e61SBarry Smith 49601a81e61SBarry Smith PetscFunctionBegin; 497e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SPAI options");CHKERRQ(ierr); 49801a81e61SBarry Smith ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr); 49901a81e61SBarry Smith if (flg) { 50001a81e61SBarry Smith ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr); 50101a81e61SBarry Smith } 50201a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr); 50301a81e61SBarry Smith if (flg) { 50401a81e61SBarry Smith ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr); 50501a81e61SBarry Smith } 50601a81e61SBarry Smith /* added 1/7/99 g.h. */ 50701a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr); 50801a81e61SBarry Smith if (flg) { 50901a81e61SBarry Smith ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr); 51001a81e61SBarry Smith } 51101a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr); 51201a81e61SBarry Smith if (flg) { 51301a81e61SBarry Smith ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr); 51401a81e61SBarry Smith } 51501a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr); 51601a81e61SBarry Smith if (flg) { 51701a81e61SBarry Smith ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr); 51801a81e61SBarry Smith } 51901a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr); 52001a81e61SBarry Smith if (flg) { 52101a81e61SBarry Smith ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr); 52201a81e61SBarry Smith } 52301a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr); 52401a81e61SBarry Smith if (flg) { 52501a81e61SBarry Smith ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr); 52601a81e61SBarry Smith } 52701a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr); 52801a81e61SBarry Smith if (flg) { 52901a81e61SBarry Smith ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr); 53001a81e61SBarry Smith } 53101a81e61SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 53201a81e61SBarry Smith PetscFunctionReturn(0); 53301a81e61SBarry Smith } 53401a81e61SBarry Smith 53501a81e61SBarry Smith /**********************************************************************/ 53601a81e61SBarry Smith 53701a81e61SBarry Smith /*MC 53801a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 53901a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 54001a81e61SBarry Smith 54101a81e61SBarry Smith Options Database Keys: 542c5ae99e3SSatish Balay + -pc_spai_epsilon <eps> - set tolerance 54301a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 54401a81e61SBarry Smith . -pc_spai_max <m> - set max 54501a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 54601a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 54701a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 54801a81e61SBarry Smith . -pc_spai_sp <m> - set sp 54901a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 55001a81e61SBarry Smith 551*95452b02SPatrick Sanan Notes: 552*95452b02SPatrick Sanan This only works with AIJ matrices. 55301a81e61SBarry Smith 55401a81e61SBarry Smith Level: beginner 55501a81e61SBarry Smith 55601a81e61SBarry Smith Concepts: approximate inverse 55701a81e61SBarry Smith 55801a81e61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 55901a81e61SBarry Smith PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(), 56001a81e61SBarry Smith PCSPAISetVerbose(), PCSPAISetSp() 56101a81e61SBarry Smith M*/ 56201a81e61SBarry Smith 5638cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc) 56401a81e61SBarry Smith { 56501a81e61SBarry Smith PC_SPAI *ispai; 56601a81e61SBarry Smith PetscErrorCode ierr; 56701a81e61SBarry Smith 56801a81e61SBarry Smith PetscFunctionBegin; 569b00a9115SJed Brown ierr = PetscNewLog(pc,&ispai);CHKERRQ(ierr); 5702a4967abSBarry Smith pc->data = ispai; 57101a81e61SBarry Smith 57201a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 57301a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 57401a81e61SBarry Smith pc->ops->applyrichardson = 0; 57501a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 57601a81e61SBarry Smith pc->ops->view = PCView_SPAI; 57701a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 57801a81e61SBarry Smith 57901a81e61SBarry Smith ispai->epsilon = .4; 58001a81e61SBarry Smith ispai->nbsteps = 5; 58101a81e61SBarry Smith ispai->max = 5000; 58201a81e61SBarry Smith ispai->maxnew = 5; 58301a81e61SBarry Smith ispai->block_size = 1; 58401a81e61SBarry Smith ispai->cache_size = 5; 58501a81e61SBarry Smith ispai->verbose = 0; 58601a81e61SBarry Smith 58701a81e61SBarry Smith ispai->sp = 1; 588ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai));CHKERRQ(ierr); 58901a81e61SBarry Smith 590bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",PCSPAISetEpsilon_SPAI);CHKERRQ(ierr); 591bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",PCSPAISetNBSteps_SPAI);CHKERRQ(ierr); 592bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",PCSPAISetMax_SPAI);CHKERRQ(ierr); 593bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",PCSPAISetMaxNew_SPAI);CHKERRQ(ierr); 594bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",PCSPAISetBlockSize_SPAI);CHKERRQ(ierr); 595bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",PCSPAISetCacheSize_SPAI);CHKERRQ(ierr); 596bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",PCSPAISetVerbose_SPAI);CHKERRQ(ierr); 597bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",PCSPAISetSp_SPAI);CHKERRQ(ierr); 59801a81e61SBarry Smith PetscFunctionReturn(0); 59901a81e61SBarry Smith } 60001a81e61SBarry Smith 60101a81e61SBarry Smith /**********************************************************************/ 60201a81e61SBarry Smith 60301a81e61SBarry Smith /* 60401a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 60501a81e61SBarry Smith */ 60601a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B) 60701a81e61SBarry Smith { 60801a81e61SBarry Smith matrix *M; 60901a81e61SBarry Smith int i,j,col; 61001a81e61SBarry Smith int row_indx; 61101a81e61SBarry Smith int len,pe,local_indx,start_indx; 61201a81e61SBarry Smith int *mapping; 61301a81e61SBarry Smith PetscErrorCode ierr; 61401a81e61SBarry Smith const int *cols; 61501a81e61SBarry Smith const double *vals; 6162122902bSSatish Balay int n,mnl,nnl,nz,rstart,rend; 6172a4967abSBarry Smith PetscMPIInt size,rank; 61801a81e61SBarry Smith struct compressed_lines *rows; 61901a81e61SBarry Smith 62001a81e61SBarry Smith PetscFunctionBegin; 62101a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 62201a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 62301a81e61SBarry Smith ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr); 62401a81e61SBarry Smith ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr); 62501a81e61SBarry Smith 62601a81e61SBarry Smith /* 62701a81e61SBarry Smith not sure why a barrier is required. commenting out 62801a81e61SBarry Smith ierr = MPI_Barrier(comm);CHKERRQ(ierr); 62901a81e61SBarry Smith */ 63001a81e61SBarry Smith 63129b2f704SSatish Balay M = new_matrix((SPAI_Comm)comm); 63201a81e61SBarry Smith 63301a81e61SBarry Smith M->n = n; 63401a81e61SBarry Smith M->bs = 1; 63501a81e61SBarry Smith M->max_block_size = 1; 63601a81e61SBarry Smith 63701a81e61SBarry Smith M->mnls = (int*)malloc(sizeof(int)*size); 63801a81e61SBarry Smith M->start_indices = (int*)malloc(sizeof(int)*size); 63901a81e61SBarry Smith M->pe = (int*)malloc(sizeof(int)*n); 64001a81e61SBarry Smith M->block_sizes = (int*)malloc(sizeof(int)*n); 64101a81e61SBarry Smith for (i=0; i<n; i++) M->block_sizes[i] = 1; 64201a81e61SBarry Smith 64301a81e61SBarry Smith ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr); 64401a81e61SBarry Smith 64501a81e61SBarry Smith M->start_indices[0] = 0; 6462fa5cd67SKarl Rupp for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1]; 64701a81e61SBarry Smith 64801a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 64901a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 65001a81e61SBarry Smith 65101a81e61SBarry Smith for (i=0; i<size; i++) { 65201a81e61SBarry Smith start_indx = M->start_indices[i]; 6532fa5cd67SKarl Rupp for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i; 65401a81e61SBarry Smith } 65501a81e61SBarry Smith 65601a81e61SBarry Smith if (AT) { 65701a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr); 65801a81e61SBarry Smith } else { 65901a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr); 66001a81e61SBarry Smith } 66101a81e61SBarry Smith 66201a81e61SBarry Smith rows = M->lines; 66301a81e61SBarry Smith 66401a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 665785e854fSJed Brown ierr = PetscMalloc1(M->n,&mapping);CHKERRQ(ierr); 66601a81e61SBarry Smith pe = 0; 66701a81e61SBarry Smith local_indx = 0; 66801a81e61SBarry Smith for (i=0; i<M->n; i++) { 66901a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 67001a81e61SBarry Smith pe++; 67101a81e61SBarry Smith local_indx = 0; 67201a81e61SBarry Smith } 67301a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 67401a81e61SBarry Smith local_indx++; 67501a81e61SBarry Smith } 67601a81e61SBarry Smith 67701a81e61SBarry Smith /*********************************************************/ 67801a81e61SBarry Smith /************** Set up the row structure *****************/ 67901a81e61SBarry Smith /*********************************************************/ 68001a81e61SBarry Smith 68101a81e61SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 68201a81e61SBarry Smith for (i=rstart; i<rend; i++) { 68301a81e61SBarry Smith row_indx = i - rstart; 68401a81e61SBarry Smith ierr = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 6852122902bSSatish Balay /* allocate buffers */ 6862122902bSSatish Balay rows->ptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 6872122902bSSatish Balay rows->A[row_indx] = (double*)malloc(nz*sizeof(double)); 6882122902bSSatish Balay /* copy the matrix */ 68901a81e61SBarry Smith for (j=0; j<nz; j++) { 69001a81e61SBarry Smith col = cols[j]; 69101a81e61SBarry Smith len = rows->len[row_indx]++; 6922fa5cd67SKarl Rupp 69301a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 69401a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 69501a81e61SBarry Smith } 69601a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 6972fa5cd67SKarl Rupp 69801a81e61SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 69901a81e61SBarry Smith } 70001a81e61SBarry Smith 70101a81e61SBarry Smith 70201a81e61SBarry Smith /************************************************************/ 70301a81e61SBarry Smith /************** Set up the column structure *****************/ 70401a81e61SBarry Smith /*********************************************************/ 70501a81e61SBarry Smith 70601a81e61SBarry Smith if (AT) { 70701a81e61SBarry Smith 70801a81e61SBarry Smith for (i=rstart; i<rend; i++) { 70901a81e61SBarry Smith row_indx = i - rstart; 71001a81e61SBarry Smith ierr = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 7112122902bSSatish Balay /* allocate buffers */ 7122122902bSSatish Balay rows->rptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 7132122902bSSatish Balay /* copy the matrix (i.e., the structure) */ 71401a81e61SBarry Smith for (j=0; j<nz; j++) { 71501a81e61SBarry Smith col = cols[j]; 71601a81e61SBarry Smith len = rows->rlen[row_indx]++; 7172fa5cd67SKarl Rupp 71801a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 71901a81e61SBarry Smith } 72001a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 72101a81e61SBarry Smith } 72201a81e61SBarry Smith } 72301a81e61SBarry Smith 72401a81e61SBarry Smith ierr = PetscFree(mapping);CHKERRQ(ierr); 72501a81e61SBarry Smith 72601a81e61SBarry Smith order_pointers(M); 72701a81e61SBarry Smith M->maxnz = calc_maxnz(M); 72801a81e61SBarry Smith *B = M; 72901a81e61SBarry Smith PetscFunctionReturn(0); 73001a81e61SBarry Smith } 73101a81e61SBarry Smith 73201a81e61SBarry Smith /**********************************************************************/ 73301a81e61SBarry Smith 73401a81e61SBarry Smith /* 73501a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 736df3898eeSBarry Smith This assumes that the SPAI matrix B is stored in 73701a81e61SBarry Smith COMPRESSED-ROW format. 73801a81e61SBarry Smith */ 73901a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB) 74001a81e61SBarry Smith { 7414b6b5fe1SBarry Smith PetscMPIInt size,rank; 74201a81e61SBarry Smith PetscErrorCode ierr; 74301a81e61SBarry Smith int m,n,M,N; 74401a81e61SBarry Smith int d_nz,o_nz; 74501a81e61SBarry Smith int *d_nnz,*o_nnz; 74601a81e61SBarry Smith int i,k,global_row,global_col,first_diag_col,last_diag_col; 74701a81e61SBarry Smith PetscScalar val; 74801a81e61SBarry Smith 74901a81e61SBarry Smith PetscFunctionBegin; 75001a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 75101a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 75201a81e61SBarry Smith 75301a81e61SBarry Smith m = n = B->mnls[rank]; 75401a81e61SBarry Smith d_nz = o_nz = 0; 75501a81e61SBarry Smith 75601a81e61SBarry Smith /* Determine preallocation for MatCreateMPIAIJ */ 757785e854fSJed Brown ierr = PetscMalloc1(m,&d_nnz);CHKERRQ(ierr); 758785e854fSJed Brown ierr = PetscMalloc1(m,&o_nnz);CHKERRQ(ierr); 75901a81e61SBarry Smith for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0; 76001a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 76101a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 76201a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 76301a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 76401a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 765db4deed7SKarl Rupp if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++; 766db4deed7SKarl Rupp else o_nnz[i]++; 76701a81e61SBarry Smith } 76801a81e61SBarry Smith } 76901a81e61SBarry Smith 77001a81e61SBarry Smith M = N = B->n; 77101a81e61SBarry Smith /* Here we only know how to create AIJ format */ 772f69a0ea3SMatthew Knepley ierr = MatCreate(comm,PB);CHKERRQ(ierr); 773f69a0ea3SMatthew Knepley ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr); 77401a81e61SBarry Smith ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr); 77501a81e61SBarry Smith ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr); 77601a81e61SBarry Smith ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr); 77701a81e61SBarry Smith 77801a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 77901a81e61SBarry Smith global_row = B->start_indices[rank]+i; 78001a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 78101a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 7822fa5cd67SKarl Rupp 78301a81e61SBarry Smith val = B->lines->A[i][k]; 78401a81e61SBarry Smith ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr); 78501a81e61SBarry Smith } 78601a81e61SBarry Smith } 78701a81e61SBarry Smith 78801a81e61SBarry Smith ierr = PetscFree(d_nnz);CHKERRQ(ierr); 78901a81e61SBarry Smith ierr = PetscFree(o_nnz);CHKERRQ(ierr); 79001a81e61SBarry Smith 79101a81e61SBarry Smith ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79201a81e61SBarry Smith ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79301a81e61SBarry Smith PetscFunctionReturn(0); 79401a81e61SBarry Smith } 79501a81e61SBarry Smith 79601a81e61SBarry Smith /**********************************************************************/ 79701a81e61SBarry Smith 79801a81e61SBarry Smith /* 79901a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 80001a81e61SBarry Smith */ 80101a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv) 80201a81e61SBarry Smith { 80301a81e61SBarry Smith PetscErrorCode ierr; 8044b6b5fe1SBarry Smith PetscMPIInt size,rank; 8054b6b5fe1SBarry Smith int m,M,i,*mnls,*start_indices,*global_indices; 80601a81e61SBarry Smith 80701a81e61SBarry Smith PetscFunctionBegin; 80801a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 80901a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 81001a81e61SBarry Smith 81101a81e61SBarry Smith m = v->mnl; 81201a81e61SBarry Smith M = v->n; 81301a81e61SBarry Smith 81401a81e61SBarry Smith 81501a81e61SBarry Smith ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr); 81601a81e61SBarry Smith 817785e854fSJed Brown ierr = PetscMalloc1(size,&mnls);CHKERRQ(ierr); 8182a4967abSBarry Smith ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRQ(ierr); 81901a81e61SBarry Smith 820785e854fSJed Brown ierr = PetscMalloc1(size,&start_indices);CHKERRQ(ierr); 8212fa5cd67SKarl Rupp 82201a81e61SBarry Smith start_indices[0] = 0; 8232fa5cd67SKarl Rupp for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1]; 82401a81e61SBarry Smith 825785e854fSJed Brown ierr = PetscMalloc1(v->mnl,&global_indices);CHKERRQ(ierr); 8262fa5cd67SKarl Rupp for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i; 82701a81e61SBarry Smith 82801a81e61SBarry Smith ierr = PetscFree(mnls);CHKERRQ(ierr); 82901a81e61SBarry Smith ierr = PetscFree(start_indices);CHKERRQ(ierr); 83001a81e61SBarry Smith 83101a81e61SBarry Smith ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr); 83201a81e61SBarry Smith ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr); 83301a81e61SBarry Smith ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr); 83401a81e61SBarry Smith 8354b6b5fe1SBarry Smith ierr = PetscFree(global_indices);CHKERRQ(ierr); 83601a81e61SBarry Smith PetscFunctionReturn(0); 83701a81e61SBarry Smith } 83801a81e61SBarry Smith 83901a81e61SBarry Smith 84001a81e61SBarry Smith 84101a81e61SBarry Smith 842