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 */ 20*762360ffSBarry Smith #define PETSC_SKIP_COMPLEX /* since spai uses I which conflicts with some complex implementations */ 2101a81e61SBarry Smith 22b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 231c7d2463SJed Brown #include <../src/ksp/pc/impls/spai/petscspai.h> 2401a81e61SBarry Smith 2501a81e61SBarry Smith /* 2601a81e61SBarry Smith These are the SPAI include files 2701a81e61SBarry Smith */ 2801a81e61SBarry Smith EXTERN_C_BEGIN 29b1bda85cSSatish Balay #define SPAI_USE_MPI /* required for setting SPAI_Comm correctly in basics.h */ 30c6db04a5SJed Brown #include <spai.h> 31c6db04a5SJed Brown #include <matrix.h> 3201a81e61SBarry Smith EXTERN_C_END 3301a81e61SBarry Smith 3409573ac7SBarry Smith extern PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**); 3509573ac7SBarry Smith extern PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix*,Mat*); 3609573ac7SBarry Smith extern PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv); 3709573ac7SBarry Smith extern PetscErrorCode MM_to_PETSC(char*,char*,char*); 3801a81e61SBarry Smith 3901a81e61SBarry Smith typedef struct { 4001a81e61SBarry Smith 4101a81e61SBarry Smith matrix *B; /* matrix in SPAI format */ 4201a81e61SBarry Smith matrix *BT; /* transpose of matrix in SPAI format */ 4301a81e61SBarry Smith matrix *M; /* the approximate inverse in SPAI format */ 4401a81e61SBarry Smith 4501a81e61SBarry Smith Mat PM; /* the approximate inverse PETSc format */ 4601a81e61SBarry Smith 4701a81e61SBarry Smith double epsilon; /* tolerance */ 4801a81e61SBarry Smith int nbsteps; /* max number of "improvement" steps per line */ 4901a81e61SBarry Smith int max; /* max dimensions of is_I, q, etc. */ 5001a81e61SBarry Smith int maxnew; /* max number of new entries per step */ 5101a81e61SBarry Smith int block_size; /* constant block size */ 5201a81e61SBarry Smith int cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */ 5301a81e61SBarry Smith int verbose; /* SPAI prints timing and statistics */ 5401a81e61SBarry Smith 5501a81e61SBarry Smith int sp; /* symmetric nonzero pattern */ 5601a81e61SBarry Smith MPI_Comm comm_spai; /* communicator to be used with spai */ 5701a81e61SBarry Smith } PC_SPAI; 5801a81e61SBarry Smith 5901a81e61SBarry Smith /**********************************************************************/ 6001a81e61SBarry Smith 6101a81e61SBarry Smith #undef __FUNCT__ 6201a81e61SBarry Smith #define __FUNCT__ "PCSetUp_SPAI" 6301a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc) 6401a81e61SBarry Smith { 6501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 6601a81e61SBarry Smith PetscErrorCode ierr; 6701a81e61SBarry Smith Mat AT; 6801a81e61SBarry Smith 6901a81e61SBarry Smith PetscFunctionBegin; 7001a81e61SBarry Smith init_SPAI(); 7101a81e61SBarry Smith 7201a81e61SBarry Smith if (ispai->sp) { 7301a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr); 7401a81e61SBarry Smith } else { 7501a81e61SBarry Smith /* Use the transpose to get the column nonzero structure. */ 76fc4dec0aSBarry Smith ierr = MatTranspose(pc->pmat,MAT_INITIAL_MATRIX,&AT);CHKERRQ(ierr); 7701a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr); 7869a3b875SJed Brown ierr = MatDestroy(&AT);CHKERRQ(ierr); 7901a81e61SBarry Smith } 8001a81e61SBarry Smith 8101a81e61SBarry Smith /* Destroy the transpose */ 8201a81e61SBarry Smith /* Don't know how to do it. PETSc developers? */ 8301a81e61SBarry Smith 8401a81e61SBarry Smith /* construct SPAI preconditioner */ 8501a81e61SBarry Smith /* FILE *messages */ /* file for warning messages */ 8601a81e61SBarry Smith /* double epsilon */ /* tolerance */ 8701a81e61SBarry Smith /* int nbsteps */ /* max number of "improvement" steps per line */ 8801a81e61SBarry Smith /* int max */ /* max dimensions of is_I, q, etc. */ 8901a81e61SBarry Smith /* int maxnew */ /* max number of new entries per step */ 9001a81e61SBarry Smith /* int block_size */ /* block_size == 1 specifies scalar elments 9101a81e61SBarry Smith block_size == n specifies nxn constant-block elements 9201a81e61SBarry Smith block_size == 0 specifies variable-block elements */ 932fa5cd67SKarl Rupp /* int cache_size */ /* one of (1,2,3,4,5,6) indicting size of cache. 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 107ce94432eSBarry Smith ierr = ConvertMatrixToMat(PetscObjectComm((PetscObject)pc),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 PetscFunctionReturn(0); 11301a81e61SBarry Smith } 11401a81e61SBarry Smith 11501a81e61SBarry Smith /**********************************************************************/ 11601a81e61SBarry Smith 11701a81e61SBarry Smith #undef __FUNCT__ 11801a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI" 11901a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y) 12001a81e61SBarry Smith { 12101a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 12201a81e61SBarry Smith PetscErrorCode ierr; 12301a81e61SBarry Smith 12401a81e61SBarry Smith PetscFunctionBegin; 12501a81e61SBarry Smith /* Now using PETSc's multiply */ 12601a81e61SBarry Smith ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr); 12701a81e61SBarry Smith PetscFunctionReturn(0); 12801a81e61SBarry Smith } 12901a81e61SBarry Smith 13001a81e61SBarry Smith /**********************************************************************/ 13101a81e61SBarry Smith 13201a81e61SBarry Smith #undef __FUNCT__ 13301a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI" 13401a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc) 13501a81e61SBarry Smith { 13601a81e61SBarry Smith PetscErrorCode ierr; 13701a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 13801a81e61SBarry Smith 13901a81e61SBarry Smith PetscFunctionBegin; 14069a3b875SJed Brown ierr = MatDestroy(&ispai->PM);CHKERRQ(ierr); 14101a81e61SBarry Smith ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr); 142c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 14301a81e61SBarry Smith PetscFunctionReturn(0); 14401a81e61SBarry Smith } 14501a81e61SBarry Smith 14601a81e61SBarry Smith /**********************************************************************/ 14701a81e61SBarry Smith 14801a81e61SBarry Smith #undef __FUNCT__ 14901a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI" 15001a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer) 15101a81e61SBarry Smith { 15201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 15301a81e61SBarry Smith PetscErrorCode ierr; 154ace3abfcSBarry Smith PetscBool iascii; 15501a81e61SBarry Smith 15601a81e61SBarry Smith PetscFunctionBegin; 157251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 15801a81e61SBarry Smith if (iascii) { 15901a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SPAI preconditioner\n");CHKERRQ(ierr); 16057622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," epsilon %g\n", (double)ispai->epsilon);CHKERRQ(ierr); 16101a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," nbsteps %d\n", ispai->nbsteps);CHKERRQ(ierr); 16201a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," max %d\n", ispai->max);CHKERRQ(ierr); 16301a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxnew %d\n", ispai->maxnew);CHKERRQ(ierr); 16401a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block_size %d\n",ispai->block_size);CHKERRQ(ierr); 16501a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," cache_size %d\n",ispai->cache_size);CHKERRQ(ierr); 16601a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," verbose %d\n", ispai->verbose);CHKERRQ(ierr); 16701a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," sp %d\n", ispai->sp);CHKERRQ(ierr); 16801a81e61SBarry Smith } 16901a81e61SBarry Smith PetscFunctionReturn(0); 17001a81e61SBarry Smith } 17101a81e61SBarry Smith 17201a81e61SBarry Smith #undef __FUNCT__ 17301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI" 174f7a08781SBarry Smith static PetscErrorCode PCSPAISetEpsilon_SPAI(PC pc,double epsilon1) 17501a81e61SBarry Smith { 17601a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1775fd66863SKarl Rupp 17801a81e61SBarry Smith PetscFunctionBegin; 17901a81e61SBarry Smith ispai->epsilon = epsilon1; 18001a81e61SBarry Smith PetscFunctionReturn(0); 18101a81e61SBarry Smith } 18201a81e61SBarry Smith 18301a81e61SBarry Smith /**********************************************************************/ 18401a81e61SBarry Smith 18501a81e61SBarry Smith #undef __FUNCT__ 18601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI" 187f7a08781SBarry Smith static PetscErrorCode PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1) 18801a81e61SBarry Smith { 18901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 1905fd66863SKarl Rupp 19101a81e61SBarry Smith PetscFunctionBegin; 19201a81e61SBarry Smith ispai->nbsteps = nbsteps1; 19301a81e61SBarry Smith PetscFunctionReturn(0); 19401a81e61SBarry Smith } 19501a81e61SBarry Smith 19601a81e61SBarry Smith /**********************************************************************/ 19701a81e61SBarry Smith 19801a81e61SBarry Smith /* added 1/7/99 g.h. */ 19901a81e61SBarry Smith #undef __FUNCT__ 20001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax_SPAI" 201f7a08781SBarry Smith static PetscErrorCode PCSPAISetMax_SPAI(PC pc,int max1) 20201a81e61SBarry Smith { 20301a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2045fd66863SKarl Rupp 20501a81e61SBarry Smith PetscFunctionBegin; 20601a81e61SBarry Smith ispai->max = max1; 20701a81e61SBarry Smith PetscFunctionReturn(0); 20801a81e61SBarry Smith } 20901a81e61SBarry Smith 21001a81e61SBarry Smith /**********************************************************************/ 21101a81e61SBarry Smith 21201a81e61SBarry Smith #undef __FUNCT__ 21301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI" 214f7a08781SBarry Smith static PetscErrorCode PCSPAISetMaxNew_SPAI(PC pc,int maxnew1) 21501a81e61SBarry Smith { 21601a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2175fd66863SKarl Rupp 21801a81e61SBarry Smith PetscFunctionBegin; 21901a81e61SBarry Smith ispai->maxnew = maxnew1; 22001a81e61SBarry Smith PetscFunctionReturn(0); 22101a81e61SBarry Smith } 22201a81e61SBarry Smith 22301a81e61SBarry Smith /**********************************************************************/ 22401a81e61SBarry Smith 22501a81e61SBarry Smith #undef __FUNCT__ 22601a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize_SPAI" 227f7a08781SBarry Smith static PetscErrorCode PCSPAISetBlockSize_SPAI(PC pc,int block_size1) 22801a81e61SBarry Smith { 22901a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2305fd66863SKarl Rupp 23101a81e61SBarry Smith PetscFunctionBegin; 23201a81e61SBarry Smith ispai->block_size = block_size1; 23301a81e61SBarry Smith PetscFunctionReturn(0); 23401a81e61SBarry Smith } 23501a81e61SBarry Smith 23601a81e61SBarry Smith /**********************************************************************/ 23701a81e61SBarry Smith 23801a81e61SBarry Smith #undef __FUNCT__ 23901a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI" 240f7a08781SBarry Smith static PetscErrorCode PCSPAISetCacheSize_SPAI(PC pc,int cache_size) 24101a81e61SBarry Smith { 24201a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2435fd66863SKarl Rupp 24401a81e61SBarry Smith PetscFunctionBegin; 24501a81e61SBarry Smith ispai->cache_size = cache_size; 24601a81e61SBarry Smith PetscFunctionReturn(0); 24701a81e61SBarry Smith } 24801a81e61SBarry Smith 24901a81e61SBarry Smith /**********************************************************************/ 25001a81e61SBarry Smith 25101a81e61SBarry Smith #undef __FUNCT__ 25201a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI" 253f7a08781SBarry Smith static PetscErrorCode PCSPAISetVerbose_SPAI(PC pc,int verbose) 25401a81e61SBarry Smith { 25501a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2565fd66863SKarl Rupp 25701a81e61SBarry Smith PetscFunctionBegin; 25801a81e61SBarry Smith ispai->verbose = verbose; 25901a81e61SBarry Smith PetscFunctionReturn(0); 26001a81e61SBarry Smith } 26101a81e61SBarry Smith 26201a81e61SBarry Smith /**********************************************************************/ 26301a81e61SBarry Smith 26401a81e61SBarry Smith #undef __FUNCT__ 26501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI" 266f7a08781SBarry Smith static PetscErrorCode PCSPAISetSp_SPAI(PC pc,int sp) 26701a81e61SBarry Smith { 26801a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 2695fd66863SKarl Rupp 27001a81e61SBarry Smith PetscFunctionBegin; 27101a81e61SBarry Smith ispai->sp = sp; 27201a81e61SBarry Smith PetscFunctionReturn(0); 27301a81e61SBarry Smith } 27401a81e61SBarry Smith 27501a81e61SBarry Smith /* -------------------------------------------------------------------*/ 27601a81e61SBarry Smith 27701a81e61SBarry Smith #undef __FUNCT__ 27801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon" 27901a81e61SBarry Smith /*@ 28001a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 28101a81e61SBarry Smith 28201a81e61SBarry Smith Input Parameters: 28301a81e61SBarry Smith + pc - the preconditioner 28401a81e61SBarry Smith - eps - epsilon (default .4) 28501a81e61SBarry Smith 28601a81e61SBarry Smith Notes: Espilon must be between 0 and 1. It controls the 28701a81e61SBarry Smith quality of the approximation of M to the inverse of 28801a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 28901a81e61SBarry Smith fill, and usually better preconditioners. In many 29001a81e61SBarry Smith cases the best choice of epsilon is the one that 29101a81e61SBarry Smith divides the total solution time equally between the 29201a81e61SBarry Smith preconditioner and the solver. 29301a81e61SBarry Smith 29401a81e61SBarry Smith Level: intermediate 29501a81e61SBarry Smith 29601a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 29701a81e61SBarry Smith @*/ 2987087cfbeSBarry Smith PetscErrorCode PCSPAISetEpsilon(PC pc,double epsilon1) 29901a81e61SBarry Smith { 3004ac538c5SBarry Smith PetscErrorCode ierr; 3015fd66863SKarl Rupp 30201a81e61SBarry Smith PetscFunctionBegin; 3034ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetEpsilon_C",(PC,double),(pc,epsilon1));CHKERRQ(ierr); 30401a81e61SBarry Smith PetscFunctionReturn(0); 30501a81e61SBarry Smith } 30601a81e61SBarry Smith 30701a81e61SBarry Smith /**********************************************************************/ 30801a81e61SBarry Smith 30901a81e61SBarry Smith #undef __FUNCT__ 31001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps" 31101a81e61SBarry Smith /*@ 31201a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 31301a81e61SBarry Smith the SPAI preconditioner 31401a81e61SBarry Smith 31501a81e61SBarry Smith Input Parameters: 31601a81e61SBarry Smith + pc - the preconditioner 31701a81e61SBarry Smith - n - number of steps (default 5) 31801a81e61SBarry Smith 31901a81e61SBarry Smith Notes: SPAI constructs to approximation to every column of 32001a81e61SBarry Smith the exact inverse of A in a series of improvement 32101a81e61SBarry Smith steps. The quality of the approximation is determined 32201a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 32301a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 32401a81e61SBarry Smith uses the best approximation constructed so far. 32501a81e61SBarry Smith 32601a81e61SBarry Smith Level: intermediate 32701a81e61SBarry Smith 32801a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew() 32901a81e61SBarry Smith @*/ 3307087cfbeSBarry Smith PetscErrorCode PCSPAISetNBSteps(PC pc,int nbsteps1) 33101a81e61SBarry Smith { 3324ac538c5SBarry Smith PetscErrorCode ierr; 3335fd66863SKarl Rupp 33401a81e61SBarry Smith PetscFunctionBegin; 3354ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetNBSteps_C",(PC,int),(pc,nbsteps1));CHKERRQ(ierr); 33601a81e61SBarry Smith PetscFunctionReturn(0); 33701a81e61SBarry Smith } 33801a81e61SBarry Smith 33901a81e61SBarry Smith /**********************************************************************/ 34001a81e61SBarry Smith 34101a81e61SBarry Smith /* added 1/7/99 g.h. */ 34201a81e61SBarry Smith #undef __FUNCT__ 34301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax" 34401a81e61SBarry Smith /*@ 34501a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 34601a81e61SBarry Smith the SPAI preconditioner 34701a81e61SBarry Smith 34801a81e61SBarry Smith Input Parameters: 34901a81e61SBarry Smith + pc - the preconditioner 35001a81e61SBarry Smith - n - size (default is 5000) 35101a81e61SBarry Smith 35201a81e61SBarry Smith Level: intermediate 35301a81e61SBarry Smith 35401a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 35501a81e61SBarry Smith @*/ 3567087cfbeSBarry Smith PetscErrorCode PCSPAISetMax(PC pc,int max1) 35701a81e61SBarry Smith { 3584ac538c5SBarry Smith PetscErrorCode ierr; 3595fd66863SKarl Rupp 36001a81e61SBarry Smith PetscFunctionBegin; 3614ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMax_C",(PC,int),(pc,max1));CHKERRQ(ierr); 36201a81e61SBarry Smith PetscFunctionReturn(0); 36301a81e61SBarry Smith } 36401a81e61SBarry Smith 36501a81e61SBarry Smith /**********************************************************************/ 36601a81e61SBarry Smith 36701a81e61SBarry Smith #undef __FUNCT__ 36801a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew" 36901a81e61SBarry Smith /*@ 37001a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 37101a81e61SBarry Smith in SPAI preconditioner 37201a81e61SBarry Smith 37301a81e61SBarry Smith Input Parameters: 37401a81e61SBarry Smith + pc - the preconditioner 37501a81e61SBarry Smith - n - maximum number (default 5) 37601a81e61SBarry Smith 37701a81e61SBarry Smith Level: intermediate 37801a81e61SBarry Smith 37901a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps() 38001a81e61SBarry Smith @*/ 3817087cfbeSBarry Smith PetscErrorCode PCSPAISetMaxNew(PC pc,int maxnew1) 38201a81e61SBarry Smith { 3834ac538c5SBarry Smith PetscErrorCode ierr; 3845fd66863SKarl Rupp 38501a81e61SBarry Smith PetscFunctionBegin; 3864ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetMaxNew_C",(PC,int),(pc,maxnew1));CHKERRQ(ierr); 38701a81e61SBarry Smith PetscFunctionReturn(0); 38801a81e61SBarry Smith } 38901a81e61SBarry Smith 39001a81e61SBarry Smith /**********************************************************************/ 39101a81e61SBarry Smith 39201a81e61SBarry Smith #undef __FUNCT__ 39301a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize" 39401a81e61SBarry Smith /*@ 39501a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 39601a81e61SBarry Smith 39701a81e61SBarry Smith Input Parameters: 39801a81e61SBarry Smith + pc - the preconditioner 39901a81e61SBarry Smith - n - block size (default 1) 40001a81e61SBarry Smith 40101a81e61SBarry Smith Notes: A block 40201a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 40301a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 40401a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 40501a81e61SBarry Smith variable sized blocks, which are determined by 40601a81e61SBarry Smith searching for dense square diagonal blocks in A. 40701a81e61SBarry Smith This can be very effective for finite-element 40801a81e61SBarry Smith matrices. 40901a81e61SBarry Smith 41001a81e61SBarry Smith SPAI will convert A to block form, use a block 41101a81e61SBarry Smith version of the preconditioner algorithm, and then 41201a81e61SBarry Smith convert the result back to scalar form. 41301a81e61SBarry Smith 41401a81e61SBarry Smith In many cases the a block-size parameter other than 1 41501a81e61SBarry Smith can lead to very significant improvement in 41601a81e61SBarry Smith performance. 41701a81e61SBarry Smith 41801a81e61SBarry Smith 41901a81e61SBarry Smith Level: intermediate 42001a81e61SBarry Smith 42101a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 42201a81e61SBarry Smith @*/ 4237087cfbeSBarry Smith PetscErrorCode PCSPAISetBlockSize(PC pc,int block_size1) 42401a81e61SBarry Smith { 4254ac538c5SBarry Smith PetscErrorCode ierr; 4265fd66863SKarl Rupp 42701a81e61SBarry Smith PetscFunctionBegin; 4284ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetBlockSize_C",(PC,int),(pc,block_size1));CHKERRQ(ierr); 42901a81e61SBarry Smith PetscFunctionReturn(0); 43001a81e61SBarry Smith } 43101a81e61SBarry Smith 43201a81e61SBarry Smith /**********************************************************************/ 43301a81e61SBarry Smith 43401a81e61SBarry Smith #undef __FUNCT__ 43501a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize" 43601a81e61SBarry Smith /*@ 43701a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 43801a81e61SBarry Smith 43901a81e61SBarry Smith Input Parameters: 44001a81e61SBarry Smith + pc - the preconditioner 44101a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 44201a81e61SBarry Smith 44301a81e61SBarry Smith Notes: SPAI uses a hash table to cache messages and avoid 44401a81e61SBarry Smith redundant communication. If suggest always using 44501a81e61SBarry Smith 5. This parameter is irrelevant in the serial 44601a81e61SBarry Smith version. 44701a81e61SBarry Smith 44801a81e61SBarry Smith Level: intermediate 44901a81e61SBarry Smith 45001a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 45101a81e61SBarry Smith @*/ 4527087cfbeSBarry Smith PetscErrorCode PCSPAISetCacheSize(PC pc,int cache_size) 45301a81e61SBarry Smith { 4544ac538c5SBarry Smith PetscErrorCode ierr; 4555fd66863SKarl Rupp 45601a81e61SBarry Smith PetscFunctionBegin; 4574ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetCacheSize_C",(PC,int),(pc,cache_size));CHKERRQ(ierr); 45801a81e61SBarry Smith PetscFunctionReturn(0); 45901a81e61SBarry Smith } 46001a81e61SBarry Smith 46101a81e61SBarry Smith /**********************************************************************/ 46201a81e61SBarry Smith 46301a81e61SBarry Smith #undef __FUNCT__ 46401a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose" 46501a81e61SBarry Smith /*@ 46601a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 46701a81e61SBarry Smith 46801a81e61SBarry Smith Input Parameters: 46901a81e61SBarry Smith + pc - the preconditioner 47001a81e61SBarry Smith - n - level (default 1) 47101a81e61SBarry Smith 47201a81e61SBarry Smith Notes: print parameters, timings and matrix statistics 47301a81e61SBarry Smith 47401a81e61SBarry Smith Level: intermediate 47501a81e61SBarry Smith 47601a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 47701a81e61SBarry Smith @*/ 4787087cfbeSBarry Smith PetscErrorCode PCSPAISetVerbose(PC pc,int verbose) 47901a81e61SBarry Smith { 4804ac538c5SBarry Smith PetscErrorCode ierr; 4815fd66863SKarl Rupp 48201a81e61SBarry Smith PetscFunctionBegin; 4834ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetVerbose_C",(PC,int),(pc,verbose));CHKERRQ(ierr); 48401a81e61SBarry Smith PetscFunctionReturn(0); 48501a81e61SBarry Smith } 48601a81e61SBarry Smith 48701a81e61SBarry Smith /**********************************************************************/ 48801a81e61SBarry Smith 48901a81e61SBarry Smith #undef __FUNCT__ 49001a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp" 49101a81e61SBarry Smith /*@ 49201a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 49301a81e61SBarry Smith 49401a81e61SBarry Smith Input Parameters: 49501a81e61SBarry Smith + pc - the preconditioner 49601a81e61SBarry Smith - n - 0 or 1 49701a81e61SBarry Smith 49801a81e61SBarry Smith Notes: If A has a symmetric nonzero pattern use -sp 1 to 49901a81e61SBarry Smith improve performance by eliminating some communication 50001a81e61SBarry Smith in the parallel version. Even if A does not have a 50101a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 50201a81e61SBarry Smith results, but the code will not follow the published 50301a81e61SBarry Smith SPAI algorithm exactly. 50401a81e61SBarry Smith 50501a81e61SBarry Smith 50601a81e61SBarry Smith Level: intermediate 50701a81e61SBarry Smith 50801a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 50901a81e61SBarry Smith @*/ 5107087cfbeSBarry Smith PetscErrorCode PCSPAISetSp(PC pc,int sp) 51101a81e61SBarry Smith { 5124ac538c5SBarry Smith PetscErrorCode ierr; 5135fd66863SKarl Rupp 51401a81e61SBarry Smith PetscFunctionBegin; 5154ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSPAISetSp_C",(PC,int),(pc,sp));CHKERRQ(ierr); 51601a81e61SBarry Smith PetscFunctionReturn(0); 51701a81e61SBarry Smith } 51801a81e61SBarry Smith 51901a81e61SBarry Smith /**********************************************************************/ 52001a81e61SBarry Smith 52101a81e61SBarry Smith /**********************************************************************/ 52201a81e61SBarry Smith 52301a81e61SBarry Smith #undef __FUNCT__ 52401a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI" 52501a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc) 52601a81e61SBarry Smith { 52701a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 52801a81e61SBarry Smith PetscErrorCode ierr; 52901a81e61SBarry Smith int nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp; 53001a81e61SBarry Smith double epsilon1; 531ace3abfcSBarry Smith PetscBool flg; 53201a81e61SBarry Smith 53301a81e61SBarry Smith PetscFunctionBegin; 53401a81e61SBarry Smith ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr); 53501a81e61SBarry Smith ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr); 53601a81e61SBarry Smith if (flg) { 53701a81e61SBarry Smith ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr); 53801a81e61SBarry Smith } 53901a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr); 54001a81e61SBarry Smith if (flg) { 54101a81e61SBarry Smith ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr); 54201a81e61SBarry Smith } 54301a81e61SBarry Smith /* added 1/7/99 g.h. */ 54401a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr); 54501a81e61SBarry Smith if (flg) { 54601a81e61SBarry Smith ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr); 54701a81e61SBarry Smith } 54801a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr); 54901a81e61SBarry Smith if (flg) { 55001a81e61SBarry Smith ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr); 55101a81e61SBarry Smith } 55201a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr); 55301a81e61SBarry Smith if (flg) { 55401a81e61SBarry Smith ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr); 55501a81e61SBarry Smith } 55601a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr); 55701a81e61SBarry Smith if (flg) { 55801a81e61SBarry Smith ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr); 55901a81e61SBarry Smith } 56001a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr); 56101a81e61SBarry Smith if (flg) { 56201a81e61SBarry Smith ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr); 56301a81e61SBarry Smith } 56401a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr); 56501a81e61SBarry Smith if (flg) { 56601a81e61SBarry Smith ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr); 56701a81e61SBarry Smith } 56801a81e61SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 56901a81e61SBarry Smith PetscFunctionReturn(0); 57001a81e61SBarry Smith } 57101a81e61SBarry Smith 57201a81e61SBarry Smith /**********************************************************************/ 57301a81e61SBarry Smith 57401a81e61SBarry Smith /*MC 57501a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 57601a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 57701a81e61SBarry Smith 57801a81e61SBarry Smith Options Database Keys: 579c5ae99e3SSatish Balay + -pc_spai_epsilon <eps> - set tolerance 58001a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 58101a81e61SBarry Smith . -pc_spai_max <m> - set max 58201a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 58301a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 58401a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 58501a81e61SBarry Smith . -pc_spai_sp <m> - set sp 58601a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 58701a81e61SBarry Smith 58801a81e61SBarry Smith Notes: This only works with AIJ matrices. 58901a81e61SBarry Smith 59001a81e61SBarry Smith Level: beginner 59101a81e61SBarry Smith 59201a81e61SBarry Smith Concepts: approximate inverse 59301a81e61SBarry Smith 59401a81e61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 59501a81e61SBarry Smith PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(), 59601a81e61SBarry Smith PCSPAISetVerbose(), PCSPAISetSp() 59701a81e61SBarry Smith M*/ 59801a81e61SBarry Smith 59901a81e61SBarry Smith #undef __FUNCT__ 60001a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI" 6018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SPAI(PC pc) 60201a81e61SBarry Smith { 60301a81e61SBarry Smith PC_SPAI *ispai; 60401a81e61SBarry Smith PetscErrorCode ierr; 60501a81e61SBarry Smith 60601a81e61SBarry Smith PetscFunctionBegin; 607b00a9115SJed Brown ierr = PetscNewLog(pc,&ispai);CHKERRQ(ierr); 6082a4967abSBarry Smith pc->data = ispai; 60901a81e61SBarry Smith 61001a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 61101a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 61201a81e61SBarry Smith pc->ops->applyrichardson = 0; 61301a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 61401a81e61SBarry Smith pc->ops->view = PCView_SPAI; 61501a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 61601a81e61SBarry Smith 61701a81e61SBarry Smith ispai->epsilon = .4; 61801a81e61SBarry Smith ispai->nbsteps = 5; 61901a81e61SBarry Smith ispai->max = 5000; 62001a81e61SBarry Smith ispai->maxnew = 5; 62101a81e61SBarry Smith ispai->block_size = 1; 62201a81e61SBarry Smith ispai->cache_size = 5; 62301a81e61SBarry Smith ispai->verbose = 0; 62401a81e61SBarry Smith 62501a81e61SBarry Smith ispai->sp = 1; 626ce94432eSBarry Smith ierr = MPI_Comm_dup(PetscObjectComm((PetscObject)pc),&(ispai->comm_spai));CHKERRQ(ierr); 62701a81e61SBarry Smith 628bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetEpsilon_C",PCSPAISetEpsilon_SPAI);CHKERRQ(ierr); 629bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetNBSteps_C",PCSPAISetNBSteps_SPAI);CHKERRQ(ierr); 630bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMax_C",PCSPAISetMax_SPAI);CHKERRQ(ierr); 631bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetMaxNew_C",PCSPAISetMaxNew_SPAI);CHKERRQ(ierr); 632bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetBlockSize_C",PCSPAISetBlockSize_SPAI);CHKERRQ(ierr); 633bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetCacheSize_C",PCSPAISetCacheSize_SPAI);CHKERRQ(ierr); 634bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetVerbose_C",PCSPAISetVerbose_SPAI);CHKERRQ(ierr); 635bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSPAISetSp_C",PCSPAISetSp_SPAI);CHKERRQ(ierr); 63601a81e61SBarry Smith PetscFunctionReturn(0); 63701a81e61SBarry Smith } 63801a81e61SBarry Smith 63901a81e61SBarry Smith /**********************************************************************/ 64001a81e61SBarry Smith 64101a81e61SBarry Smith /* 64201a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 64301a81e61SBarry Smith */ 64401a81e61SBarry Smith #undef __FUNCT__ 64501a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix" 64601a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B) 64701a81e61SBarry Smith { 64801a81e61SBarry Smith matrix *M; 64901a81e61SBarry Smith int i,j,col; 65001a81e61SBarry Smith int row_indx; 65101a81e61SBarry Smith int len,pe,local_indx,start_indx; 65201a81e61SBarry Smith int *mapping; 65301a81e61SBarry Smith PetscErrorCode ierr; 65401a81e61SBarry Smith const int *cols; 65501a81e61SBarry Smith const double *vals; 6562122902bSSatish Balay int n,mnl,nnl,nz,rstart,rend; 6572a4967abSBarry Smith PetscMPIInt size,rank; 65801a81e61SBarry Smith struct compressed_lines *rows; 65901a81e61SBarry Smith 66001a81e61SBarry Smith PetscFunctionBegin; 66101a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 66201a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 66301a81e61SBarry Smith ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr); 66401a81e61SBarry Smith ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr); 66501a81e61SBarry Smith 66601a81e61SBarry Smith /* 66701a81e61SBarry Smith not sure why a barrier is required. commenting out 66801a81e61SBarry Smith ierr = MPI_Barrier(comm);CHKERRQ(ierr); 66901a81e61SBarry Smith */ 67001a81e61SBarry Smith 67129b2f704SSatish Balay M = new_matrix((SPAI_Comm)comm); 67201a81e61SBarry Smith 67301a81e61SBarry Smith M->n = n; 67401a81e61SBarry Smith M->bs = 1; 67501a81e61SBarry Smith M->max_block_size = 1; 67601a81e61SBarry Smith 67701a81e61SBarry Smith M->mnls = (int*)malloc(sizeof(int)*size); 67801a81e61SBarry Smith M->start_indices = (int*)malloc(sizeof(int)*size); 67901a81e61SBarry Smith M->pe = (int*)malloc(sizeof(int)*n); 68001a81e61SBarry Smith M->block_sizes = (int*)malloc(sizeof(int)*n); 68101a81e61SBarry Smith for (i=0; i<n; i++) M->block_sizes[i] = 1; 68201a81e61SBarry Smith 68301a81e61SBarry Smith ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr); 68401a81e61SBarry Smith 68501a81e61SBarry Smith M->start_indices[0] = 0; 6862fa5cd67SKarl Rupp for (i=1; i<size; i++) M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1]; 68701a81e61SBarry Smith 68801a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 68901a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 69001a81e61SBarry Smith 69101a81e61SBarry Smith for (i=0; i<size; i++) { 69201a81e61SBarry Smith start_indx = M->start_indices[i]; 6932fa5cd67SKarl Rupp for (j=0; j<M->mnls[i]; j++) M->pe[start_indx+j] = i; 69401a81e61SBarry Smith } 69501a81e61SBarry Smith 69601a81e61SBarry Smith if (AT) { 69701a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr); 69801a81e61SBarry Smith } else { 69901a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr); 70001a81e61SBarry Smith } 70101a81e61SBarry Smith 70201a81e61SBarry Smith rows = M->lines; 70301a81e61SBarry Smith 70401a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 705785e854fSJed Brown ierr = PetscMalloc1(M->n,&mapping);CHKERRQ(ierr); 70601a81e61SBarry Smith pe = 0; 70701a81e61SBarry Smith local_indx = 0; 70801a81e61SBarry Smith for (i=0; i<M->n; i++) { 70901a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 71001a81e61SBarry Smith pe++; 71101a81e61SBarry Smith local_indx = 0; 71201a81e61SBarry Smith } 71301a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 71401a81e61SBarry Smith local_indx++; 71501a81e61SBarry Smith } 71601a81e61SBarry Smith 71701a81e61SBarry Smith /*********************************************************/ 71801a81e61SBarry Smith /************** Set up the row structure *****************/ 71901a81e61SBarry Smith /*********************************************************/ 72001a81e61SBarry Smith 72101a81e61SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 72201a81e61SBarry Smith for (i=rstart; i<rend; i++) { 72301a81e61SBarry Smith row_indx = i - rstart; 72401a81e61SBarry Smith ierr = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 7252122902bSSatish Balay /* allocate buffers */ 7262122902bSSatish Balay rows->ptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 7272122902bSSatish Balay rows->A[row_indx] = (double*)malloc(nz*sizeof(double)); 7282122902bSSatish Balay /* copy the matrix */ 72901a81e61SBarry Smith for (j=0; j<nz; j++) { 73001a81e61SBarry Smith col = cols[j]; 73101a81e61SBarry Smith len = rows->len[row_indx]++; 7322fa5cd67SKarl Rupp 73301a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 73401a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 73501a81e61SBarry Smith } 73601a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 7372fa5cd67SKarl Rupp 73801a81e61SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 73901a81e61SBarry Smith } 74001a81e61SBarry Smith 74101a81e61SBarry Smith 74201a81e61SBarry Smith /************************************************************/ 74301a81e61SBarry Smith /************** Set up the column structure *****************/ 74401a81e61SBarry Smith /*********************************************************/ 74501a81e61SBarry Smith 74601a81e61SBarry Smith if (AT) { 74701a81e61SBarry Smith 74801a81e61SBarry Smith for (i=rstart; i<rend; i++) { 74901a81e61SBarry Smith row_indx = i - rstart; 75001a81e61SBarry Smith ierr = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 7512122902bSSatish Balay /* allocate buffers */ 7522122902bSSatish Balay rows->rptrs[row_indx] = (int*)malloc(nz*sizeof(int)); 7532122902bSSatish Balay /* copy the matrix (i.e., the structure) */ 75401a81e61SBarry Smith for (j=0; j<nz; j++) { 75501a81e61SBarry Smith col = cols[j]; 75601a81e61SBarry Smith len = rows->rlen[row_indx]++; 7572fa5cd67SKarl Rupp 75801a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 75901a81e61SBarry Smith } 76001a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 76101a81e61SBarry Smith } 76201a81e61SBarry Smith } 76301a81e61SBarry Smith 76401a81e61SBarry Smith ierr = PetscFree(mapping);CHKERRQ(ierr); 76501a81e61SBarry Smith 76601a81e61SBarry Smith order_pointers(M); 76701a81e61SBarry Smith M->maxnz = calc_maxnz(M); 76801a81e61SBarry Smith *B = M; 76901a81e61SBarry Smith PetscFunctionReturn(0); 77001a81e61SBarry Smith } 77101a81e61SBarry Smith 77201a81e61SBarry Smith /**********************************************************************/ 77301a81e61SBarry Smith 77401a81e61SBarry Smith /* 77501a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 776df3898eeSBarry Smith This assumes that the SPAI matrix B is stored in 77701a81e61SBarry Smith COMPRESSED-ROW format. 77801a81e61SBarry Smith */ 77901a81e61SBarry Smith #undef __FUNCT__ 78001a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat" 78101a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB) 78201a81e61SBarry Smith { 7834b6b5fe1SBarry Smith PetscMPIInt size,rank; 78401a81e61SBarry Smith PetscErrorCode ierr; 78501a81e61SBarry Smith int m,n,M,N; 78601a81e61SBarry Smith int d_nz,o_nz; 78701a81e61SBarry Smith int *d_nnz,*o_nnz; 78801a81e61SBarry Smith int i,k,global_row,global_col,first_diag_col,last_diag_col; 78901a81e61SBarry Smith PetscScalar val; 79001a81e61SBarry Smith 79101a81e61SBarry Smith PetscFunctionBegin; 79201a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 79301a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 79401a81e61SBarry Smith 79501a81e61SBarry Smith m = n = B->mnls[rank]; 79601a81e61SBarry Smith d_nz = o_nz = 0; 79701a81e61SBarry Smith 79801a81e61SBarry Smith /* Determine preallocation for MatCreateMPIAIJ */ 799785e854fSJed Brown ierr = PetscMalloc1(m,&d_nnz);CHKERRQ(ierr); 800785e854fSJed Brown ierr = PetscMalloc1(m,&o_nnz);CHKERRQ(ierr); 80101a81e61SBarry Smith for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0; 80201a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 80301a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 80401a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 80501a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 80601a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 807db4deed7SKarl Rupp if ((global_col >= first_diag_col) && (global_col < last_diag_col)) d_nnz[i]++; 808db4deed7SKarl Rupp else o_nnz[i]++; 80901a81e61SBarry Smith } 81001a81e61SBarry Smith } 81101a81e61SBarry Smith 81201a81e61SBarry Smith M = N = B->n; 81301a81e61SBarry Smith /* Here we only know how to create AIJ format */ 814f69a0ea3SMatthew Knepley ierr = MatCreate(comm,PB);CHKERRQ(ierr); 815f69a0ea3SMatthew Knepley ierr = MatSetSizes(*PB,m,n,M,N);CHKERRQ(ierr); 81601a81e61SBarry Smith ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr); 81701a81e61SBarry Smith ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr); 81801a81e61SBarry Smith ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr); 81901a81e61SBarry Smith 82001a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 82101a81e61SBarry Smith global_row = B->start_indices[rank]+i; 82201a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 82301a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 8242fa5cd67SKarl Rupp 82501a81e61SBarry Smith val = B->lines->A[i][k]; 82601a81e61SBarry Smith ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr); 82701a81e61SBarry Smith } 82801a81e61SBarry Smith } 82901a81e61SBarry Smith 83001a81e61SBarry Smith ierr = PetscFree(d_nnz);CHKERRQ(ierr); 83101a81e61SBarry Smith ierr = PetscFree(o_nnz);CHKERRQ(ierr); 83201a81e61SBarry Smith 83301a81e61SBarry Smith ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83401a81e61SBarry Smith ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 83501a81e61SBarry Smith PetscFunctionReturn(0); 83601a81e61SBarry Smith } 83701a81e61SBarry Smith 83801a81e61SBarry Smith /**********************************************************************/ 83901a81e61SBarry Smith 84001a81e61SBarry Smith /* 84101a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 84201a81e61SBarry Smith */ 84301a81e61SBarry Smith #undef __FUNCT__ 84401a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec" 84501a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv) 84601a81e61SBarry Smith { 84701a81e61SBarry Smith PetscErrorCode ierr; 8484b6b5fe1SBarry Smith PetscMPIInt size,rank; 8494b6b5fe1SBarry Smith int m,M,i,*mnls,*start_indices,*global_indices; 85001a81e61SBarry Smith 85101a81e61SBarry Smith PetscFunctionBegin; 85201a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 85301a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 85401a81e61SBarry Smith 85501a81e61SBarry Smith m = v->mnl; 85601a81e61SBarry Smith M = v->n; 85701a81e61SBarry Smith 85801a81e61SBarry Smith 85901a81e61SBarry Smith ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr); 86001a81e61SBarry Smith 861785e854fSJed Brown ierr = PetscMalloc1(size,&mnls);CHKERRQ(ierr); 8622a4967abSBarry Smith ierr = MPI_Allgather(&v->mnl,1,MPI_INT,mnls,1,MPI_INT,comm);CHKERRQ(ierr); 86301a81e61SBarry Smith 864785e854fSJed Brown ierr = PetscMalloc1(size,&start_indices);CHKERRQ(ierr); 8652fa5cd67SKarl Rupp 86601a81e61SBarry Smith start_indices[0] = 0; 8672fa5cd67SKarl Rupp for (i=1; i<size; i++) start_indices[i] = start_indices[i-1] +mnls[i-1]; 86801a81e61SBarry Smith 869785e854fSJed Brown ierr = PetscMalloc1(v->mnl,&global_indices);CHKERRQ(ierr); 8702fa5cd67SKarl Rupp for (i=0; i<v->mnl; i++) global_indices[i] = start_indices[rank] + i; 87101a81e61SBarry Smith 87201a81e61SBarry Smith ierr = PetscFree(mnls);CHKERRQ(ierr); 87301a81e61SBarry Smith ierr = PetscFree(start_indices);CHKERRQ(ierr); 87401a81e61SBarry Smith 87501a81e61SBarry Smith ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr); 87601a81e61SBarry Smith ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr); 87701a81e61SBarry Smith ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr); 87801a81e61SBarry Smith 8794b6b5fe1SBarry Smith ierr = PetscFree(global_indices);CHKERRQ(ierr); 88001a81e61SBarry Smith PetscFunctionReturn(0); 88101a81e61SBarry Smith } 88201a81e61SBarry Smith 88301a81e61SBarry Smith 88401a81e61SBarry Smith 88501a81e61SBarry Smith 886