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 *); 36b03515a0SUmberto Zerbinati extern PetscErrorCode ConvertVectorToVec(MPI_Comm, vector *, Vec *); 3709573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char *, char *, char *); 3801a81e61SBarry Smith 3901a81e61SBarry Smith typedef struct { 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 60*9371c9d4SSatish Balay static PetscErrorCode PCSetUp_SPAI(PC pc) { 6101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 6201a81e61SBarry Smith Mat AT; 6301a81e61SBarry Smith 6401a81e61SBarry Smith PetscFunctionBegin; 6501a81e61SBarry Smith init_SPAI(); 6601a81e61SBarry Smith 6701a81e61SBarry Smith if (ispai->sp) { 689566063dSJacob Faibussowitsch PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, pc->pmat, &ispai->B)); 6901a81e61SBarry Smith } else { 7001a81e61SBarry Smith /* Use the transpose to get the column nonzero structure. */ 719566063dSJacob Faibussowitsch PetscCall(MatTranspose(pc->pmat, MAT_INITIAL_MATRIX, &AT)); 729566063dSJacob Faibussowitsch PetscCall(ConvertMatToMatrix(ispai->comm_spai, pc->pmat, AT, &ispai->B)); 739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&AT)); 7401a81e61SBarry Smith } 7501a81e61SBarry Smith 7601a81e61SBarry Smith /* Destroy the transpose */ 7701a81e61SBarry Smith /* Don't know how to do it. PETSc developers? */ 7801a81e61SBarry Smith 7901a81e61SBarry Smith /* construct SPAI preconditioner */ 8001a81e61SBarry Smith /* FILE *messages */ /* file for warning messages */ 8101a81e61SBarry Smith /* double epsilon */ /* tolerance */ 8201a81e61SBarry Smith /* int nbsteps */ /* max number of "improvement" steps per line */ 8301a81e61SBarry Smith /* int max */ /* max dimensions of is_I, q, etc. */ 8401a81e61SBarry Smith /* int maxnew */ /* max number of new entries per step */ 8501a81e61SBarry Smith /* int block_size */ /* block_size == 1 specifies scalar elments 8601a81e61SBarry Smith block_size == n specifies nxn constant-block elements 8701a81e61SBarry Smith block_size == 0 specifies variable-block elements */ 882fa5cd67SKarl Rupp /* int cache_size */ /* one of (1,2,3,4,5,6) indicting size of cache. cache_size == 0 indicates no caching */ 8901a81e61SBarry Smith /* int verbose */ /* verbose == 0 specifies that SPAI is silent 9001a81e61SBarry Smith verbose == 1 prints timing and matrix statistics */ 9101a81e61SBarry Smith 92*9371c9d4SSatish Balay PetscCall(bspai(ispai->B, &ispai->M, stdout, ispai->epsilon, ispai->nbsteps, ispai->max, ispai->maxnew, ispai->block_size, ispai->cache_size, ispai->verbose)); 9301a81e61SBarry Smith 949566063dSJacob Faibussowitsch PetscCall(ConvertMatrixToMat(PetscObjectComm((PetscObject)pc), ispai->M, &ispai->PM)); 9501a81e61SBarry Smith 9601a81e61SBarry Smith /* free the SPAI matrices */ 9701a81e61SBarry Smith sp_free_matrix(ispai->B); 9801a81e61SBarry Smith sp_free_matrix(ispai->M); 9901a81e61SBarry Smith PetscFunctionReturn(0); 10001a81e61SBarry Smith } 10101a81e61SBarry Smith 10201a81e61SBarry Smith /**********************************************************************/ 10301a81e61SBarry Smith 104*9371c9d4SSatish Balay static PetscErrorCode PCApply_SPAI(PC pc, Vec xx, Vec y) { 10501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 10601a81e61SBarry Smith 10701a81e61SBarry Smith PetscFunctionBegin; 10801a81e61SBarry Smith /* Now using PETSc's multiply */ 1099566063dSJacob Faibussowitsch PetscCall(MatMult(ispai->PM, xx, y)); 11001a81e61SBarry Smith PetscFunctionReturn(0); 11101a81e61SBarry Smith } 11201a81e61SBarry Smith 113*9371c9d4SSatish Balay static PetscErrorCode PCMatApply_SPAI(PC pc, Mat X, Mat Y) { 11448e38a8aSPierre Jolivet PC_SPAI *ispai = (PC_SPAI *)pc->data; 11548e38a8aSPierre Jolivet 11648e38a8aSPierre Jolivet PetscFunctionBegin; 11748e38a8aSPierre Jolivet /* Now using PETSc's multiply */ 1189566063dSJacob Faibussowitsch PetscCall(MatMatMult(ispai->PM, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 11948e38a8aSPierre Jolivet PetscFunctionReturn(0); 12048e38a8aSPierre Jolivet } 12148e38a8aSPierre Jolivet 12201a81e61SBarry Smith /**********************************************************************/ 12301a81e61SBarry Smith 124*9371c9d4SSatish Balay static PetscErrorCode PCDestroy_SPAI(PC pc) { 12501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 12601a81e61SBarry Smith 12701a81e61SBarry Smith PetscFunctionBegin; 1289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ispai->PM)); 1299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&(ispai->comm_spai))); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", NULL)); 1322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", NULL)); 1332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", NULL)); 1342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", NULL)); 1352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", NULL)); 1362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", NULL)); 1372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", NULL)); 1382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", NULL)); 13901a81e61SBarry Smith PetscFunctionReturn(0); 14001a81e61SBarry Smith } 14101a81e61SBarry Smith 14201a81e61SBarry Smith /**********************************************************************/ 14301a81e61SBarry Smith 144*9371c9d4SSatish Balay static PetscErrorCode PCView_SPAI(PC pc, PetscViewer viewer) { 14501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 146ace3abfcSBarry Smith PetscBool iascii; 14701a81e61SBarry Smith 14801a81e61SBarry Smith PetscFunctionBegin; 1499566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 15001a81e61SBarry Smith if (iascii) { 151b03515a0SUmberto Zerbinati PetscCall(PetscViewerASCIIPrintf(viewer, " epsilon %g\n", ispai->epsilon)); 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " nbsteps %d\n", ispai->nbsteps)); 1539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " max %d\n", ispai->max)); 1549566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " maxnew %d\n", ispai->maxnew)); 1559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " block_size %d\n", ispai->block_size)); 1569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " cache_size %d\n", ispai->cache_size)); 1579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " verbose %d\n", ispai->verbose)); 1589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " sp %d\n", ispai->sp)); 15901a81e61SBarry Smith } 16001a81e61SBarry Smith PetscFunctionReturn(0); 16101a81e61SBarry Smith } 16201a81e61SBarry Smith 163*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc, PetscReal epsilon1) { 16401a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 1655fd66863SKarl Rupp 16601a81e61SBarry Smith PetscFunctionBegin; 167b03515a0SUmberto Zerbinati ispai->epsilon = (double)epsilon1; 16801a81e61SBarry Smith PetscFunctionReturn(0); 16901a81e61SBarry Smith } 17001a81e61SBarry Smith 17101a81e61SBarry Smith /**********************************************************************/ 17201a81e61SBarry Smith 173*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc, PetscInt nbsteps1) { 17401a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 1755fd66863SKarl Rupp 17601a81e61SBarry Smith PetscFunctionBegin; 177b03515a0SUmberto Zerbinati ispai->nbsteps = (int)nbsteps1; 17801a81e61SBarry Smith PetscFunctionReturn(0); 17901a81e61SBarry Smith } 18001a81e61SBarry Smith 18101a81e61SBarry Smith /**********************************************************************/ 18201a81e61SBarry Smith 18301a81e61SBarry Smith /* added 1/7/99 g.h. */ 184*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetMax_SPAI(PC pc, PetscInt max1) { 18501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 1865fd66863SKarl Rupp 18701a81e61SBarry Smith PetscFunctionBegin; 188b03515a0SUmberto Zerbinati ispai->max = (int)max1; 18901a81e61SBarry Smith PetscFunctionReturn(0); 19001a81e61SBarry Smith } 19101a81e61SBarry Smith 19201a81e61SBarry Smith /**********************************************************************/ 19301a81e61SBarry Smith 194*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc, PetscInt maxnew1) { 19501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 1965fd66863SKarl Rupp 19701a81e61SBarry Smith PetscFunctionBegin; 198b03515a0SUmberto Zerbinati ispai->maxnew = (int)maxnew1; 19901a81e61SBarry Smith PetscFunctionReturn(0); 20001a81e61SBarry Smith } 20101a81e61SBarry Smith 20201a81e61SBarry Smith /**********************************************************************/ 20301a81e61SBarry Smith 204*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc, PetscInt block_size1) { 20501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 2065fd66863SKarl Rupp 20701a81e61SBarry Smith PetscFunctionBegin; 208b03515a0SUmberto Zerbinati ispai->block_size = (int)block_size1; 20901a81e61SBarry Smith PetscFunctionReturn(0); 21001a81e61SBarry Smith } 21101a81e61SBarry Smith 21201a81e61SBarry Smith /**********************************************************************/ 21301a81e61SBarry Smith 214*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc, PetscInt cache_size) { 21501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 2165fd66863SKarl Rupp 21701a81e61SBarry Smith PetscFunctionBegin; 218b03515a0SUmberto Zerbinati ispai->cache_size = (int)cache_size; 21901a81e61SBarry Smith PetscFunctionReturn(0); 22001a81e61SBarry Smith } 22101a81e61SBarry Smith 22201a81e61SBarry Smith /**********************************************************************/ 22301a81e61SBarry Smith 224*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc, PetscInt verbose) { 22501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 2265fd66863SKarl Rupp 22701a81e61SBarry Smith PetscFunctionBegin; 228b03515a0SUmberto Zerbinati ispai->verbose = (int)verbose; 22901a81e61SBarry Smith PetscFunctionReturn(0); 23001a81e61SBarry Smith } 23101a81e61SBarry Smith 23201a81e61SBarry Smith /**********************************************************************/ 23301a81e61SBarry Smith 234*9371c9d4SSatish Balay static PetscErrorCode PCSPAISetSp_SPAI(PC pc, PetscInt sp) { 23501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 2365fd66863SKarl Rupp 23701a81e61SBarry Smith PetscFunctionBegin; 238b03515a0SUmberto Zerbinati ispai->sp = (int)sp; 23901a81e61SBarry Smith PetscFunctionReturn(0); 24001a81e61SBarry Smith } 24101a81e61SBarry Smith 24201a81e61SBarry Smith /* -------------------------------------------------------------------*/ 24301a81e61SBarry Smith 24401a81e61SBarry Smith /*@ 24501a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 24601a81e61SBarry Smith 24701a81e61SBarry Smith Input Parameters: 24801a81e61SBarry Smith + pc - the preconditioner 24901a81e61SBarry Smith - eps - epsilon (default .4) 25001a81e61SBarry Smith 25195452b02SPatrick Sanan Notes: 25295452b02SPatrick Sanan Espilon must be between 0 and 1. It controls the 25301a81e61SBarry Smith quality of the approximation of M to the inverse of 25401a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 25501a81e61SBarry Smith fill, and usually better preconditioners. In many 25601a81e61SBarry Smith cases the best choice of epsilon is the one that 25701a81e61SBarry Smith divides the total solution time equally between the 25801a81e61SBarry Smith preconditioner and the solver. 25901a81e61SBarry Smith 26001a81e61SBarry Smith Level: intermediate 26101a81e61SBarry Smith 262db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 26301a81e61SBarry Smith @*/ 264*9371c9d4SSatish Balay PetscErrorCode PCSPAISetEpsilon(PC pc, PetscReal epsilon1) { 26501a81e61SBarry Smith PetscFunctionBegin; 266b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetEpsilon_C", (PC, PetscReal), (pc, epsilon1)); 26701a81e61SBarry Smith PetscFunctionReturn(0); 26801a81e61SBarry Smith } 26901a81e61SBarry Smith 27001a81e61SBarry Smith /**********************************************************************/ 27101a81e61SBarry Smith 27201a81e61SBarry Smith /*@ 27301a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 27401a81e61SBarry Smith the SPAI preconditioner 27501a81e61SBarry Smith 27601a81e61SBarry Smith Input Parameters: 27701a81e61SBarry Smith + pc - the preconditioner 27801a81e61SBarry Smith - n - number of steps (default 5) 27901a81e61SBarry Smith 28095452b02SPatrick Sanan Notes: 28195452b02SPatrick Sanan SPAI constructs to approximation to every column of 28201a81e61SBarry Smith the exact inverse of A in a series of improvement 28301a81e61SBarry Smith steps. The quality of the approximation is determined 28401a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 28501a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 28601a81e61SBarry Smith uses the best approximation constructed so far. 28701a81e61SBarry Smith 28801a81e61SBarry Smith Level: intermediate 28901a81e61SBarry Smith 290db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetMaxNew()` 29101a81e61SBarry Smith @*/ 292*9371c9d4SSatish Balay PetscErrorCode PCSPAISetNBSteps(PC pc, PetscInt nbsteps1) { 29301a81e61SBarry Smith PetscFunctionBegin; 294b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetNBSteps_C", (PC, PetscInt), (pc, nbsteps1)); 29501a81e61SBarry Smith PetscFunctionReturn(0); 29601a81e61SBarry Smith } 29701a81e61SBarry Smith 29801a81e61SBarry Smith /**********************************************************************/ 29901a81e61SBarry Smith 30001a81e61SBarry Smith /* added 1/7/99 g.h. */ 30101a81e61SBarry Smith /*@ 30201a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 30301a81e61SBarry Smith the SPAI preconditioner 30401a81e61SBarry Smith 30501a81e61SBarry Smith Input Parameters: 30601a81e61SBarry Smith + pc - the preconditioner 30701a81e61SBarry Smith - n - size (default is 5000) 30801a81e61SBarry Smith 30901a81e61SBarry Smith Level: intermediate 31001a81e61SBarry Smith 311db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 31201a81e61SBarry Smith @*/ 313*9371c9d4SSatish Balay PetscErrorCode PCSPAISetMax(PC pc, PetscInt max1) { 31401a81e61SBarry Smith PetscFunctionBegin; 315b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetMax_C", (PC, PetscInt), (pc, max1)); 31601a81e61SBarry Smith PetscFunctionReturn(0); 31701a81e61SBarry Smith } 31801a81e61SBarry Smith 31901a81e61SBarry Smith /**********************************************************************/ 32001a81e61SBarry Smith 32101a81e61SBarry Smith /*@ 32201a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 32301a81e61SBarry Smith in SPAI preconditioner 32401a81e61SBarry Smith 32501a81e61SBarry Smith Input Parameters: 32601a81e61SBarry Smith + pc - the preconditioner 32701a81e61SBarry Smith - n - maximum number (default 5) 32801a81e61SBarry Smith 32901a81e61SBarry Smith Level: intermediate 33001a81e61SBarry Smith 331db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()`, `PCSPAISetNBSteps()` 33201a81e61SBarry Smith @*/ 333*9371c9d4SSatish Balay PetscErrorCode PCSPAISetMaxNew(PC pc, PetscInt maxnew1) { 33401a81e61SBarry Smith PetscFunctionBegin; 335b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetMaxNew_C", (PC, PetscInt), (pc, maxnew1)); 33601a81e61SBarry Smith PetscFunctionReturn(0); 33701a81e61SBarry Smith } 33801a81e61SBarry Smith 33901a81e61SBarry Smith /**********************************************************************/ 34001a81e61SBarry Smith 34101a81e61SBarry Smith /*@ 34201a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 34301a81e61SBarry Smith 34401a81e61SBarry Smith Input Parameters: 34501a81e61SBarry Smith + pc - the preconditioner 34601a81e61SBarry Smith - n - block size (default 1) 34701a81e61SBarry Smith 34895452b02SPatrick Sanan Notes: 34995452b02SPatrick Sanan A block 35001a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 35101a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 35201a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 35301a81e61SBarry Smith variable sized blocks, which are determined by 35401a81e61SBarry Smith searching for dense square diagonal blocks in A. 35501a81e61SBarry Smith This can be very effective for finite-element 35601a81e61SBarry Smith matrices. 35701a81e61SBarry Smith 35801a81e61SBarry Smith SPAI will convert A to block form, use a block 35901a81e61SBarry Smith version of the preconditioner algorithm, and then 36001a81e61SBarry Smith convert the result back to scalar form. 36101a81e61SBarry Smith 36201a81e61SBarry Smith In many cases the a block-size parameter other than 1 36301a81e61SBarry Smith can lead to very significant improvement in 36401a81e61SBarry Smith performance. 36501a81e61SBarry Smith 36601a81e61SBarry Smith Level: intermediate 36701a81e61SBarry Smith 368db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 36901a81e61SBarry Smith @*/ 370*9371c9d4SSatish Balay PetscErrorCode PCSPAISetBlockSize(PC pc, PetscInt block_size1) { 37101a81e61SBarry Smith PetscFunctionBegin; 372b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetBlockSize_C", (PC, PetscInt), (pc, block_size1)); 37301a81e61SBarry Smith PetscFunctionReturn(0); 37401a81e61SBarry Smith } 37501a81e61SBarry Smith 37601a81e61SBarry Smith /**********************************************************************/ 37701a81e61SBarry Smith 37801a81e61SBarry Smith /*@ 37901a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 38001a81e61SBarry Smith 38101a81e61SBarry Smith Input Parameters: 38201a81e61SBarry Smith + pc - the preconditioner 38301a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 38401a81e61SBarry Smith 38595452b02SPatrick Sanan Notes: 38695452b02SPatrick Sanan SPAI uses a hash table to cache messages and avoid 38701a81e61SBarry Smith redundant communication. If suggest always using 38801a81e61SBarry Smith 5. This parameter is irrelevant in the serial 38901a81e61SBarry Smith version. 39001a81e61SBarry Smith 39101a81e61SBarry Smith Level: intermediate 39201a81e61SBarry Smith 393db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 39401a81e61SBarry Smith @*/ 395*9371c9d4SSatish Balay PetscErrorCode PCSPAISetCacheSize(PC pc, PetscInt cache_size) { 39601a81e61SBarry Smith PetscFunctionBegin; 397b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetCacheSize_C", (PC, PetscInt), (pc, cache_size)); 39801a81e61SBarry Smith PetscFunctionReturn(0); 39901a81e61SBarry Smith } 40001a81e61SBarry Smith 40101a81e61SBarry Smith /**********************************************************************/ 40201a81e61SBarry Smith 40301a81e61SBarry Smith /*@ 40401a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 40501a81e61SBarry Smith 40601a81e61SBarry Smith Input Parameters: 40701a81e61SBarry Smith + pc - the preconditioner 40801a81e61SBarry Smith - n - level (default 1) 40901a81e61SBarry Smith 41095452b02SPatrick Sanan Notes: 41195452b02SPatrick Sanan print parameters, timings and matrix statistics 41201a81e61SBarry Smith 41301a81e61SBarry Smith Level: intermediate 41401a81e61SBarry Smith 415db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 41601a81e61SBarry Smith @*/ 417*9371c9d4SSatish Balay PetscErrorCode PCSPAISetVerbose(PC pc, PetscInt verbose) { 41801a81e61SBarry Smith PetscFunctionBegin; 419b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetVerbose_C", (PC, PetscInt), (pc, verbose)); 42001a81e61SBarry Smith PetscFunctionReturn(0); 42101a81e61SBarry Smith } 42201a81e61SBarry Smith 42301a81e61SBarry Smith /**********************************************************************/ 42401a81e61SBarry Smith 42501a81e61SBarry Smith /*@ 42601a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 42701a81e61SBarry Smith 42801a81e61SBarry Smith Input Parameters: 42901a81e61SBarry Smith + pc - the preconditioner 43001a81e61SBarry Smith - n - 0 or 1 43101a81e61SBarry Smith 43295452b02SPatrick Sanan Notes: 43395452b02SPatrick Sanan If A has a symmetric nonzero pattern use -sp 1 to 43401a81e61SBarry Smith improve performance by eliminating some communication 43501a81e61SBarry Smith in the parallel version. Even if A does not have a 43601a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 43701a81e61SBarry Smith results, but the code will not follow the published 43801a81e61SBarry Smith SPAI algorithm exactly. 43901a81e61SBarry Smith 44001a81e61SBarry Smith Level: intermediate 44101a81e61SBarry Smith 442db781477SPatrick Sanan .seealso: `PCSPAI`, `PCSetType()` 44301a81e61SBarry Smith @*/ 444*9371c9d4SSatish Balay PetscErrorCode PCSPAISetSp(PC pc, PetscInt sp) { 44501a81e61SBarry Smith PetscFunctionBegin; 446b03515a0SUmberto Zerbinati PetscTryMethod(pc, "PCSPAISetSp_C", (PC, PetscInt), (pc, sp)); 44701a81e61SBarry Smith PetscFunctionReturn(0); 44801a81e61SBarry Smith } 44901a81e61SBarry Smith 45001a81e61SBarry Smith /**********************************************************************/ 45101a81e61SBarry Smith 45201a81e61SBarry Smith /**********************************************************************/ 45301a81e61SBarry Smith 454*9371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_SPAI(PC pc, PetscOptionItems *PetscOptionsObject) { 45501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI *)pc->data; 45601a81e61SBarry Smith int nbsteps1, max1, maxnew1, block_size1, cache_size, verbose, sp; 45701a81e61SBarry Smith double epsilon1; 458ace3abfcSBarry Smith PetscBool flg; 45901a81e61SBarry Smith 46001a81e61SBarry Smith PetscFunctionBegin; 461d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SPAI options"); 4629566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_spai_epsilon", "", "PCSPAISetEpsilon", ispai->epsilon, &epsilon1, &flg)); 4631baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetEpsilon(pc, epsilon1)); 4649566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_nbsteps", "", "PCSPAISetNBSteps", ispai->nbsteps, &nbsteps1, &flg)); 4651baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetNBSteps(pc, nbsteps1)); 46601a81e61SBarry Smith /* added 1/7/99 g.h. */ 4679566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_max", "", "PCSPAISetMax", ispai->max, &max1, &flg)); 4681baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetMax(pc, max1)); 4699566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_maxnew", "", "PCSPAISetMaxNew", ispai->maxnew, &maxnew1, &flg)); 4701baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetMaxNew(pc, maxnew1)); 4719566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_block_size", "", "PCSPAISetBlockSize", ispai->block_size, &block_size1, &flg)); 4721baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetBlockSize(pc, block_size1)); 4739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_cache_size", "", "PCSPAISetCacheSize", ispai->cache_size, &cache_size, &flg)); 4741baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetCacheSize(pc, cache_size)); 4759566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_verbose", "", "PCSPAISetVerbose", ispai->verbose, &verbose, &flg)); 4761baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetVerbose(pc, verbose)); 4779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_spai_sp", "", "PCSPAISetSp", ispai->sp, &sp, &flg)); 4781baa6e33SBarry Smith if (flg) PetscCall(PCSPAISetSp(pc, sp)); 479d0609cedSBarry Smith PetscOptionsHeadEnd(); 48001a81e61SBarry Smith PetscFunctionReturn(0); 48101a81e61SBarry Smith } 48201a81e61SBarry Smith 48301a81e61SBarry Smith /**********************************************************************/ 48401a81e61SBarry Smith 48501a81e61SBarry Smith /*MC 48601a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 48701a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 48801a81e61SBarry Smith 48901a81e61SBarry Smith Options Database Keys: 490c5ae99e3SSatish Balay + -pc_spai_epsilon <eps> - set tolerance 49101a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 49201a81e61SBarry Smith . -pc_spai_max <m> - set max 49301a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 49401a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 49501a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 49601a81e61SBarry Smith . -pc_spai_sp <m> - set sp 49701a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 49801a81e61SBarry Smith 49995452b02SPatrick Sanan Notes: 50095452b02SPatrick Sanan This only works with AIJ matrices. 50101a81e61SBarry Smith 50201a81e61SBarry Smith Level: beginner 50301a81e61SBarry Smith 504db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 505db781477SPatrick Sanan `PCSPAISetEpsilon()`, `PCSPAISetMax()`, `PCSPAISetMaxNew()`, `PCSPAISetBlockSize()`, 506db781477SPatrick Sanan `PCSPAISetVerbose()`, `PCSPAISetSp()` 50701a81e61SBarry Smith M*/ 50801a81e61SBarry Smith 509*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc) { 51001a81e61SBarry Smith PC_SPAI *ispai; 51101a81e61SBarry Smith 51201a81e61SBarry Smith PetscFunctionBegin; 5139566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &ispai)); 5142a4967abSBarry Smith pc->data = ispai; 51501a81e61SBarry Smith 51601a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 51701a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 51848e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_SPAI; 51901a81e61SBarry Smith pc->ops->applyrichardson = 0; 52001a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 52101a81e61SBarry Smith pc->ops->view = PCView_SPAI; 52201a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 52301a81e61SBarry Smith 52401a81e61SBarry Smith ispai->epsilon = .4; 52501a81e61SBarry Smith ispai->nbsteps = 5; 52601a81e61SBarry Smith ispai->max = 5000; 52701a81e61SBarry Smith ispai->maxnew = 5; 52801a81e61SBarry Smith ispai->block_size = 1; 52901a81e61SBarry Smith ispai->cache_size = 5; 53001a81e61SBarry Smith ispai->verbose = 0; 53101a81e61SBarry Smith 53201a81e61SBarry Smith ispai->sp = 1; 5339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_dup(PetscObjectComm((PetscObject)pc), &(ispai->comm_spai))); 53401a81e61SBarry Smith 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetEpsilon_C", PCSPAISetEpsilon_SPAI)); 5369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetNBSteps_C", PCSPAISetNBSteps_SPAI)); 5379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMax_C", PCSPAISetMax_SPAI)); 5389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetMaxNew_C", PCSPAISetMaxNew_SPAI)); 5399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetBlockSize_C", PCSPAISetBlockSize_SPAI)); 5409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetCacheSize_C", PCSPAISetCacheSize_SPAI)); 5419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetVerbose_C", PCSPAISetVerbose_SPAI)); 5429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSPAISetSp_C", PCSPAISetSp_SPAI)); 54301a81e61SBarry Smith PetscFunctionReturn(0); 54401a81e61SBarry Smith } 54501a81e61SBarry Smith 54601a81e61SBarry Smith /**********************************************************************/ 54701a81e61SBarry Smith 54801a81e61SBarry Smith /* 54901a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 55001a81e61SBarry Smith */ 551*9371c9d4SSatish Balay PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A, Mat AT, matrix **B) { 55201a81e61SBarry Smith matrix *M; 55301a81e61SBarry Smith int i, j, col; 55401a81e61SBarry Smith int row_indx; 55501a81e61SBarry Smith int len, pe, local_indx, start_indx; 55601a81e61SBarry Smith int *mapping; 55701a81e61SBarry Smith const int *cols; 55801a81e61SBarry Smith const double *vals; 5592122902bSSatish Balay int n, mnl, nnl, nz, rstart, rend; 5602a4967abSBarry Smith PetscMPIInt size, rank; 56101a81e61SBarry Smith struct compressed_lines *rows; 56201a81e61SBarry Smith 56301a81e61SBarry Smith PetscFunctionBegin; 5649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 5659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 5669566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &n, &n)); 5679566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &mnl, &nnl)); 56801a81e61SBarry Smith 56901a81e61SBarry Smith /* 57001a81e61SBarry Smith not sure why a barrier is required. commenting out 5719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(comm)); 57201a81e61SBarry Smith */ 57301a81e61SBarry Smith 57429b2f704SSatish Balay M = new_matrix((SPAI_Comm)comm); 57501a81e61SBarry Smith 57601a81e61SBarry Smith M->n = n; 57701a81e61SBarry Smith M->bs = 1; 57801a81e61SBarry Smith M->max_block_size = 1; 57901a81e61SBarry Smith 58001a81e61SBarry Smith M->mnls = (int *)malloc(sizeof(int) * size); 58101a81e61SBarry Smith M->start_indices = (int *)malloc(sizeof(int) * size); 58201a81e61SBarry Smith M->pe = (int *)malloc(sizeof(int) * n); 58301a81e61SBarry Smith M->block_sizes = (int *)malloc(sizeof(int) * n); 58401a81e61SBarry Smith for (i = 0; i < n; i++) M->block_sizes[i] = 1; 58501a81e61SBarry Smith 5869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&mnl, 1, MPI_INT, M->mnls, 1, MPI_INT, comm)); 58701a81e61SBarry Smith 58801a81e61SBarry Smith M->start_indices[0] = 0; 5892fa5cd67SKarl Rupp for (i = 1; i < size; i++) M->start_indices[i] = M->start_indices[i - 1] + M->mnls[i - 1]; 59001a81e61SBarry Smith 59101a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 59201a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 59301a81e61SBarry Smith 59401a81e61SBarry Smith for (i = 0; i < size; i++) { 59501a81e61SBarry Smith start_indx = M->start_indices[i]; 5962fa5cd67SKarl Rupp for (j = 0; j < M->mnls[i]; j++) M->pe[start_indx + j] = i; 59701a81e61SBarry Smith } 59801a81e61SBarry Smith 59901a81e61SBarry Smith if (AT) { 6002f613bf5SBarry Smith M->lines = new_compressed_lines(M->mnls[rank], 1); 60101a81e61SBarry Smith } else { 6022f613bf5SBarry Smith M->lines = new_compressed_lines(M->mnls[rank], 0); 60301a81e61SBarry Smith } 60401a81e61SBarry Smith 60501a81e61SBarry Smith rows = M->lines; 60601a81e61SBarry Smith 60701a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 6089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M->n, &mapping)); 60901a81e61SBarry Smith pe = 0; 61001a81e61SBarry Smith local_indx = 0; 61101a81e61SBarry Smith for (i = 0; i < M->n; i++) { 61201a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 61301a81e61SBarry Smith pe++; 61401a81e61SBarry Smith local_indx = 0; 61501a81e61SBarry Smith } 61601a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 61701a81e61SBarry Smith local_indx++; 61801a81e61SBarry Smith } 61901a81e61SBarry Smith 62001a81e61SBarry Smith /*********************************************************/ 62101a81e61SBarry Smith /************** Set up the row structure *****************/ 62201a81e61SBarry Smith /*********************************************************/ 62301a81e61SBarry Smith 6249566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 62501a81e61SBarry Smith for (i = rstart; i < rend; i++) { 62601a81e61SBarry Smith row_indx = i - rstart; 6279566063dSJacob Faibussowitsch PetscCall(MatGetRow(A, i, &nz, &cols, &vals)); 6282122902bSSatish Balay /* allocate buffers */ 6292122902bSSatish Balay rows->ptrs[row_indx] = (int *)malloc(nz * sizeof(int)); 6302122902bSSatish Balay rows->A[row_indx] = (double *)malloc(nz * sizeof(double)); 6312122902bSSatish Balay /* copy the matrix */ 63201a81e61SBarry Smith for (j = 0; j < nz; j++) { 63301a81e61SBarry Smith col = cols[j]; 63401a81e61SBarry Smith len = rows->len[row_indx]++; 6352fa5cd67SKarl Rupp 63601a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 63701a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 63801a81e61SBarry Smith } 63901a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 6402fa5cd67SKarl Rupp 6419566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A, i, &nz, &cols, &vals)); 64201a81e61SBarry Smith } 64301a81e61SBarry Smith 64401a81e61SBarry Smith /************************************************************/ 64501a81e61SBarry Smith /************** Set up the column structure *****************/ 64601a81e61SBarry Smith /*********************************************************/ 64701a81e61SBarry Smith 64801a81e61SBarry Smith if (AT) { 64901a81e61SBarry Smith for (i = rstart; i < rend; i++) { 65001a81e61SBarry Smith row_indx = i - rstart; 6519566063dSJacob Faibussowitsch PetscCall(MatGetRow(AT, i, &nz, &cols, &vals)); 6522122902bSSatish Balay /* allocate buffers */ 6532122902bSSatish Balay rows->rptrs[row_indx] = (int *)malloc(nz * sizeof(int)); 6542122902bSSatish Balay /* copy the matrix (i.e., the structure) */ 65501a81e61SBarry Smith for (j = 0; j < nz; j++) { 65601a81e61SBarry Smith col = cols[j]; 65701a81e61SBarry Smith len = rows->rlen[row_indx]++; 6582fa5cd67SKarl Rupp 65901a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 66001a81e61SBarry Smith } 6619566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(AT, i, &nz, &cols, &vals)); 66201a81e61SBarry Smith } 66301a81e61SBarry Smith } 66401a81e61SBarry Smith 6659566063dSJacob Faibussowitsch PetscCall(PetscFree(mapping)); 66601a81e61SBarry Smith 66701a81e61SBarry Smith order_pointers(M); 66801a81e61SBarry Smith M->maxnz = calc_maxnz(M); 66901a81e61SBarry Smith *B = M; 67001a81e61SBarry Smith PetscFunctionReturn(0); 67101a81e61SBarry Smith } 67201a81e61SBarry Smith 67301a81e61SBarry Smith /**********************************************************************/ 67401a81e61SBarry Smith 67501a81e61SBarry Smith /* 67601a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 677df3898eeSBarry Smith This assumes that the SPAI matrix B is stored in 67801a81e61SBarry Smith COMPRESSED-ROW format. 67901a81e61SBarry Smith */ 680*9371c9d4SSatish Balay PetscErrorCode ConvertMatrixToMat(MPI_Comm comm, matrix *B, Mat *PB) { 6814b6b5fe1SBarry Smith PetscMPIInt size, rank; 68201a81e61SBarry Smith int m, n, M, N; 68301a81e61SBarry Smith int d_nz, o_nz; 68401a81e61SBarry Smith int *d_nnz, *o_nnz; 68501a81e61SBarry Smith int i, k, global_row, global_col, first_diag_col, last_diag_col; 68601a81e61SBarry Smith PetscScalar val; 68701a81e61SBarry Smith 68801a81e61SBarry Smith PetscFunctionBegin; 6899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 6909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 69101a81e61SBarry Smith 69201a81e61SBarry Smith m = n = B->mnls[rank]; 69301a81e61SBarry Smith d_nz = o_nz = 0; 69401a81e61SBarry Smith 69506946f3aSJose E. Roman /* Determine preallocation for MatCreateAIJ */ 6969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &d_nnz)); 6979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &o_nnz)); 69801a81e61SBarry Smith for (i = 0; i < m; i++) d_nnz[i] = o_nnz[i] = 0; 69901a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 70001a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 70101a81e61SBarry Smith for (i = 0; i < B->mnls[rank]; i++) { 70201a81e61SBarry Smith for (k = 0; k < B->lines->len[i]; k++) { 70301a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 704db4deed7SKarl Rupp if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++; 705db4deed7SKarl Rupp else o_nnz[i]++; 70601a81e61SBarry Smith } 70701a81e61SBarry Smith } 70801a81e61SBarry Smith 70901a81e61SBarry Smith M = N = B->n; 71001a81e61SBarry Smith /* Here we only know how to create AIJ format */ 7119566063dSJacob Faibussowitsch PetscCall(MatCreate(comm, PB)); 7129566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*PB, m, n, M, N)); 7139566063dSJacob Faibussowitsch PetscCall(MatSetType(*PB, MATAIJ)); 7149566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*PB, d_nz, d_nnz)); 7159566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*PB, d_nz, d_nnz, o_nz, o_nnz)); 71601a81e61SBarry Smith 71701a81e61SBarry Smith for (i = 0; i < B->mnls[rank]; i++) { 71801a81e61SBarry Smith global_row = B->start_indices[rank] + i; 71901a81e61SBarry Smith for (k = 0; k < B->lines->len[i]; k++) { 72001a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 7212fa5cd67SKarl Rupp 72201a81e61SBarry Smith val = B->lines->A[i][k]; 7239566063dSJacob Faibussowitsch PetscCall(MatSetValues(*PB, 1, &global_row, 1, &global_col, &val, ADD_VALUES)); 72401a81e61SBarry Smith } 72501a81e61SBarry Smith } 72601a81e61SBarry Smith 7279566063dSJacob Faibussowitsch PetscCall(PetscFree(d_nnz)); 7289566063dSJacob Faibussowitsch PetscCall(PetscFree(o_nnz)); 72901a81e61SBarry Smith 7309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(*PB, MAT_FINAL_ASSEMBLY)); 7319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(*PB, MAT_FINAL_ASSEMBLY)); 73201a81e61SBarry Smith PetscFunctionReturn(0); 73301a81e61SBarry Smith } 73401a81e61SBarry Smith 73501a81e61SBarry Smith /**********************************************************************/ 73601a81e61SBarry Smith 73701a81e61SBarry Smith /* 73801a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 73901a81e61SBarry Smith */ 740*9371c9d4SSatish Balay PetscErrorCode ConvertVectorToVec(MPI_Comm comm, vector *v, Vec *Pv) { 7414b6b5fe1SBarry Smith PetscMPIInt size, rank; 7424b6b5fe1SBarry Smith int m, M, i, *mnls, *start_indices, *global_indices; 74301a81e61SBarry Smith 74401a81e61SBarry Smith PetscFunctionBegin; 7459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 7469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 74701a81e61SBarry Smith 74801a81e61SBarry Smith m = v->mnl; 74901a81e61SBarry Smith M = v->n; 75001a81e61SBarry Smith 7519566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm, m, M, Pv)); 75201a81e61SBarry Smith 7539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &mnls)); 7549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allgather(&v->mnl, 1, MPI_INT, mnls, 1, MPI_INT, comm)); 75501a81e61SBarry Smith 7569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(size, &start_indices)); 7572fa5cd67SKarl Rupp 75801a81e61SBarry Smith start_indices[0] = 0; 7592fa5cd67SKarl Rupp for (i = 1; i < size; i++) start_indices[i] = start_indices[i - 1] + mnls[i - 1]; 76001a81e61SBarry Smith 7619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(v->mnl, &global_indices)); 7622fa5cd67SKarl Rupp for (i = 0; i < v->mnl; i++) global_indices[i] = start_indices[rank] + i; 76301a81e61SBarry Smith 7649566063dSJacob Faibussowitsch PetscCall(PetscFree(mnls)); 7659566063dSJacob Faibussowitsch PetscCall(PetscFree(start_indices)); 76601a81e61SBarry Smith 7679566063dSJacob Faibussowitsch PetscCall(VecSetValues(*Pv, v->mnl, global_indices, v->v, INSERT_VALUES)); 7689566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(*Pv)); 7699566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(*Pv)); 77001a81e61SBarry Smith 7719566063dSJacob Faibussowitsch PetscCall(PetscFree(global_indices)); 77201a81e61SBarry Smith PetscFunctionReturn(0); 77301a81e61SBarry Smith } 774