101a81e61SBarry Smith #define PETSCKSP_DLL 201a81e61SBarry Smith 301a81e61SBarry Smith /* 401a81e61SBarry Smith 3/99 Modified by Stephen Barnard to support SPAI version 3.0 501a81e61SBarry Smith */ 601a81e61SBarry Smith 701a81e61SBarry Smith /* 801a81e61SBarry Smith Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner 901a81e61SBarry Smith Code written by Stephen Barnard. 1001a81e61SBarry Smith 1101a81e61SBarry Smith Note: there is some BAD memory bleeding below! 1201a81e61SBarry Smith 1301a81e61SBarry Smith This code needs work 1401a81e61SBarry Smith 1501a81e61SBarry Smith 1) get rid of all memory bleeding 1601a81e61SBarry Smith 2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix 1701a81e61SBarry Smith rather than having the sp flag for PC_SPAI 1801a81e61SBarry Smith 1901a81e61SBarry Smith */ 2001a81e61SBarry Smith 216356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 2275031707SSatish Balay #include "petscspai.h" 2301a81e61SBarry Smith 2401a81e61SBarry Smith /* 2501a81e61SBarry Smith These are the SPAI include files 2601a81e61SBarry Smith */ 2701a81e61SBarry Smith EXTERN_C_BEGIN 2801a81e61SBarry Smith #include "spai.h" 2901a81e61SBarry Smith #include "matrix.h" 3001a81e61SBarry Smith EXTERN_C_END 3101a81e61SBarry Smith 3201a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**); 3301a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix *,Mat *); 3401a81e61SBarry Smith EXTERN PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv); 3501a81e61SBarry Smith EXTERN PetscErrorCode MM_to_PETSC(char *,char *,char *); 3601a81e61SBarry Smith 3701a81e61SBarry Smith typedef struct { 3801a81e61SBarry Smith 3901a81e61SBarry Smith matrix *B; /* matrix in SPAI format */ 4001a81e61SBarry Smith matrix *BT; /* transpose of matrix in SPAI format */ 4101a81e61SBarry Smith matrix *M; /* the approximate inverse in SPAI format */ 4201a81e61SBarry Smith 4301a81e61SBarry Smith Mat PM; /* the approximate inverse PETSc format */ 4401a81e61SBarry Smith 4501a81e61SBarry Smith double epsilon; /* tolerance */ 4601a81e61SBarry Smith int nbsteps; /* max number of "improvement" steps per line */ 4701a81e61SBarry Smith int max; /* max dimensions of is_I, q, etc. */ 4801a81e61SBarry Smith int maxnew; /* max number of new entries per step */ 4901a81e61SBarry Smith int block_size; /* constant block size */ 5001a81e61SBarry Smith int cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */ 5101a81e61SBarry Smith int verbose; /* SPAI prints timing and statistics */ 5201a81e61SBarry Smith 5301a81e61SBarry Smith int sp; /* symmetric nonzero pattern */ 5401a81e61SBarry Smith MPI_Comm comm_spai; /* communicator to be used with spai */ 5501a81e61SBarry Smith } PC_SPAI; 5601a81e61SBarry Smith 5701a81e61SBarry Smith /**********************************************************************/ 5801a81e61SBarry Smith 5901a81e61SBarry Smith #undef __FUNCT__ 6001a81e61SBarry Smith #define __FUNCT__ "PCSetUp_SPAI" 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 6901a81e61SBarry Smith init_SPAI(); 7001a81e61SBarry Smith 7101a81e61SBarry Smith if (ispai->sp) { 7201a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr); 7301a81e61SBarry Smith } else { 7401a81e61SBarry Smith /* Use the transpose to get the column nonzero structure. */ 7501a81e61SBarry Smith ierr = MatTranspose(pc->pmat,&AT);CHKERRQ(ierr); 7601a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr); 7701a81e61SBarry Smith ierr = MatDestroy(AT);CHKERRQ(ierr); 7801a81e61SBarry Smith } 7901a81e61SBarry Smith 8001a81e61SBarry Smith /* Destroy the transpose */ 8101a81e61SBarry Smith /* Don't know how to do it. PETSc developers? */ 8201a81e61SBarry Smith 8301a81e61SBarry Smith /* construct SPAI preconditioner */ 8401a81e61SBarry Smith /* FILE *messages */ /* file for warning messages */ 8501a81e61SBarry Smith /* double epsilon */ /* tolerance */ 8601a81e61SBarry Smith /* int nbsteps */ /* max number of "improvement" steps per line */ 8701a81e61SBarry Smith /* int max */ /* max dimensions of is_I, q, etc. */ 8801a81e61SBarry Smith /* int maxnew */ /* max number of new entries per step */ 8901a81e61SBarry Smith /* int block_size */ /* block_size == 1 specifies scalar elments 9001a81e61SBarry Smith block_size == n specifies nxn constant-block elements 9101a81e61SBarry Smith block_size == 0 specifies variable-block elements */ 9201a81e61SBarry Smith /* int cache_size */ /* one of (1,2,3,4,5,6) indicting size of cache */ 9301a81e61SBarry Smith /* cache_size == 0 indicates no caching */ 9401a81e61SBarry Smith /* int verbose */ /* verbose == 0 specifies that SPAI is silent 9501a81e61SBarry Smith verbose == 1 prints timing and matrix statistics */ 9601a81e61SBarry Smith 9701a81e61SBarry Smith ierr = bspai(ispai->B,&ispai->M, 9801a81e61SBarry Smith stdout, 9901a81e61SBarry Smith ispai->epsilon, 10001a81e61SBarry Smith ispai->nbsteps, 10101a81e61SBarry Smith ispai->max, 10201a81e61SBarry Smith ispai->maxnew, 10301a81e61SBarry Smith ispai->block_size, 10401a81e61SBarry Smith ispai->cache_size, 10501a81e61SBarry Smith ispai->verbose); CHKERRQ(ierr); 10601a81e61SBarry Smith 10701a81e61SBarry Smith ierr = ConvertMatrixToMat(pc->comm,ispai->M,&ispai->PM);CHKERRQ(ierr); 10801a81e61SBarry Smith 10901a81e61SBarry Smith /* free the SPAI matrices */ 11001a81e61SBarry Smith sp_free_matrix(ispai->B); 11101a81e61SBarry Smith sp_free_matrix(ispai->M); 11201a81e61SBarry Smith 11301a81e61SBarry Smith PetscFunctionReturn(0); 11401a81e61SBarry Smith } 11501a81e61SBarry Smith 11601a81e61SBarry Smith /**********************************************************************/ 11701a81e61SBarry Smith 11801a81e61SBarry Smith #undef __FUNCT__ 11901a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI" 12001a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y) 12101a81e61SBarry Smith { 12201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 12301a81e61SBarry Smith PetscErrorCode ierr; 12401a81e61SBarry Smith 12501a81e61SBarry Smith PetscFunctionBegin; 12601a81e61SBarry Smith /* Now using PETSc's multiply */ 12701a81e61SBarry Smith ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr); 12801a81e61SBarry Smith PetscFunctionReturn(0); 12901a81e61SBarry Smith } 13001a81e61SBarry Smith 13101a81e61SBarry Smith /**********************************************************************/ 13201a81e61SBarry Smith 13301a81e61SBarry Smith #undef __FUNCT__ 13401a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI" 13501a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc) 13601a81e61SBarry Smith { 13701a81e61SBarry Smith PetscErrorCode ierr; 13801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 13901a81e61SBarry Smith 14001a81e61SBarry Smith PetscFunctionBegin; 14101a81e61SBarry Smith if (ispai->PM) {ierr = MatDestroy(ispai->PM);CHKERRQ(ierr);} 14201a81e61SBarry Smith ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr); 14301a81e61SBarry Smith ierr = PetscFree(ispai);CHKERRQ(ierr); 14401a81e61SBarry Smith PetscFunctionReturn(0); 14501a81e61SBarry Smith } 14601a81e61SBarry Smith 14701a81e61SBarry Smith /**********************************************************************/ 14801a81e61SBarry Smith 14901a81e61SBarry Smith #undef __FUNCT__ 15001a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI" 15101a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer) 15201a81e61SBarry Smith { 15301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 15401a81e61SBarry Smith PetscErrorCode ierr; 15501a81e61SBarry Smith PetscTruth iascii; 15601a81e61SBarry Smith 15701a81e61SBarry Smith PetscFunctionBegin; 15801a81e61SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 15901a81e61SBarry Smith if (iascii) { 16001a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SPAI preconditioner\n");CHKERRQ(ierr); 161a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," epsilon %G\n", ispai->epsilon);CHKERRQ(ierr); 16201a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," nbsteps %d\n", ispai->nbsteps);CHKERRQ(ierr); 16301a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," max %d\n", ispai->max);CHKERRQ(ierr); 16401a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxnew %d\n", ispai->maxnew);CHKERRQ(ierr); 16501a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block_size %d\n",ispai->block_size);CHKERRQ(ierr); 16601a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," cache_size %d\n",ispai->cache_size);CHKERRQ(ierr); 16701a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," verbose %d\n", ispai->verbose);CHKERRQ(ierr); 16801a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," sp %d\n", ispai->sp);CHKERRQ(ierr); 16901a81e61SBarry Smith 17001a81e61SBarry Smith } 17101a81e61SBarry Smith PetscFunctionReturn(0); 17201a81e61SBarry Smith } 17301a81e61SBarry Smith 17401a81e61SBarry Smith EXTERN_C_BEGIN 17501a81e61SBarry Smith #undef __FUNCT__ 17601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI" 17701a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon_SPAI(PC pc,double epsilon1) 17801a81e61SBarry Smith { 17901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 18001a81e61SBarry Smith PetscFunctionBegin; 18101a81e61SBarry Smith ispai->epsilon = epsilon1; 18201a81e61SBarry Smith PetscFunctionReturn(0); 18301a81e61SBarry Smith } 18401a81e61SBarry Smith EXTERN_C_END 18501a81e61SBarry Smith 18601a81e61SBarry Smith /**********************************************************************/ 18701a81e61SBarry Smith 18801a81e61SBarry Smith EXTERN_C_BEGIN 18901a81e61SBarry Smith #undef __FUNCT__ 19001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI" 19101a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1) 19201a81e61SBarry Smith { 19301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 19401a81e61SBarry Smith PetscFunctionBegin; 19501a81e61SBarry Smith ispai->nbsteps = nbsteps1; 19601a81e61SBarry Smith PetscFunctionReturn(0); 19701a81e61SBarry Smith } 19801a81e61SBarry Smith EXTERN_C_END 19901a81e61SBarry Smith 20001a81e61SBarry Smith /**********************************************************************/ 20101a81e61SBarry Smith 20201a81e61SBarry Smith /* added 1/7/99 g.h. */ 20301a81e61SBarry Smith EXTERN_C_BEGIN 20401a81e61SBarry Smith #undef __FUNCT__ 20501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax_SPAI" 20601a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax_SPAI(PC pc,int max1) 20701a81e61SBarry Smith { 20801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 20901a81e61SBarry Smith PetscFunctionBegin; 21001a81e61SBarry Smith ispai->max = max1; 21101a81e61SBarry Smith PetscFunctionReturn(0); 21201a81e61SBarry Smith } 21301a81e61SBarry Smith EXTERN_C_END 21401a81e61SBarry Smith 21501a81e61SBarry Smith /**********************************************************************/ 21601a81e61SBarry Smith 21701a81e61SBarry Smith EXTERN_C_BEGIN 21801a81e61SBarry Smith #undef __FUNCT__ 21901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI" 22001a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew_SPAI(PC pc,int maxnew1) 22101a81e61SBarry Smith { 22201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 22301a81e61SBarry Smith PetscFunctionBegin; 22401a81e61SBarry Smith ispai->maxnew = maxnew1; 22501a81e61SBarry Smith PetscFunctionReturn(0); 22601a81e61SBarry Smith } 22701a81e61SBarry Smith EXTERN_C_END 22801a81e61SBarry Smith 22901a81e61SBarry Smith /**********************************************************************/ 23001a81e61SBarry Smith 23101a81e61SBarry Smith EXTERN_C_BEGIN 23201a81e61SBarry Smith #undef __FUNCT__ 23301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize_SPAI" 23401a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize_SPAI(PC pc,int block_size1) 23501a81e61SBarry Smith { 23601a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 23701a81e61SBarry Smith PetscFunctionBegin; 23801a81e61SBarry Smith ispai->block_size = block_size1; 23901a81e61SBarry Smith PetscFunctionReturn(0); 24001a81e61SBarry Smith } 24101a81e61SBarry Smith EXTERN_C_END 24201a81e61SBarry Smith 24301a81e61SBarry Smith /**********************************************************************/ 24401a81e61SBarry Smith 24501a81e61SBarry Smith EXTERN_C_BEGIN 24601a81e61SBarry Smith #undef __FUNCT__ 24701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI" 24801a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize_SPAI(PC pc,int cache_size) 24901a81e61SBarry Smith { 25001a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 25101a81e61SBarry Smith PetscFunctionBegin; 25201a81e61SBarry Smith ispai->cache_size = cache_size; 25301a81e61SBarry Smith PetscFunctionReturn(0); 25401a81e61SBarry Smith } 25501a81e61SBarry Smith EXTERN_C_END 25601a81e61SBarry Smith 25701a81e61SBarry Smith /**********************************************************************/ 25801a81e61SBarry Smith 25901a81e61SBarry Smith EXTERN_C_BEGIN 26001a81e61SBarry Smith #undef __FUNCT__ 26101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI" 26201a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose_SPAI(PC pc,int verbose) 26301a81e61SBarry Smith { 26401a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 26501a81e61SBarry Smith PetscFunctionBegin; 26601a81e61SBarry Smith ispai->verbose = verbose; 26701a81e61SBarry Smith PetscFunctionReturn(0); 26801a81e61SBarry Smith } 26901a81e61SBarry Smith EXTERN_C_END 27001a81e61SBarry Smith 27101a81e61SBarry Smith /**********************************************************************/ 27201a81e61SBarry Smith 27301a81e61SBarry Smith EXTERN_C_BEGIN 27401a81e61SBarry Smith #undef __FUNCT__ 27501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI" 27601a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp_SPAI(PC pc,int sp) 27701a81e61SBarry Smith { 27801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 27901a81e61SBarry Smith PetscFunctionBegin; 28001a81e61SBarry Smith ispai->sp = sp; 28101a81e61SBarry Smith PetscFunctionReturn(0); 28201a81e61SBarry Smith } 28301a81e61SBarry Smith EXTERN_C_END 28401a81e61SBarry Smith 28501a81e61SBarry Smith /* -------------------------------------------------------------------*/ 28601a81e61SBarry Smith 28701a81e61SBarry Smith #undef __FUNCT__ 28801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon" 28901a81e61SBarry Smith /*@ 29001a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 29101a81e61SBarry Smith 29201a81e61SBarry Smith Input Parameters: 29301a81e61SBarry Smith + pc - the preconditioner 29401a81e61SBarry Smith - eps - epsilon (default .4) 29501a81e61SBarry Smith 29601a81e61SBarry Smith Notes: Espilon must be between 0 and 1. It controls the 29701a81e61SBarry Smith quality of the approximation of M to the inverse of 29801a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 29901a81e61SBarry Smith fill, and usually better preconditioners. In many 30001a81e61SBarry Smith cases the best choice of epsilon is the one that 30101a81e61SBarry Smith divides the total solution time equally between the 30201a81e61SBarry Smith preconditioner and the solver. 30301a81e61SBarry Smith 30401a81e61SBarry Smith Level: intermediate 30501a81e61SBarry Smith 30601a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 30701a81e61SBarry Smith @*/ 30801a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC pc,double epsilon1) 30901a81e61SBarry Smith { 31001a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,double); 31101a81e61SBarry Smith PetscFunctionBegin; 31201a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetEpsilon_C",(void (**)(void))&f);CHKERRQ(ierr); 31301a81e61SBarry Smith if (f) { 31401a81e61SBarry Smith ierr = (*f)(pc,epsilon1);CHKERRQ(ierr); 31501a81e61SBarry Smith } 31601a81e61SBarry Smith PetscFunctionReturn(0); 31701a81e61SBarry Smith } 31801a81e61SBarry Smith 31901a81e61SBarry Smith /**********************************************************************/ 32001a81e61SBarry Smith 32101a81e61SBarry Smith #undef __FUNCT__ 32201a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps" 32301a81e61SBarry Smith /*@ 32401a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 32501a81e61SBarry Smith the SPAI preconditioner 32601a81e61SBarry Smith 32701a81e61SBarry Smith Input Parameters: 32801a81e61SBarry Smith + pc - the preconditioner 32901a81e61SBarry Smith - n - number of steps (default 5) 33001a81e61SBarry Smith 33101a81e61SBarry Smith Notes: SPAI constructs to approximation to every column of 33201a81e61SBarry Smith the exact inverse of A in a series of improvement 33301a81e61SBarry Smith steps. The quality of the approximation is determined 33401a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 33501a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 33601a81e61SBarry Smith uses the best approximation constructed so far. 33701a81e61SBarry Smith 33801a81e61SBarry Smith Level: intermediate 33901a81e61SBarry Smith 34001a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew() 34101a81e61SBarry Smith @*/ 34201a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC pc,int nbsteps1) 34301a81e61SBarry Smith { 34401a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 34501a81e61SBarry Smith PetscFunctionBegin; 34601a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetNBSteps_C",(void (**)(void))&f);CHKERRQ(ierr); 34701a81e61SBarry Smith if (f) { 34801a81e61SBarry Smith ierr = (*f)(pc,nbsteps1);CHKERRQ(ierr); 34901a81e61SBarry Smith } 35001a81e61SBarry Smith PetscFunctionReturn(0); 35101a81e61SBarry Smith } 35201a81e61SBarry Smith 35301a81e61SBarry Smith /**********************************************************************/ 35401a81e61SBarry Smith 35501a81e61SBarry Smith /* added 1/7/99 g.h. */ 35601a81e61SBarry Smith #undef __FUNCT__ 35701a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax" 35801a81e61SBarry Smith /*@ 35901a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 36001a81e61SBarry Smith the SPAI preconditioner 36101a81e61SBarry Smith 36201a81e61SBarry Smith Input Parameters: 36301a81e61SBarry Smith + pc - the preconditioner 36401a81e61SBarry Smith - n - size (default is 5000) 36501a81e61SBarry Smith 36601a81e61SBarry Smith Level: intermediate 36701a81e61SBarry Smith 36801a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 36901a81e61SBarry Smith @*/ 37001a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC pc,int max1) 37101a81e61SBarry Smith { 37201a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 37301a81e61SBarry Smith PetscFunctionBegin; 37401a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMax_C",(void (**)(void))&f);CHKERRQ(ierr); 37501a81e61SBarry Smith if (f) { 37601a81e61SBarry Smith ierr = (*f)(pc,max1);CHKERRQ(ierr); 37701a81e61SBarry Smith } 37801a81e61SBarry Smith PetscFunctionReturn(0); 37901a81e61SBarry Smith } 38001a81e61SBarry Smith 38101a81e61SBarry Smith /**********************************************************************/ 38201a81e61SBarry Smith 38301a81e61SBarry Smith #undef __FUNCT__ 38401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew" 38501a81e61SBarry Smith /*@ 38601a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 38701a81e61SBarry Smith in SPAI preconditioner 38801a81e61SBarry Smith 38901a81e61SBarry Smith Input Parameters: 39001a81e61SBarry Smith + pc - the preconditioner 39101a81e61SBarry Smith - n - maximum number (default 5) 39201a81e61SBarry Smith 39301a81e61SBarry Smith Level: intermediate 39401a81e61SBarry Smith 39501a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps() 39601a81e61SBarry Smith @*/ 39701a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC pc,int maxnew1) 39801a81e61SBarry Smith { 39901a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 40001a81e61SBarry Smith PetscFunctionBegin; 40101a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMaxNew_C",(void (**)(void))&f);CHKERRQ(ierr); 40201a81e61SBarry Smith if (f) { 40301a81e61SBarry Smith ierr = (*f)(pc,maxnew1);CHKERRQ(ierr); 40401a81e61SBarry Smith } 40501a81e61SBarry Smith PetscFunctionReturn(0); 40601a81e61SBarry Smith } 40701a81e61SBarry Smith 40801a81e61SBarry Smith /**********************************************************************/ 40901a81e61SBarry Smith 41001a81e61SBarry Smith #undef __FUNCT__ 41101a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize" 41201a81e61SBarry Smith /*@ 41301a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 41401a81e61SBarry Smith 41501a81e61SBarry Smith Input Parameters: 41601a81e61SBarry Smith + pc - the preconditioner 41701a81e61SBarry Smith - n - block size (default 1) 41801a81e61SBarry Smith 41901a81e61SBarry Smith Notes: A block 42001a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 42101a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 42201a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 42301a81e61SBarry Smith variable sized blocks, which are determined by 42401a81e61SBarry Smith searching for dense square diagonal blocks in A. 42501a81e61SBarry Smith This can be very effective for finite-element 42601a81e61SBarry Smith matrices. 42701a81e61SBarry Smith 42801a81e61SBarry Smith SPAI will convert A to block form, use a block 42901a81e61SBarry Smith version of the preconditioner algorithm, and then 43001a81e61SBarry Smith convert the result back to scalar form. 43101a81e61SBarry Smith 43201a81e61SBarry Smith In many cases the a block-size parameter other than 1 43301a81e61SBarry Smith can lead to very significant improvement in 43401a81e61SBarry Smith performance. 43501a81e61SBarry Smith 43601a81e61SBarry Smith 43701a81e61SBarry Smith Level: intermediate 43801a81e61SBarry Smith 43901a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 44001a81e61SBarry Smith @*/ 44101a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC pc,int block_size1) 44201a81e61SBarry Smith { 44301a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 44401a81e61SBarry Smith PetscFunctionBegin; 44501a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 44601a81e61SBarry Smith if (f) { 44701a81e61SBarry Smith ierr = (*f)(pc,block_size1);CHKERRQ(ierr); 44801a81e61SBarry Smith } 44901a81e61SBarry Smith PetscFunctionReturn(0); 45001a81e61SBarry Smith } 45101a81e61SBarry Smith 45201a81e61SBarry Smith /**********************************************************************/ 45301a81e61SBarry Smith 45401a81e61SBarry Smith #undef __FUNCT__ 45501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize" 45601a81e61SBarry Smith /*@ 45701a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 45801a81e61SBarry Smith 45901a81e61SBarry Smith Input Parameters: 46001a81e61SBarry Smith + pc - the preconditioner 46101a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 46201a81e61SBarry Smith 46301a81e61SBarry Smith Notes: SPAI uses a hash table to cache messages and avoid 46401a81e61SBarry Smith redundant communication. If suggest always using 46501a81e61SBarry Smith 5. This parameter is irrelevant in the serial 46601a81e61SBarry Smith version. 46701a81e61SBarry Smith 46801a81e61SBarry Smith Level: intermediate 46901a81e61SBarry Smith 47001a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 47101a81e61SBarry Smith @*/ 47201a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC pc,int cache_size) 47301a81e61SBarry Smith { 47401a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 47501a81e61SBarry Smith PetscFunctionBegin; 47601a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetCacheSize_C",(void (**)(void))&f);CHKERRQ(ierr); 47701a81e61SBarry Smith if (f) { 47801a81e61SBarry Smith ierr = (*f)(pc,cache_size);CHKERRQ(ierr); 47901a81e61SBarry Smith } 48001a81e61SBarry Smith PetscFunctionReturn(0); 48101a81e61SBarry Smith } 48201a81e61SBarry Smith 48301a81e61SBarry Smith /**********************************************************************/ 48401a81e61SBarry Smith 48501a81e61SBarry Smith #undef __FUNCT__ 48601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose" 48701a81e61SBarry Smith /*@ 48801a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 48901a81e61SBarry Smith 49001a81e61SBarry Smith Input Parameters: 49101a81e61SBarry Smith + pc - the preconditioner 49201a81e61SBarry Smith - n - level (default 1) 49301a81e61SBarry Smith 49401a81e61SBarry Smith Notes: print parameters, timings and matrix statistics 49501a81e61SBarry Smith 49601a81e61SBarry Smith Level: intermediate 49701a81e61SBarry Smith 49801a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 49901a81e61SBarry Smith @*/ 50001a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC pc,int verbose) 50101a81e61SBarry Smith { 50201a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 50301a81e61SBarry Smith PetscFunctionBegin; 50401a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetVerbose_C",(void (**)(void))&f);CHKERRQ(ierr); 50501a81e61SBarry Smith if (f) { 50601a81e61SBarry Smith ierr = (*f)(pc,verbose);CHKERRQ(ierr); 50701a81e61SBarry Smith } 50801a81e61SBarry Smith PetscFunctionReturn(0); 50901a81e61SBarry Smith } 51001a81e61SBarry Smith 51101a81e61SBarry Smith /**********************************************************************/ 51201a81e61SBarry Smith 51301a81e61SBarry Smith #undef __FUNCT__ 51401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp" 51501a81e61SBarry Smith /*@ 51601a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 51701a81e61SBarry Smith 51801a81e61SBarry Smith Input Parameters: 51901a81e61SBarry Smith + pc - the preconditioner 52001a81e61SBarry Smith - n - 0 or 1 52101a81e61SBarry Smith 52201a81e61SBarry Smith Notes: If A has a symmetric nonzero pattern use -sp 1 to 52301a81e61SBarry Smith improve performance by eliminating some communication 52401a81e61SBarry Smith in the parallel version. Even if A does not have a 52501a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 52601a81e61SBarry Smith results, but the code will not follow the published 52701a81e61SBarry Smith SPAI algorithm exactly. 52801a81e61SBarry Smith 52901a81e61SBarry Smith 53001a81e61SBarry Smith Level: intermediate 53101a81e61SBarry Smith 53201a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 53301a81e61SBarry Smith @*/ 53401a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC pc,int sp) 53501a81e61SBarry Smith { 53601a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 53701a81e61SBarry Smith PetscFunctionBegin; 53801a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetSp_C",(void (**)(void))&f);CHKERRQ(ierr); 53901a81e61SBarry Smith if (f) { 54001a81e61SBarry Smith ierr = (*f)(pc,sp);CHKERRQ(ierr); 54101a81e61SBarry Smith } 54201a81e61SBarry Smith PetscFunctionReturn(0); 54301a81e61SBarry Smith } 54401a81e61SBarry Smith 54501a81e61SBarry Smith /**********************************************************************/ 54601a81e61SBarry Smith 54701a81e61SBarry Smith /**********************************************************************/ 54801a81e61SBarry Smith 54901a81e61SBarry Smith #undef __FUNCT__ 55001a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI" 55101a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc) 55201a81e61SBarry Smith { 55301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 55401a81e61SBarry Smith PetscErrorCode ierr; 55501a81e61SBarry Smith int nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp; 55601a81e61SBarry Smith double epsilon1; 55701a81e61SBarry Smith PetscTruth flg; 55801a81e61SBarry Smith 55901a81e61SBarry Smith PetscFunctionBegin; 56001a81e61SBarry Smith ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr); 56101a81e61SBarry Smith ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr); 56201a81e61SBarry Smith if (flg) { 56301a81e61SBarry Smith ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr); 56401a81e61SBarry Smith } 56501a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr); 56601a81e61SBarry Smith if (flg) { 56701a81e61SBarry Smith ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr); 56801a81e61SBarry Smith } 56901a81e61SBarry Smith /* added 1/7/99 g.h. */ 57001a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr); 57101a81e61SBarry Smith if (flg) { 57201a81e61SBarry Smith ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr); 57301a81e61SBarry Smith } 57401a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr); 57501a81e61SBarry Smith if (flg) { 57601a81e61SBarry Smith ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr); 57701a81e61SBarry Smith } 57801a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr); 57901a81e61SBarry Smith if (flg) { 58001a81e61SBarry Smith ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr); 58101a81e61SBarry Smith } 58201a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr); 58301a81e61SBarry Smith if (flg) { 58401a81e61SBarry Smith ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr); 58501a81e61SBarry Smith } 58601a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr); 58701a81e61SBarry Smith if (flg) { 58801a81e61SBarry Smith ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr); 58901a81e61SBarry Smith } 59001a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr); 59101a81e61SBarry Smith if (flg) { 59201a81e61SBarry Smith ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr); 59301a81e61SBarry Smith } 59401a81e61SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 59501a81e61SBarry Smith PetscFunctionReturn(0); 59601a81e61SBarry Smith } 59701a81e61SBarry Smith 59801a81e61SBarry Smith /**********************************************************************/ 59901a81e61SBarry Smith 60001a81e61SBarry Smith /*MC 60101a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 60201a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 60301a81e61SBarry Smith 60401a81e61SBarry Smith Options Database Keys: 60501a81e61SBarry Smith + -pc_spai_set_epsilon <eps> - set tolerance 60601a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 60701a81e61SBarry Smith . -pc_spai_max <m> - set max 60801a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 60901a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 61001a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 61101a81e61SBarry Smith . -pc_spai_sp <m> - set sp 61201a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 61301a81e61SBarry Smith 61401a81e61SBarry Smith Notes: This only works with AIJ matrices. 61501a81e61SBarry Smith 61601a81e61SBarry Smith Level: beginner 61701a81e61SBarry Smith 61801a81e61SBarry Smith Concepts: approximate inverse 61901a81e61SBarry Smith 62001a81e61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 62101a81e61SBarry Smith PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(), 62201a81e61SBarry Smith PCSPAISetVerbose(), PCSPAISetSp() 62301a81e61SBarry Smith M*/ 62401a81e61SBarry Smith 62501a81e61SBarry Smith EXTERN_C_BEGIN 62601a81e61SBarry Smith #undef __FUNCT__ 62701a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI" 62801a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SPAI(PC pc) 62901a81e61SBarry Smith { 63001a81e61SBarry Smith PC_SPAI *ispai; 63101a81e61SBarry Smith PetscErrorCode ierr; 63201a81e61SBarry Smith 63301a81e61SBarry Smith PetscFunctionBegin; 634*38f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_SPAI,&ispai);CHKERRQ(ierr); 63501a81e61SBarry Smith pc->data = (void*)ispai; 63601a81e61SBarry Smith 63701a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 63801a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 63901a81e61SBarry Smith pc->ops->applyrichardson = 0; 64001a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 64101a81e61SBarry Smith pc->ops->view = PCView_SPAI; 64201a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 64301a81e61SBarry Smith 64401a81e61SBarry Smith pc->name = 0; 64501a81e61SBarry Smith ispai->epsilon = .4; 64601a81e61SBarry Smith ispai->nbsteps = 5; 64701a81e61SBarry Smith ispai->max = 5000; 64801a81e61SBarry Smith ispai->maxnew = 5; 64901a81e61SBarry Smith ispai->block_size = 1; 65001a81e61SBarry Smith ispai->cache_size = 5; 65101a81e61SBarry Smith ispai->verbose = 0; 65201a81e61SBarry Smith 65301a81e61SBarry Smith ispai->sp = 1; 65401a81e61SBarry Smith ierr = MPI_Comm_dup(pc->comm,&(ispai->comm_spai));CHKERRQ(ierr); 65501a81e61SBarry Smith 65601a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C", 65701a81e61SBarry Smith "PCSPAISetEpsilon_SPAI", 65801a81e61SBarry Smith PCSPAISetEpsilon_SPAI);CHKERRQ(ierr); 65901a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C", 66001a81e61SBarry Smith "PCSPAISetNBSteps_SPAI", 66101a81e61SBarry Smith PCSPAISetNBSteps_SPAI);CHKERRQ(ierr); 66201a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C", 66301a81e61SBarry Smith "PCSPAISetMax_SPAI", 66401a81e61SBarry Smith PCSPAISetMax_SPAI);CHKERRQ(ierr); 66501a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC", 66601a81e61SBarry Smith "PCSPAISetMaxNew_SPAI", 66701a81e61SBarry Smith PCSPAISetMaxNew_SPAI);CHKERRQ(ierr); 66801a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C", 66901a81e61SBarry Smith "PCSPAISetBlockSize_SPAI", 67001a81e61SBarry Smith PCSPAISetBlockSize_SPAI);CHKERRQ(ierr); 67101a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C", 67201a81e61SBarry Smith "PCSPAISetCacheSize_SPAI", 67301a81e61SBarry Smith PCSPAISetCacheSize_SPAI);CHKERRQ(ierr); 67401a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C", 67501a81e61SBarry Smith "PCSPAISetVerbose_SPAI", 67601a81e61SBarry Smith PCSPAISetVerbose_SPAI);CHKERRQ(ierr); 67701a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C", 67801a81e61SBarry Smith "PCSPAISetSp_SPAI", 67901a81e61SBarry Smith PCSPAISetSp_SPAI);CHKERRQ(ierr); 68001a81e61SBarry Smith 68101a81e61SBarry Smith PetscFunctionReturn(0); 68201a81e61SBarry Smith } 68301a81e61SBarry Smith EXTERN_C_END 68401a81e61SBarry Smith 68501a81e61SBarry Smith /**********************************************************************/ 68601a81e61SBarry Smith 68701a81e61SBarry Smith /* 68801a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 68901a81e61SBarry Smith */ 69001a81e61SBarry Smith #undef __FUNCT__ 69101a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix" 69201a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B) 69301a81e61SBarry Smith { 69401a81e61SBarry Smith matrix *M; 69501a81e61SBarry Smith int i,j,col; 69601a81e61SBarry Smith int row_indx; 69701a81e61SBarry Smith int len,pe,local_indx,start_indx; 69801a81e61SBarry Smith int *mapping; 69901a81e61SBarry Smith PetscErrorCode ierr; 70001a81e61SBarry Smith const int *cols; 70101a81e61SBarry Smith const double *vals; 70201a81e61SBarry Smith int *num_ptr,n,mnl,nnl,rank,size,nz,rstart,rend; 70301a81e61SBarry Smith struct compressed_lines *rows; 70401a81e61SBarry Smith 70501a81e61SBarry Smith PetscFunctionBegin; 70601a81e61SBarry Smith 70701a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 70801a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 70901a81e61SBarry Smith ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr); 71001a81e61SBarry Smith ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr); 71101a81e61SBarry Smith 71201a81e61SBarry Smith /* 71301a81e61SBarry Smith not sure why a barrier is required. commenting out 71401a81e61SBarry Smith ierr = MPI_Barrier(comm);CHKERRQ(ierr); 71501a81e61SBarry Smith */ 71601a81e61SBarry Smith 71701a81e61SBarry Smith M = new_matrix((void*)comm); 71801a81e61SBarry Smith 71901a81e61SBarry Smith M->n = n; 72001a81e61SBarry Smith M->bs = 1; 72101a81e61SBarry Smith M->max_block_size = 1; 72201a81e61SBarry Smith 72301a81e61SBarry Smith M->mnls = (int*)malloc(sizeof(int)*size); 72401a81e61SBarry Smith M->start_indices = (int*)malloc(sizeof(int)*size); 72501a81e61SBarry Smith M->pe = (int*)malloc(sizeof(int)*n); 72601a81e61SBarry Smith M->block_sizes = (int*)malloc(sizeof(int)*n); 72701a81e61SBarry Smith for (i=0; i<n; i++) M->block_sizes[i] = 1; 72801a81e61SBarry Smith 72901a81e61SBarry Smith ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr); 73001a81e61SBarry Smith 73101a81e61SBarry Smith M->start_indices[0] = 0; 73201a81e61SBarry Smith for (i=1; i<size; i++) { 73301a81e61SBarry Smith M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1]; 73401a81e61SBarry Smith } 73501a81e61SBarry Smith 73601a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 73701a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 73801a81e61SBarry Smith 73901a81e61SBarry Smith for (i=0; i<size; i++) { 74001a81e61SBarry Smith start_indx = M->start_indices[i]; 74101a81e61SBarry Smith for (j=0; j<M->mnls[i]; j++) 74201a81e61SBarry Smith M->pe[start_indx+j] = i; 74301a81e61SBarry Smith } 74401a81e61SBarry Smith 74501a81e61SBarry Smith if (AT) { 74601a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr); 74701a81e61SBarry Smith } else { 74801a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr); 74901a81e61SBarry Smith } 75001a81e61SBarry Smith 75101a81e61SBarry Smith rows = M->lines; 75201a81e61SBarry Smith 75301a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 75401a81e61SBarry Smith ierr = PetscMalloc(M->n*sizeof(int),&mapping);CHKERRQ(ierr); 75501a81e61SBarry Smith pe = 0; 75601a81e61SBarry Smith local_indx = 0; 75701a81e61SBarry Smith for (i=0; i<M->n; i++) { 75801a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 75901a81e61SBarry Smith pe++; 76001a81e61SBarry Smith local_indx = 0; 76101a81e61SBarry Smith } 76201a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 76301a81e61SBarry Smith local_indx++; 76401a81e61SBarry Smith } 76501a81e61SBarry Smith 76601a81e61SBarry Smith 76701a81e61SBarry Smith ierr = PetscMalloc(mnl*sizeof(int),&num_ptr);CHKERRQ(ierr); 76801a81e61SBarry Smith 76901a81e61SBarry Smith /*********************************************************/ 77001a81e61SBarry Smith /************** Set up the row structure *****************/ 77101a81e61SBarry Smith /*********************************************************/ 77201a81e61SBarry Smith 77301a81e61SBarry Smith /* count number of nonzeros in every row */ 77401a81e61SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 77501a81e61SBarry Smith for (i=rstart; i<rend; i++) { 77601a81e61SBarry Smith ierr = MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 77701a81e61SBarry Smith ierr = MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 77801a81e61SBarry Smith } 77901a81e61SBarry Smith 78001a81e61SBarry Smith /* allocate buffers */ 78101a81e61SBarry Smith len = 0; 78201a81e61SBarry Smith for (i=0; i<mnl; i++) { 78301a81e61SBarry Smith if (len < num_ptr[i]) len = num_ptr[i]; 78401a81e61SBarry Smith } 78501a81e61SBarry Smith 78601a81e61SBarry Smith for (i=rstart; i<rend; i++) { 78701a81e61SBarry Smith row_indx = i-rstart; 78801a81e61SBarry Smith len = num_ptr[row_indx]; 78901a81e61SBarry Smith rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int)); 79001a81e61SBarry Smith rows->A[row_indx] = (double*)malloc(len*sizeof(double)); 79101a81e61SBarry Smith } 79201a81e61SBarry Smith 79301a81e61SBarry Smith /* copy the matrix */ 79401a81e61SBarry Smith for (i=rstart; i<rend; i++) { 79501a81e61SBarry Smith row_indx = i - rstart; 79601a81e61SBarry Smith ierr = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 79701a81e61SBarry Smith for (j=0; j<nz; j++) { 79801a81e61SBarry Smith col = cols[j]; 79901a81e61SBarry Smith len = rows->len[row_indx]++; 80001a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 80101a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 80201a81e61SBarry Smith } 80301a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 80401a81e61SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 80501a81e61SBarry Smith } 80601a81e61SBarry Smith 80701a81e61SBarry Smith 80801a81e61SBarry Smith /************************************************************/ 80901a81e61SBarry Smith /************** Set up the column structure *****************/ 81001a81e61SBarry Smith /*********************************************************/ 81101a81e61SBarry Smith 81201a81e61SBarry Smith if (AT) { 81301a81e61SBarry Smith 81401a81e61SBarry Smith /* count number of nonzeros in every column */ 81501a81e61SBarry Smith for (i=rstart; i<rend; i++) { 81601a81e61SBarry Smith ierr = MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 81701a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 81801a81e61SBarry Smith } 81901a81e61SBarry Smith 82001a81e61SBarry Smith /* allocate buffers */ 82101a81e61SBarry Smith len = 0; 82201a81e61SBarry Smith for (i=0; i<mnl; i++) { 82301a81e61SBarry Smith if (len < num_ptr[i]) len = num_ptr[i]; 82401a81e61SBarry Smith } 82501a81e61SBarry Smith 82601a81e61SBarry Smith for (i=rstart; i<rend; i++) { 82701a81e61SBarry Smith row_indx = i-rstart; 82801a81e61SBarry Smith len = num_ptr[row_indx]; 82901a81e61SBarry Smith rows->rptrs[row_indx] = (int*)malloc(len*sizeof(int)); 83001a81e61SBarry Smith } 83101a81e61SBarry Smith 83201a81e61SBarry Smith /* copy the matrix (i.e., the structure) */ 83301a81e61SBarry Smith for (i=rstart; i<rend; i++) { 83401a81e61SBarry Smith row_indx = i - rstart; 83501a81e61SBarry Smith ierr = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 83601a81e61SBarry Smith for (j=0; j<nz; j++) { 83701a81e61SBarry Smith col = cols[j]; 83801a81e61SBarry Smith len = rows->rlen[row_indx]++; 83901a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 84001a81e61SBarry Smith } 84101a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 84201a81e61SBarry Smith } 84301a81e61SBarry Smith } 84401a81e61SBarry Smith 84501a81e61SBarry Smith ierr = PetscFree(num_ptr);CHKERRQ(ierr); 84601a81e61SBarry Smith ierr = PetscFree(mapping);CHKERRQ(ierr); 84701a81e61SBarry Smith 84801a81e61SBarry Smith order_pointers(M); 84901a81e61SBarry Smith M->maxnz = calc_maxnz(M); 85001a81e61SBarry Smith 85101a81e61SBarry Smith *B = M; 85201a81e61SBarry Smith 85301a81e61SBarry Smith PetscFunctionReturn(0); 85401a81e61SBarry Smith } 85501a81e61SBarry Smith 85601a81e61SBarry Smith /**********************************************************************/ 85701a81e61SBarry Smith 85801a81e61SBarry Smith /* 85901a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 86001a81e61SBarry Smith This assumes that the the SPAI matrix B is stored in 86101a81e61SBarry Smith COMPRESSED-ROW format. 86201a81e61SBarry Smith */ 86301a81e61SBarry Smith #undef __FUNCT__ 86401a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat" 86501a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB) 86601a81e61SBarry Smith { 86701a81e61SBarry Smith int size,rank; 86801a81e61SBarry Smith PetscErrorCode ierr; 86901a81e61SBarry Smith int m,n,M,N; 87001a81e61SBarry Smith int d_nz,o_nz; 87101a81e61SBarry Smith int *d_nnz,*o_nnz; 87201a81e61SBarry Smith int i,k,global_row,global_col,first_diag_col,last_diag_col; 87301a81e61SBarry Smith PetscScalar val; 87401a81e61SBarry Smith 87501a81e61SBarry Smith PetscFunctionBegin; 87601a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 87701a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 87801a81e61SBarry Smith 87901a81e61SBarry Smith m = n = B->mnls[rank]; 88001a81e61SBarry Smith d_nz = o_nz = 0; 88101a81e61SBarry Smith 88201a81e61SBarry Smith /* Determine preallocation for MatCreateMPIAIJ */ 88301a81e61SBarry Smith ierr = PetscMalloc(m*sizeof(int),&d_nnz);CHKERRQ(ierr); 88401a81e61SBarry Smith ierr = PetscMalloc(m*sizeof(int),&o_nnz);CHKERRQ(ierr); 88501a81e61SBarry Smith for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0; 88601a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 88701a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 88801a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 88901a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 89001a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 89101a81e61SBarry Smith if ((global_col >= first_diag_col) && (global_col <= last_diag_col)) 89201a81e61SBarry Smith d_nnz[i]++; 89301a81e61SBarry Smith else 89401a81e61SBarry Smith o_nnz[i]++; 89501a81e61SBarry Smith } 89601a81e61SBarry Smith } 89701a81e61SBarry Smith 89801a81e61SBarry Smith M = N = B->n; 89901a81e61SBarry Smith /* Here we only know how to create AIJ format */ 900f69a0ea3SMatthew Knepley ierr = MatCreate(comm,PB);CHKERRQ(ierr); 901f69a0ea3SMatthew Knepley ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr); 90201a81e61SBarry Smith ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr); 90301a81e61SBarry Smith ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr); 90401a81e61SBarry Smith ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr); 90501a81e61SBarry Smith 90601a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 90701a81e61SBarry Smith global_row = B->start_indices[rank]+i; 90801a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 90901a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 91001a81e61SBarry Smith val = B->lines->A[i][k]; 91101a81e61SBarry Smith ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr); 91201a81e61SBarry Smith } 91301a81e61SBarry Smith } 91401a81e61SBarry Smith 91501a81e61SBarry Smith ierr = PetscFree(d_nnz);CHKERRQ(ierr); 91601a81e61SBarry Smith ierr = PetscFree(o_nnz);CHKERRQ(ierr); 91701a81e61SBarry Smith 91801a81e61SBarry Smith ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91901a81e61SBarry Smith ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 92001a81e61SBarry Smith 92101a81e61SBarry Smith PetscFunctionReturn(0); 92201a81e61SBarry Smith } 92301a81e61SBarry Smith 92401a81e61SBarry Smith /**********************************************************************/ 92501a81e61SBarry Smith 92601a81e61SBarry Smith /* 92701a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 92801a81e61SBarry Smith */ 92901a81e61SBarry Smith #undef __FUNCT__ 93001a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec" 93101a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv) 93201a81e61SBarry Smith { 93301a81e61SBarry Smith PetscErrorCode ierr; 93401a81e61SBarry Smith int size,rank,m,M,i,*mnls,*start_indices,*global_indices; 93501a81e61SBarry Smith 93601a81e61SBarry Smith PetscFunctionBegin; 93701a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 93801a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 93901a81e61SBarry Smith 94001a81e61SBarry Smith m = v->mnl; 94101a81e61SBarry Smith M = v->n; 94201a81e61SBarry Smith 94301a81e61SBarry Smith 94401a81e61SBarry Smith ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr); 94501a81e61SBarry Smith 94601a81e61SBarry Smith ierr = PetscMalloc(size*sizeof(int),&mnls);CHKERRQ(ierr); 94701a81e61SBarry Smith ierr = MPI_Allgather((void*)&v->mnl,1,MPI_INT,(void*)mnls,1,MPI_INT,comm);CHKERRQ(ierr); 94801a81e61SBarry Smith 94901a81e61SBarry Smith ierr = PetscMalloc(size*sizeof(int),&start_indices);CHKERRQ(ierr); 95001a81e61SBarry Smith start_indices[0] = 0; 95101a81e61SBarry Smith for (i=1; i<size; i++) 95201a81e61SBarry Smith start_indices[i] = start_indices[i-1] +mnls[i-1]; 95301a81e61SBarry Smith 95401a81e61SBarry Smith ierr = PetscMalloc(v->mnl*sizeof(int),&global_indices);CHKERRQ(ierr); 95501a81e61SBarry Smith for (i=0; i<v->mnl; i++) 95601a81e61SBarry Smith global_indices[i] = start_indices[rank] + i; 95701a81e61SBarry Smith 95801a81e61SBarry Smith ierr = PetscFree(mnls);CHKERRQ(ierr); 95901a81e61SBarry Smith ierr = PetscFree(start_indices);CHKERRQ(ierr); 96001a81e61SBarry Smith 96101a81e61SBarry Smith ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr); 96201a81e61SBarry Smith 96301a81e61SBarry Smith ierr = PetscFree(global_indices);CHKERRQ(ierr); 96401a81e61SBarry Smith 96501a81e61SBarry Smith ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr); 96601a81e61SBarry Smith ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr); 96701a81e61SBarry Smith 96801a81e61SBarry Smith PetscFunctionReturn(0); 96901a81e61SBarry Smith } 97001a81e61SBarry Smith 97101a81e61SBarry Smith 97201a81e61SBarry Smith 97301a81e61SBarry Smith 974