1*01a81e61SBarry Smith #define PETSCKSP_DLL 2*01a81e61SBarry Smith 3*01a81e61SBarry Smith /* 4*01a81e61SBarry Smith 3/99 Modified by Stephen Barnard to support SPAI version 3.0 5*01a81e61SBarry Smith */ 6*01a81e61SBarry Smith 7*01a81e61SBarry Smith /* 8*01a81e61SBarry Smith Provides an interface to the SPAI Sparse Approximate Inverse Preconditioner 9*01a81e61SBarry Smith Code written by Stephen Barnard. 10*01a81e61SBarry Smith 11*01a81e61SBarry Smith Note: there is some BAD memory bleeding below! 12*01a81e61SBarry Smith 13*01a81e61SBarry Smith This code needs work 14*01a81e61SBarry Smith 15*01a81e61SBarry Smith 1) get rid of all memory bleeding 16*01a81e61SBarry Smith 2) fix PETSc/interface so that it gets if the matrix is symmetric from the matrix 17*01a81e61SBarry Smith rather than having the sp flag for PC_SPAI 18*01a81e61SBarry Smith 19*01a81e61SBarry Smith */ 20*01a81e61SBarry Smith 21*01a81e61SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 22*01a81e61SBarry Smith 23*01a81e61SBarry Smith /* 24*01a81e61SBarry Smith These are the SPAI include files 25*01a81e61SBarry Smith */ 26*01a81e61SBarry Smith EXTERN_C_BEGIN 27*01a81e61SBarry Smith #include "spai.h" 28*01a81e61SBarry Smith #include "matrix.h" 29*01a81e61SBarry Smith EXTERN_C_END 30*01a81e61SBarry Smith 31*01a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatToMatrix(MPI_Comm,Mat,Mat,matrix**); 32*01a81e61SBarry Smith EXTERN PetscErrorCode ConvertMatrixToMat(MPI_Comm,matrix *,Mat *); 33*01a81e61SBarry Smith EXTERN PetscErrorCode ConvertVectorToVec(MPI_Comm,vector *v,Vec *Pv); 34*01a81e61SBarry Smith EXTERN PetscErrorCode MM_to_PETSC(char *,char *,char *); 35*01a81e61SBarry Smith 36*01a81e61SBarry Smith typedef struct { 37*01a81e61SBarry Smith 38*01a81e61SBarry Smith matrix *B; /* matrix in SPAI format */ 39*01a81e61SBarry Smith matrix *BT; /* transpose of matrix in SPAI format */ 40*01a81e61SBarry Smith matrix *M; /* the approximate inverse in SPAI format */ 41*01a81e61SBarry Smith 42*01a81e61SBarry Smith Mat PM; /* the approximate inverse PETSc format */ 43*01a81e61SBarry Smith 44*01a81e61SBarry Smith double epsilon; /* tolerance */ 45*01a81e61SBarry Smith int nbsteps; /* max number of "improvement" steps per line */ 46*01a81e61SBarry Smith int max; /* max dimensions of is_I, q, etc. */ 47*01a81e61SBarry Smith int maxnew; /* max number of new entries per step */ 48*01a81e61SBarry Smith int block_size; /* constant block size */ 49*01a81e61SBarry Smith int cache_size; /* one of (1,2,3,4,5,6) indicting size of cache */ 50*01a81e61SBarry Smith int verbose; /* SPAI prints timing and statistics */ 51*01a81e61SBarry Smith 52*01a81e61SBarry Smith int sp; /* symmetric nonzero pattern */ 53*01a81e61SBarry Smith MPI_Comm comm_spai; /* communicator to be used with spai */ 54*01a81e61SBarry Smith } PC_SPAI; 55*01a81e61SBarry Smith 56*01a81e61SBarry Smith /**********************************************************************/ 57*01a81e61SBarry Smith 58*01a81e61SBarry Smith #undef __FUNCT__ 59*01a81e61SBarry Smith #define __FUNCT__ "PCSetUp_SPAI" 60*01a81e61SBarry Smith static PetscErrorCode PCSetUp_SPAI(PC pc) 61*01a81e61SBarry Smith { 62*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 63*01a81e61SBarry Smith PetscErrorCode ierr; 64*01a81e61SBarry Smith Mat AT; 65*01a81e61SBarry Smith 66*01a81e61SBarry Smith PetscFunctionBegin; 67*01a81e61SBarry Smith 68*01a81e61SBarry Smith init_SPAI(); 69*01a81e61SBarry Smith 70*01a81e61SBarry Smith if (ispai->sp) { 71*01a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,pc->pmat,&ispai->B);CHKERRQ(ierr); 72*01a81e61SBarry Smith } else { 73*01a81e61SBarry Smith /* Use the transpose to get the column nonzero structure. */ 74*01a81e61SBarry Smith ierr = MatTranspose(pc->pmat,&AT);CHKERRQ(ierr); 75*01a81e61SBarry Smith ierr = ConvertMatToMatrix(ispai->comm_spai,pc->pmat,AT,&ispai->B);CHKERRQ(ierr); 76*01a81e61SBarry Smith ierr = MatDestroy(AT);CHKERRQ(ierr); 77*01a81e61SBarry Smith } 78*01a81e61SBarry Smith 79*01a81e61SBarry Smith /* Destroy the transpose */ 80*01a81e61SBarry Smith /* Don't know how to do it. PETSc developers? */ 81*01a81e61SBarry Smith 82*01a81e61SBarry Smith /* construct SPAI preconditioner */ 83*01a81e61SBarry Smith /* FILE *messages */ /* file for warning messages */ 84*01a81e61SBarry Smith /* double epsilon */ /* tolerance */ 85*01a81e61SBarry Smith /* int nbsteps */ /* max number of "improvement" steps per line */ 86*01a81e61SBarry Smith /* int max */ /* max dimensions of is_I, q, etc. */ 87*01a81e61SBarry Smith /* int maxnew */ /* max number of new entries per step */ 88*01a81e61SBarry Smith /* int block_size */ /* block_size == 1 specifies scalar elments 89*01a81e61SBarry Smith block_size == n specifies nxn constant-block elements 90*01a81e61SBarry Smith block_size == 0 specifies variable-block elements */ 91*01a81e61SBarry Smith /* int cache_size */ /* one of (1,2,3,4,5,6) indicting size of cache */ 92*01a81e61SBarry Smith /* cache_size == 0 indicates no caching */ 93*01a81e61SBarry Smith /* int verbose */ /* verbose == 0 specifies that SPAI is silent 94*01a81e61SBarry Smith verbose == 1 prints timing and matrix statistics */ 95*01a81e61SBarry Smith 96*01a81e61SBarry Smith ierr = bspai(ispai->B,&ispai->M, 97*01a81e61SBarry Smith stdout, 98*01a81e61SBarry Smith ispai->epsilon, 99*01a81e61SBarry Smith ispai->nbsteps, 100*01a81e61SBarry Smith ispai->max, 101*01a81e61SBarry Smith ispai->maxnew, 102*01a81e61SBarry Smith ispai->block_size, 103*01a81e61SBarry Smith ispai->cache_size, 104*01a81e61SBarry Smith ispai->verbose); CHKERRQ(ierr); 105*01a81e61SBarry Smith 106*01a81e61SBarry Smith ierr = ConvertMatrixToMat(pc->comm,ispai->M,&ispai->PM);CHKERRQ(ierr); 107*01a81e61SBarry Smith 108*01a81e61SBarry Smith /* free the SPAI matrices */ 109*01a81e61SBarry Smith sp_free_matrix(ispai->B); 110*01a81e61SBarry Smith sp_free_matrix(ispai->M); 111*01a81e61SBarry Smith 112*01a81e61SBarry Smith PetscFunctionReturn(0); 113*01a81e61SBarry Smith } 114*01a81e61SBarry Smith 115*01a81e61SBarry Smith /**********************************************************************/ 116*01a81e61SBarry Smith 117*01a81e61SBarry Smith #undef __FUNCT__ 118*01a81e61SBarry Smith #define __FUNCT__ "PCApply_SPAI" 119*01a81e61SBarry Smith static PetscErrorCode PCApply_SPAI(PC pc,Vec xx,Vec y) 120*01a81e61SBarry Smith { 121*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 122*01a81e61SBarry Smith PetscErrorCode ierr; 123*01a81e61SBarry Smith 124*01a81e61SBarry Smith PetscFunctionBegin; 125*01a81e61SBarry Smith /* Now using PETSc's multiply */ 126*01a81e61SBarry Smith ierr = MatMult(ispai->PM,xx,y);CHKERRQ(ierr); 127*01a81e61SBarry Smith PetscFunctionReturn(0); 128*01a81e61SBarry Smith } 129*01a81e61SBarry Smith 130*01a81e61SBarry Smith /**********************************************************************/ 131*01a81e61SBarry Smith 132*01a81e61SBarry Smith #undef __FUNCT__ 133*01a81e61SBarry Smith #define __FUNCT__ "PCDestroy_SPAI" 134*01a81e61SBarry Smith static PetscErrorCode PCDestroy_SPAI(PC pc) 135*01a81e61SBarry Smith { 136*01a81e61SBarry Smith PetscErrorCode ierr; 137*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 138*01a81e61SBarry Smith 139*01a81e61SBarry Smith PetscFunctionBegin; 140*01a81e61SBarry Smith if (ispai->PM) {ierr = MatDestroy(ispai->PM);CHKERRQ(ierr);} 141*01a81e61SBarry Smith ierr = MPI_Comm_free(&(ispai->comm_spai));CHKERRQ(ierr); 142*01a81e61SBarry Smith ierr = PetscFree(ispai);CHKERRQ(ierr); 143*01a81e61SBarry Smith PetscFunctionReturn(0); 144*01a81e61SBarry Smith } 145*01a81e61SBarry Smith 146*01a81e61SBarry Smith /**********************************************************************/ 147*01a81e61SBarry Smith 148*01a81e61SBarry Smith #undef __FUNCT__ 149*01a81e61SBarry Smith #define __FUNCT__ "PCView_SPAI" 150*01a81e61SBarry Smith static PetscErrorCode PCView_SPAI(PC pc,PetscViewer viewer) 151*01a81e61SBarry Smith { 152*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 153*01a81e61SBarry Smith PetscErrorCode ierr; 154*01a81e61SBarry Smith PetscTruth iascii; 155*01a81e61SBarry Smith 156*01a81e61SBarry Smith PetscFunctionBegin; 157*01a81e61SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 158*01a81e61SBarry Smith if (iascii) { 159*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SPAI preconditioner\n");CHKERRQ(ierr); 160*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," epsilon %g\n", ispai->epsilon);CHKERRQ(ierr); 161*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," nbsteps %d\n", ispai->nbsteps);CHKERRQ(ierr); 162*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," max %d\n", ispai->max);CHKERRQ(ierr); 163*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maxnew %d\n", ispai->maxnew);CHKERRQ(ierr); 164*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," block_size %d\n",ispai->block_size);CHKERRQ(ierr); 165*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," cache_size %d\n",ispai->cache_size);CHKERRQ(ierr); 166*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," verbose %d\n", ispai->verbose);CHKERRQ(ierr); 167*01a81e61SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," sp %d\n", ispai->sp);CHKERRQ(ierr); 168*01a81e61SBarry Smith 169*01a81e61SBarry Smith } 170*01a81e61SBarry Smith PetscFunctionReturn(0); 171*01a81e61SBarry Smith } 172*01a81e61SBarry Smith 173*01a81e61SBarry Smith EXTERN_C_BEGIN 174*01a81e61SBarry Smith #undef __FUNCT__ 175*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon_SPAI" 176*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon_SPAI(PC pc,double epsilon1) 177*01a81e61SBarry Smith { 178*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 179*01a81e61SBarry Smith PetscFunctionBegin; 180*01a81e61SBarry Smith ispai->epsilon = epsilon1; 181*01a81e61SBarry Smith PetscFunctionReturn(0); 182*01a81e61SBarry Smith } 183*01a81e61SBarry Smith EXTERN_C_END 184*01a81e61SBarry Smith 185*01a81e61SBarry Smith /**********************************************************************/ 186*01a81e61SBarry Smith 187*01a81e61SBarry Smith EXTERN_C_BEGIN 188*01a81e61SBarry Smith #undef __FUNCT__ 189*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps_SPAI" 190*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps_SPAI(PC pc,int nbsteps1) 191*01a81e61SBarry Smith { 192*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 193*01a81e61SBarry Smith PetscFunctionBegin; 194*01a81e61SBarry Smith ispai->nbsteps = nbsteps1; 195*01a81e61SBarry Smith PetscFunctionReturn(0); 196*01a81e61SBarry Smith } 197*01a81e61SBarry Smith EXTERN_C_END 198*01a81e61SBarry Smith 199*01a81e61SBarry Smith /**********************************************************************/ 200*01a81e61SBarry Smith 201*01a81e61SBarry Smith /* added 1/7/99 g.h. */ 202*01a81e61SBarry Smith EXTERN_C_BEGIN 203*01a81e61SBarry Smith #undef __FUNCT__ 204*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax_SPAI" 205*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax_SPAI(PC pc,int max1) 206*01a81e61SBarry Smith { 207*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 208*01a81e61SBarry Smith PetscFunctionBegin; 209*01a81e61SBarry Smith ispai->max = max1; 210*01a81e61SBarry Smith PetscFunctionReturn(0); 211*01a81e61SBarry Smith } 212*01a81e61SBarry Smith EXTERN_C_END 213*01a81e61SBarry Smith 214*01a81e61SBarry Smith /**********************************************************************/ 215*01a81e61SBarry Smith 216*01a81e61SBarry Smith EXTERN_C_BEGIN 217*01a81e61SBarry Smith #undef __FUNCT__ 218*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew_SPAI" 219*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew_SPAI(PC pc,int maxnew1) 220*01a81e61SBarry Smith { 221*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 222*01a81e61SBarry Smith PetscFunctionBegin; 223*01a81e61SBarry Smith ispai->maxnew = maxnew1; 224*01a81e61SBarry Smith PetscFunctionReturn(0); 225*01a81e61SBarry Smith } 226*01a81e61SBarry Smith EXTERN_C_END 227*01a81e61SBarry Smith 228*01a81e61SBarry Smith /**********************************************************************/ 229*01a81e61SBarry Smith 230*01a81e61SBarry Smith EXTERN_C_BEGIN 231*01a81e61SBarry Smith #undef __FUNCT__ 232*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize_SPAI" 233*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize_SPAI(PC pc,int block_size1) 234*01a81e61SBarry Smith { 235*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 236*01a81e61SBarry Smith PetscFunctionBegin; 237*01a81e61SBarry Smith ispai->block_size = block_size1; 238*01a81e61SBarry Smith PetscFunctionReturn(0); 239*01a81e61SBarry Smith } 240*01a81e61SBarry Smith EXTERN_C_END 241*01a81e61SBarry Smith 242*01a81e61SBarry Smith /**********************************************************************/ 243*01a81e61SBarry Smith 244*01a81e61SBarry Smith EXTERN_C_BEGIN 245*01a81e61SBarry Smith #undef __FUNCT__ 246*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize_SPAI" 247*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize_SPAI(PC pc,int cache_size) 248*01a81e61SBarry Smith { 249*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 250*01a81e61SBarry Smith PetscFunctionBegin; 251*01a81e61SBarry Smith ispai->cache_size = cache_size; 252*01a81e61SBarry Smith PetscFunctionReturn(0); 253*01a81e61SBarry Smith } 254*01a81e61SBarry Smith EXTERN_C_END 255*01a81e61SBarry Smith 256*01a81e61SBarry Smith /**********************************************************************/ 257*01a81e61SBarry Smith 258*01a81e61SBarry Smith EXTERN_C_BEGIN 259*01a81e61SBarry Smith #undef __FUNCT__ 260*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose_SPAI" 261*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose_SPAI(PC pc,int verbose) 262*01a81e61SBarry Smith { 263*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 264*01a81e61SBarry Smith PetscFunctionBegin; 265*01a81e61SBarry Smith ispai->verbose = verbose; 266*01a81e61SBarry Smith PetscFunctionReturn(0); 267*01a81e61SBarry Smith } 268*01a81e61SBarry Smith EXTERN_C_END 269*01a81e61SBarry Smith 270*01a81e61SBarry Smith /**********************************************************************/ 271*01a81e61SBarry Smith 272*01a81e61SBarry Smith EXTERN_C_BEGIN 273*01a81e61SBarry Smith #undef __FUNCT__ 274*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp_SPAI" 275*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp_SPAI(PC pc,int sp) 276*01a81e61SBarry Smith { 277*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 278*01a81e61SBarry Smith PetscFunctionBegin; 279*01a81e61SBarry Smith ispai->sp = sp; 280*01a81e61SBarry Smith PetscFunctionReturn(0); 281*01a81e61SBarry Smith } 282*01a81e61SBarry Smith EXTERN_C_END 283*01a81e61SBarry Smith 284*01a81e61SBarry Smith /* -------------------------------------------------------------------*/ 285*01a81e61SBarry Smith 286*01a81e61SBarry Smith #undef __FUNCT__ 287*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetEpsilon" 288*01a81e61SBarry Smith /*@ 289*01a81e61SBarry Smith PCSPAISetEpsilon -- Set the tolerance for the SPAI preconditioner 290*01a81e61SBarry Smith 291*01a81e61SBarry Smith Input Parameters: 292*01a81e61SBarry Smith + pc - the preconditioner 293*01a81e61SBarry Smith - eps - epsilon (default .4) 294*01a81e61SBarry Smith 295*01a81e61SBarry Smith Notes: Espilon must be between 0 and 1. It controls the 296*01a81e61SBarry Smith quality of the approximation of M to the inverse of 297*01a81e61SBarry Smith A. Higher values of epsilon lead to more work, more 298*01a81e61SBarry Smith fill, and usually better preconditioners. In many 299*01a81e61SBarry Smith cases the best choice of epsilon is the one that 300*01a81e61SBarry Smith divides the total solution time equally between the 301*01a81e61SBarry Smith preconditioner and the solver. 302*01a81e61SBarry Smith 303*01a81e61SBarry Smith Level: intermediate 304*01a81e61SBarry Smith 305*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 306*01a81e61SBarry Smith @*/ 307*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetEpsilon(PC pc,double epsilon1) 308*01a81e61SBarry Smith { 309*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,double); 310*01a81e61SBarry Smith PetscFunctionBegin; 311*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetEpsilon_C",(void (**)(void))&f);CHKERRQ(ierr); 312*01a81e61SBarry Smith if (f) { 313*01a81e61SBarry Smith ierr = (*f)(pc,epsilon1);CHKERRQ(ierr); 314*01a81e61SBarry Smith } 315*01a81e61SBarry Smith PetscFunctionReturn(0); 316*01a81e61SBarry Smith } 317*01a81e61SBarry Smith 318*01a81e61SBarry Smith /**********************************************************************/ 319*01a81e61SBarry Smith 320*01a81e61SBarry Smith #undef __FUNCT__ 321*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetNBSteps" 322*01a81e61SBarry Smith /*@ 323*01a81e61SBarry Smith PCSPAISetNBSteps - set maximum number of improvement steps per row in 324*01a81e61SBarry Smith the SPAI preconditioner 325*01a81e61SBarry Smith 326*01a81e61SBarry Smith Input Parameters: 327*01a81e61SBarry Smith + pc - the preconditioner 328*01a81e61SBarry Smith - n - number of steps (default 5) 329*01a81e61SBarry Smith 330*01a81e61SBarry Smith Notes: SPAI constructs to approximation to every column of 331*01a81e61SBarry Smith the exact inverse of A in a series of improvement 332*01a81e61SBarry Smith steps. The quality of the approximation is determined 333*01a81e61SBarry Smith by epsilon. If an approximation achieving an accuracy 334*01a81e61SBarry Smith of epsilon is not obtained after ns steps, SPAI simply 335*01a81e61SBarry Smith uses the best approximation constructed so far. 336*01a81e61SBarry Smith 337*01a81e61SBarry Smith Level: intermediate 338*01a81e61SBarry Smith 339*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetMaxNew() 340*01a81e61SBarry Smith @*/ 341*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetNBSteps(PC pc,int nbsteps1) 342*01a81e61SBarry Smith { 343*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 344*01a81e61SBarry Smith PetscFunctionBegin; 345*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetNBSteps_C",(void (**)(void))&f);CHKERRQ(ierr); 346*01a81e61SBarry Smith if (f) { 347*01a81e61SBarry Smith ierr = (*f)(pc,nbsteps1);CHKERRQ(ierr); 348*01a81e61SBarry Smith } 349*01a81e61SBarry Smith PetscFunctionReturn(0); 350*01a81e61SBarry Smith } 351*01a81e61SBarry Smith 352*01a81e61SBarry Smith /**********************************************************************/ 353*01a81e61SBarry Smith 354*01a81e61SBarry Smith /* added 1/7/99 g.h. */ 355*01a81e61SBarry Smith #undef __FUNCT__ 356*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMax" 357*01a81e61SBarry Smith /*@ 358*01a81e61SBarry Smith PCSPAISetMax - set the size of various working buffers in 359*01a81e61SBarry Smith the SPAI preconditioner 360*01a81e61SBarry Smith 361*01a81e61SBarry Smith Input Parameters: 362*01a81e61SBarry Smith + pc - the preconditioner 363*01a81e61SBarry Smith - n - size (default is 5000) 364*01a81e61SBarry Smith 365*01a81e61SBarry Smith Level: intermediate 366*01a81e61SBarry Smith 367*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 368*01a81e61SBarry Smith @*/ 369*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMax(PC pc,int max1) 370*01a81e61SBarry Smith { 371*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 372*01a81e61SBarry Smith PetscFunctionBegin; 373*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMax_C",(void (**)(void))&f);CHKERRQ(ierr); 374*01a81e61SBarry Smith if (f) { 375*01a81e61SBarry Smith ierr = (*f)(pc,max1);CHKERRQ(ierr); 376*01a81e61SBarry Smith } 377*01a81e61SBarry Smith PetscFunctionReturn(0); 378*01a81e61SBarry Smith } 379*01a81e61SBarry Smith 380*01a81e61SBarry Smith /**********************************************************************/ 381*01a81e61SBarry Smith 382*01a81e61SBarry Smith #undef __FUNCT__ 383*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetMaxNew" 384*01a81e61SBarry Smith /*@ 385*01a81e61SBarry Smith PCSPAISetMaxNew - set maximum number of new nonzero candidates per step 386*01a81e61SBarry Smith in SPAI preconditioner 387*01a81e61SBarry Smith 388*01a81e61SBarry Smith Input Parameters: 389*01a81e61SBarry Smith + pc - the preconditioner 390*01a81e61SBarry Smith - n - maximum number (default 5) 391*01a81e61SBarry Smith 392*01a81e61SBarry Smith Level: intermediate 393*01a81e61SBarry Smith 394*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType(), PCSPAISetNBSteps() 395*01a81e61SBarry Smith @*/ 396*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetMaxNew(PC pc,int maxnew1) 397*01a81e61SBarry Smith { 398*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 399*01a81e61SBarry Smith PetscFunctionBegin; 400*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetMaxNew_C",(void (**)(void))&f);CHKERRQ(ierr); 401*01a81e61SBarry Smith if (f) { 402*01a81e61SBarry Smith ierr = (*f)(pc,maxnew1);CHKERRQ(ierr); 403*01a81e61SBarry Smith } 404*01a81e61SBarry Smith PetscFunctionReturn(0); 405*01a81e61SBarry Smith } 406*01a81e61SBarry Smith 407*01a81e61SBarry Smith /**********************************************************************/ 408*01a81e61SBarry Smith 409*01a81e61SBarry Smith #undef __FUNCT__ 410*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetBlockSize" 411*01a81e61SBarry Smith /*@ 412*01a81e61SBarry Smith PCSPAISetBlockSize - set the block size for the SPAI preconditioner 413*01a81e61SBarry Smith 414*01a81e61SBarry Smith Input Parameters: 415*01a81e61SBarry Smith + pc - the preconditioner 416*01a81e61SBarry Smith - n - block size (default 1) 417*01a81e61SBarry Smith 418*01a81e61SBarry Smith Notes: A block 419*01a81e61SBarry Smith size of 1 treats A as a matrix of scalar elements. A 420*01a81e61SBarry Smith block size of s > 1 treats A as a matrix of sxs 421*01a81e61SBarry Smith blocks. A block size of 0 treats A as a matrix with 422*01a81e61SBarry Smith variable sized blocks, which are determined by 423*01a81e61SBarry Smith searching for dense square diagonal blocks in A. 424*01a81e61SBarry Smith This can be very effective for finite-element 425*01a81e61SBarry Smith matrices. 426*01a81e61SBarry Smith 427*01a81e61SBarry Smith SPAI will convert A to block form, use a block 428*01a81e61SBarry Smith version of the preconditioner algorithm, and then 429*01a81e61SBarry Smith convert the result back to scalar form. 430*01a81e61SBarry Smith 431*01a81e61SBarry Smith In many cases the a block-size parameter other than 1 432*01a81e61SBarry Smith can lead to very significant improvement in 433*01a81e61SBarry Smith performance. 434*01a81e61SBarry Smith 435*01a81e61SBarry Smith 436*01a81e61SBarry Smith Level: intermediate 437*01a81e61SBarry Smith 438*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 439*01a81e61SBarry Smith @*/ 440*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetBlockSize(PC pc,int block_size1) 441*01a81e61SBarry Smith { 442*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 443*01a81e61SBarry Smith PetscFunctionBegin; 444*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetBlockSize_C",(void (**)(void))&f);CHKERRQ(ierr); 445*01a81e61SBarry Smith if (f) { 446*01a81e61SBarry Smith ierr = (*f)(pc,block_size1);CHKERRQ(ierr); 447*01a81e61SBarry Smith } 448*01a81e61SBarry Smith PetscFunctionReturn(0); 449*01a81e61SBarry Smith } 450*01a81e61SBarry Smith 451*01a81e61SBarry Smith /**********************************************************************/ 452*01a81e61SBarry Smith 453*01a81e61SBarry Smith #undef __FUNCT__ 454*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetCacheSize" 455*01a81e61SBarry Smith /*@ 456*01a81e61SBarry Smith PCSPAISetCacheSize - specify cache size in the SPAI preconditioner 457*01a81e61SBarry Smith 458*01a81e61SBarry Smith Input Parameters: 459*01a81e61SBarry Smith + pc - the preconditioner 460*01a81e61SBarry Smith - n - cache size {0,1,2,3,4,5} (default 5) 461*01a81e61SBarry Smith 462*01a81e61SBarry Smith Notes: SPAI uses a hash table to cache messages and avoid 463*01a81e61SBarry Smith redundant communication. If suggest always using 464*01a81e61SBarry Smith 5. This parameter is irrelevant in the serial 465*01a81e61SBarry Smith version. 466*01a81e61SBarry Smith 467*01a81e61SBarry Smith Level: intermediate 468*01a81e61SBarry Smith 469*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 470*01a81e61SBarry Smith @*/ 471*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetCacheSize(PC pc,int cache_size) 472*01a81e61SBarry Smith { 473*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 474*01a81e61SBarry Smith PetscFunctionBegin; 475*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetCacheSize_C",(void (**)(void))&f);CHKERRQ(ierr); 476*01a81e61SBarry Smith if (f) { 477*01a81e61SBarry Smith ierr = (*f)(pc,cache_size);CHKERRQ(ierr); 478*01a81e61SBarry Smith } 479*01a81e61SBarry Smith PetscFunctionReturn(0); 480*01a81e61SBarry Smith } 481*01a81e61SBarry Smith 482*01a81e61SBarry Smith /**********************************************************************/ 483*01a81e61SBarry Smith 484*01a81e61SBarry Smith #undef __FUNCT__ 485*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetVerbose" 486*01a81e61SBarry Smith /*@ 487*01a81e61SBarry Smith PCSPAISetVerbose - verbosity level for the SPAI preconditioner 488*01a81e61SBarry Smith 489*01a81e61SBarry Smith Input Parameters: 490*01a81e61SBarry Smith + pc - the preconditioner 491*01a81e61SBarry Smith - n - level (default 1) 492*01a81e61SBarry Smith 493*01a81e61SBarry Smith Notes: print parameters, timings and matrix statistics 494*01a81e61SBarry Smith 495*01a81e61SBarry Smith Level: intermediate 496*01a81e61SBarry Smith 497*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 498*01a81e61SBarry Smith @*/ 499*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetVerbose(PC pc,int verbose) 500*01a81e61SBarry Smith { 501*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 502*01a81e61SBarry Smith PetscFunctionBegin; 503*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetVerbose_C",(void (**)(void))&f);CHKERRQ(ierr); 504*01a81e61SBarry Smith if (f) { 505*01a81e61SBarry Smith ierr = (*f)(pc,verbose);CHKERRQ(ierr); 506*01a81e61SBarry Smith } 507*01a81e61SBarry Smith PetscFunctionReturn(0); 508*01a81e61SBarry Smith } 509*01a81e61SBarry Smith 510*01a81e61SBarry Smith /**********************************************************************/ 511*01a81e61SBarry Smith 512*01a81e61SBarry Smith #undef __FUNCT__ 513*01a81e61SBarry Smith #define __FUNCT__ "PCSPAISetSp" 514*01a81e61SBarry Smith /*@ 515*01a81e61SBarry Smith PCSPAISetSp - specify a symmetric matrix sparsity pattern in the SPAI preconditioner 516*01a81e61SBarry Smith 517*01a81e61SBarry Smith Input Parameters: 518*01a81e61SBarry Smith + pc - the preconditioner 519*01a81e61SBarry Smith - n - 0 or 1 520*01a81e61SBarry Smith 521*01a81e61SBarry Smith Notes: If A has a symmetric nonzero pattern use -sp 1 to 522*01a81e61SBarry Smith improve performance by eliminating some communication 523*01a81e61SBarry Smith in the parallel version. Even if A does not have a 524*01a81e61SBarry Smith symmetric nonzero pattern -sp 1 may well lead to good 525*01a81e61SBarry Smith results, but the code will not follow the published 526*01a81e61SBarry Smith SPAI algorithm exactly. 527*01a81e61SBarry Smith 528*01a81e61SBarry Smith 529*01a81e61SBarry Smith Level: intermediate 530*01a81e61SBarry Smith 531*01a81e61SBarry Smith .seealso: PCSPAI, PCSetType() 532*01a81e61SBarry Smith @*/ 533*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCSPAISetSp(PC pc,int sp) 534*01a81e61SBarry Smith { 535*01a81e61SBarry Smith PetscErrorCode ierr,(*f)(PC,int); 536*01a81e61SBarry Smith PetscFunctionBegin; 537*01a81e61SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSPAISetSp_C",(void (**)(void))&f);CHKERRQ(ierr); 538*01a81e61SBarry Smith if (f) { 539*01a81e61SBarry Smith ierr = (*f)(pc,sp);CHKERRQ(ierr); 540*01a81e61SBarry Smith } 541*01a81e61SBarry Smith PetscFunctionReturn(0); 542*01a81e61SBarry Smith } 543*01a81e61SBarry Smith 544*01a81e61SBarry Smith /**********************************************************************/ 545*01a81e61SBarry Smith 546*01a81e61SBarry Smith /**********************************************************************/ 547*01a81e61SBarry Smith 548*01a81e61SBarry Smith #undef __FUNCT__ 549*01a81e61SBarry Smith #define __FUNCT__ "PCSetFromOptions_SPAI" 550*01a81e61SBarry Smith static PetscErrorCode PCSetFromOptions_SPAI(PC pc) 551*01a81e61SBarry Smith { 552*01a81e61SBarry Smith PC_SPAI *ispai = (PC_SPAI*)pc->data; 553*01a81e61SBarry Smith PetscErrorCode ierr; 554*01a81e61SBarry Smith int nbsteps1,max1,maxnew1,block_size1,cache_size,verbose,sp; 555*01a81e61SBarry Smith double epsilon1; 556*01a81e61SBarry Smith PetscTruth flg; 557*01a81e61SBarry Smith 558*01a81e61SBarry Smith PetscFunctionBegin; 559*01a81e61SBarry Smith ierr = PetscOptionsHead("SPAI options");CHKERRQ(ierr); 560*01a81e61SBarry Smith ierr = PetscOptionsReal("-pc_spai_epsilon","","PCSPAISetEpsilon",ispai->epsilon,&epsilon1,&flg);CHKERRQ(ierr); 561*01a81e61SBarry Smith if (flg) { 562*01a81e61SBarry Smith ierr = PCSPAISetEpsilon(pc,epsilon1);CHKERRQ(ierr); 563*01a81e61SBarry Smith } 564*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_nbsteps","","PCSPAISetNBSteps",ispai->nbsteps,&nbsteps1,&flg);CHKERRQ(ierr); 565*01a81e61SBarry Smith if (flg) { 566*01a81e61SBarry Smith ierr = PCSPAISetNBSteps(pc,nbsteps1);CHKERRQ(ierr); 567*01a81e61SBarry Smith } 568*01a81e61SBarry Smith /* added 1/7/99 g.h. */ 569*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_max","","PCSPAISetMax",ispai->max,&max1,&flg);CHKERRQ(ierr); 570*01a81e61SBarry Smith if (flg) { 571*01a81e61SBarry Smith ierr = PCSPAISetMax(pc,max1);CHKERRQ(ierr); 572*01a81e61SBarry Smith } 573*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_maxnew","","PCSPAISetMaxNew",ispai->maxnew,&maxnew1,&flg);CHKERRQ(ierr); 574*01a81e61SBarry Smith if (flg) { 575*01a81e61SBarry Smith ierr = PCSPAISetMaxNew(pc,maxnew1);CHKERRQ(ierr); 576*01a81e61SBarry Smith } 577*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_block_size","","PCSPAISetBlockSize",ispai->block_size,&block_size1,&flg);CHKERRQ(ierr); 578*01a81e61SBarry Smith if (flg) { 579*01a81e61SBarry Smith ierr = PCSPAISetBlockSize(pc,block_size1);CHKERRQ(ierr); 580*01a81e61SBarry Smith } 581*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_cache_size","","PCSPAISetCacheSize",ispai->cache_size,&cache_size,&flg);CHKERRQ(ierr); 582*01a81e61SBarry Smith if (flg) { 583*01a81e61SBarry Smith ierr = PCSPAISetCacheSize(pc,cache_size);CHKERRQ(ierr); 584*01a81e61SBarry Smith } 585*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_verbose","","PCSPAISetVerbose",ispai->verbose,&verbose,&flg);CHKERRQ(ierr); 586*01a81e61SBarry Smith if (flg) { 587*01a81e61SBarry Smith ierr = PCSPAISetVerbose(pc,verbose);CHKERRQ(ierr); 588*01a81e61SBarry Smith } 589*01a81e61SBarry Smith ierr = PetscOptionsInt("-pc_spai_sp","","PCSPAISetSp",ispai->sp,&sp,&flg);CHKERRQ(ierr); 590*01a81e61SBarry Smith if (flg) { 591*01a81e61SBarry Smith ierr = PCSPAISetSp(pc,sp);CHKERRQ(ierr); 592*01a81e61SBarry Smith } 593*01a81e61SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 594*01a81e61SBarry Smith PetscFunctionReturn(0); 595*01a81e61SBarry Smith } 596*01a81e61SBarry Smith 597*01a81e61SBarry Smith /**********************************************************************/ 598*01a81e61SBarry Smith 599*01a81e61SBarry Smith /*MC 600*01a81e61SBarry Smith PCSPAI - Use the Sparse Approximate Inverse method of Grote and Barnard 601*01a81e61SBarry Smith as a preconditioner (SIAM J. Sci. Comput.; vol 18, nr 3) 602*01a81e61SBarry Smith 603*01a81e61SBarry Smith Options Database Keys: 604*01a81e61SBarry Smith + -pc_spai_set_epsilon <eps> - set tolerance 605*01a81e61SBarry Smith . -pc_spai_nbstep <n> - set nbsteps 606*01a81e61SBarry Smith . -pc_spai_max <m> - set max 607*01a81e61SBarry Smith . -pc_spai_max_new <m> - set maxnew 608*01a81e61SBarry Smith . -pc_spai_block_size <n> - set block size 609*01a81e61SBarry Smith . -pc_spai_cache_size <n> - set cache size 610*01a81e61SBarry Smith . -pc_spai_sp <m> - set sp 611*01a81e61SBarry Smith - -pc_spai_set_verbose <true,false> - verbose output 612*01a81e61SBarry Smith 613*01a81e61SBarry Smith Notes: This only works with AIJ matrices. 614*01a81e61SBarry Smith 615*01a81e61SBarry Smith Level: beginner 616*01a81e61SBarry Smith 617*01a81e61SBarry Smith Concepts: approximate inverse 618*01a81e61SBarry Smith 619*01a81e61SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 620*01a81e61SBarry Smith PCSPAISetEpsilon(), PCSPAISetMax(), PCSPAISetMaxNew(), PCSPAISetBlockSize(), 621*01a81e61SBarry Smith PCSPAISetVerbose(), PCSPAISetSp() 622*01a81e61SBarry Smith M*/ 623*01a81e61SBarry Smith 624*01a81e61SBarry Smith EXTERN_C_BEGIN 625*01a81e61SBarry Smith #undef __FUNCT__ 626*01a81e61SBarry Smith #define __FUNCT__ "PCCreate_SPAI" 627*01a81e61SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SPAI(PC pc) 628*01a81e61SBarry Smith { 629*01a81e61SBarry Smith PC_SPAI *ispai; 630*01a81e61SBarry Smith PetscErrorCode ierr; 631*01a81e61SBarry Smith 632*01a81e61SBarry Smith PetscFunctionBegin; 633*01a81e61SBarry Smith ierr = PetscNew(PC_SPAI,&ispai);CHKERRQ(ierr); 634*01a81e61SBarry Smith pc->data = (void*)ispai; 635*01a81e61SBarry Smith 636*01a81e61SBarry Smith pc->ops->destroy = PCDestroy_SPAI; 637*01a81e61SBarry Smith pc->ops->apply = PCApply_SPAI; 638*01a81e61SBarry Smith pc->ops->applyrichardson = 0; 639*01a81e61SBarry Smith pc->ops->setup = PCSetUp_SPAI; 640*01a81e61SBarry Smith pc->ops->view = PCView_SPAI; 641*01a81e61SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SPAI; 642*01a81e61SBarry Smith 643*01a81e61SBarry Smith pc->name = 0; 644*01a81e61SBarry Smith ispai->epsilon = .4; 645*01a81e61SBarry Smith ispai->nbsteps = 5; 646*01a81e61SBarry Smith ispai->max = 5000; 647*01a81e61SBarry Smith ispai->maxnew = 5; 648*01a81e61SBarry Smith ispai->block_size = 1; 649*01a81e61SBarry Smith ispai->cache_size = 5; 650*01a81e61SBarry Smith ispai->verbose = 0; 651*01a81e61SBarry Smith 652*01a81e61SBarry Smith ispai->sp = 1; 653*01a81e61SBarry Smith ierr = MPI_Comm_dup(pc->comm,&(ispai->comm_spai));CHKERRQ(ierr); 654*01a81e61SBarry Smith 655*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetEpsilon_C", 656*01a81e61SBarry Smith "PCSPAISetEpsilon_SPAI", 657*01a81e61SBarry Smith PCSPAISetEpsilon_SPAI);CHKERRQ(ierr); 658*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetNBSteps_C", 659*01a81e61SBarry Smith "PCSPAISetNBSteps_SPAI", 660*01a81e61SBarry Smith PCSPAISetNBSteps_SPAI);CHKERRQ(ierr); 661*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMax_C", 662*01a81e61SBarry Smith "PCSPAISetMax_SPAI", 663*01a81e61SBarry Smith PCSPAISetMax_SPAI);CHKERRQ(ierr); 664*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetMaxNew_CC", 665*01a81e61SBarry Smith "PCSPAISetMaxNew_SPAI", 666*01a81e61SBarry Smith PCSPAISetMaxNew_SPAI);CHKERRQ(ierr); 667*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetBlockSize_C", 668*01a81e61SBarry Smith "PCSPAISetBlockSize_SPAI", 669*01a81e61SBarry Smith PCSPAISetBlockSize_SPAI);CHKERRQ(ierr); 670*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetCacheSize_C", 671*01a81e61SBarry Smith "PCSPAISetCacheSize_SPAI", 672*01a81e61SBarry Smith PCSPAISetCacheSize_SPAI);CHKERRQ(ierr); 673*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetVerbose_C", 674*01a81e61SBarry Smith "PCSPAISetVerbose_SPAI", 675*01a81e61SBarry Smith PCSPAISetVerbose_SPAI);CHKERRQ(ierr); 676*01a81e61SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSPAISetSp_C", 677*01a81e61SBarry Smith "PCSPAISetSp_SPAI", 678*01a81e61SBarry Smith PCSPAISetSp_SPAI);CHKERRQ(ierr); 679*01a81e61SBarry Smith 680*01a81e61SBarry Smith PetscFunctionReturn(0); 681*01a81e61SBarry Smith } 682*01a81e61SBarry Smith EXTERN_C_END 683*01a81e61SBarry Smith 684*01a81e61SBarry Smith /**********************************************************************/ 685*01a81e61SBarry Smith 686*01a81e61SBarry Smith /* 687*01a81e61SBarry Smith Converts from a PETSc matrix to an SPAI matrix 688*01a81e61SBarry Smith */ 689*01a81e61SBarry Smith #undef __FUNCT__ 690*01a81e61SBarry Smith #define __FUNCT__ "ConvertMatToMatrix" 691*01a81e61SBarry Smith PetscErrorCode ConvertMatToMatrix(MPI_Comm comm, Mat A,Mat AT,matrix **B) 692*01a81e61SBarry Smith { 693*01a81e61SBarry Smith matrix *M; 694*01a81e61SBarry Smith int i,j,col; 695*01a81e61SBarry Smith int row_indx; 696*01a81e61SBarry Smith int len,pe,local_indx,start_indx; 697*01a81e61SBarry Smith int *mapping; 698*01a81e61SBarry Smith PetscErrorCode ierr; 699*01a81e61SBarry Smith const int *cols; 700*01a81e61SBarry Smith const double *vals; 701*01a81e61SBarry Smith int *num_ptr,n,mnl,nnl,rank,size,nz,rstart,rend; 702*01a81e61SBarry Smith struct compressed_lines *rows; 703*01a81e61SBarry Smith 704*01a81e61SBarry Smith PetscFunctionBegin; 705*01a81e61SBarry Smith 706*01a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 707*01a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 708*01a81e61SBarry Smith ierr = MatGetSize(A,&n,&n);CHKERRQ(ierr); 709*01a81e61SBarry Smith ierr = MatGetLocalSize(A,&mnl,&nnl);CHKERRQ(ierr); 710*01a81e61SBarry Smith 711*01a81e61SBarry Smith /* 712*01a81e61SBarry Smith not sure why a barrier is required. commenting out 713*01a81e61SBarry Smith ierr = MPI_Barrier(comm);CHKERRQ(ierr); 714*01a81e61SBarry Smith */ 715*01a81e61SBarry Smith 716*01a81e61SBarry Smith M = new_matrix((void*)comm); 717*01a81e61SBarry Smith 718*01a81e61SBarry Smith M->n = n; 719*01a81e61SBarry Smith M->bs = 1; 720*01a81e61SBarry Smith M->max_block_size = 1; 721*01a81e61SBarry Smith 722*01a81e61SBarry Smith M->mnls = (int*)malloc(sizeof(int)*size); 723*01a81e61SBarry Smith M->start_indices = (int*)malloc(sizeof(int)*size); 724*01a81e61SBarry Smith M->pe = (int*)malloc(sizeof(int)*n); 725*01a81e61SBarry Smith M->block_sizes = (int*)malloc(sizeof(int)*n); 726*01a81e61SBarry Smith for (i=0; i<n; i++) M->block_sizes[i] = 1; 727*01a81e61SBarry Smith 728*01a81e61SBarry Smith ierr = MPI_Allgather(&mnl,1,MPI_INT,M->mnls,1,MPI_INT,comm);CHKERRQ(ierr); 729*01a81e61SBarry Smith 730*01a81e61SBarry Smith M->start_indices[0] = 0; 731*01a81e61SBarry Smith for (i=1; i<size; i++) { 732*01a81e61SBarry Smith M->start_indices[i] = M->start_indices[i-1] + M->mnls[i-1]; 733*01a81e61SBarry Smith } 734*01a81e61SBarry Smith 735*01a81e61SBarry Smith M->mnl = M->mnls[M->myid]; 736*01a81e61SBarry Smith M->my_start_index = M->start_indices[M->myid]; 737*01a81e61SBarry Smith 738*01a81e61SBarry Smith for (i=0; i<size; i++) { 739*01a81e61SBarry Smith start_indx = M->start_indices[i]; 740*01a81e61SBarry Smith for (j=0; j<M->mnls[i]; j++) 741*01a81e61SBarry Smith M->pe[start_indx+j] = i; 742*01a81e61SBarry Smith } 743*01a81e61SBarry Smith 744*01a81e61SBarry Smith if (AT) { 745*01a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],1);CHKERRQ(ierr); 746*01a81e61SBarry Smith } else { 747*01a81e61SBarry Smith M->lines = new_compressed_lines(M->mnls[rank],0);CHKERRQ(ierr); 748*01a81e61SBarry Smith } 749*01a81e61SBarry Smith 750*01a81e61SBarry Smith rows = M->lines; 751*01a81e61SBarry Smith 752*01a81e61SBarry Smith /* Determine the mapping from global indices to pointers */ 753*01a81e61SBarry Smith ierr = PetscMalloc(M->n*sizeof(int),&mapping);CHKERRQ(ierr); 754*01a81e61SBarry Smith pe = 0; 755*01a81e61SBarry Smith local_indx = 0; 756*01a81e61SBarry Smith for (i=0; i<M->n; i++) { 757*01a81e61SBarry Smith if (local_indx >= M->mnls[pe]) { 758*01a81e61SBarry Smith pe++; 759*01a81e61SBarry Smith local_indx = 0; 760*01a81e61SBarry Smith } 761*01a81e61SBarry Smith mapping[i] = local_indx + M->start_indices[pe]; 762*01a81e61SBarry Smith local_indx++; 763*01a81e61SBarry Smith } 764*01a81e61SBarry Smith 765*01a81e61SBarry Smith 766*01a81e61SBarry Smith ierr = PetscMalloc(mnl*sizeof(int),&num_ptr);CHKERRQ(ierr); 767*01a81e61SBarry Smith 768*01a81e61SBarry Smith /*********************************************************/ 769*01a81e61SBarry Smith /************** Set up the row structure *****************/ 770*01a81e61SBarry Smith /*********************************************************/ 771*01a81e61SBarry Smith 772*01a81e61SBarry Smith /* count number of nonzeros in every row */ 773*01a81e61SBarry Smith ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 774*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 775*01a81e61SBarry Smith ierr = MatGetRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 776*01a81e61SBarry Smith ierr = MatRestoreRow(A,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 777*01a81e61SBarry Smith } 778*01a81e61SBarry Smith 779*01a81e61SBarry Smith /* allocate buffers */ 780*01a81e61SBarry Smith len = 0; 781*01a81e61SBarry Smith for (i=0; i<mnl; i++) { 782*01a81e61SBarry Smith if (len < num_ptr[i]) len = num_ptr[i]; 783*01a81e61SBarry Smith } 784*01a81e61SBarry Smith 785*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 786*01a81e61SBarry Smith row_indx = i-rstart; 787*01a81e61SBarry Smith len = num_ptr[row_indx]; 788*01a81e61SBarry Smith rows->ptrs[row_indx] = (int*)malloc(len*sizeof(int)); 789*01a81e61SBarry Smith rows->A[row_indx] = (double*)malloc(len*sizeof(double)); 790*01a81e61SBarry Smith } 791*01a81e61SBarry Smith 792*01a81e61SBarry Smith /* copy the matrix */ 793*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 794*01a81e61SBarry Smith row_indx = i - rstart; 795*01a81e61SBarry Smith ierr = MatGetRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 796*01a81e61SBarry Smith for (j=0; j<nz; j++) { 797*01a81e61SBarry Smith col = cols[j]; 798*01a81e61SBarry Smith len = rows->len[row_indx]++; 799*01a81e61SBarry Smith rows->ptrs[row_indx][len] = mapping[col]; 800*01a81e61SBarry Smith rows->A[row_indx][len] = vals[j]; 801*01a81e61SBarry Smith } 802*01a81e61SBarry Smith rows->slen[row_indx] = rows->len[row_indx]; 803*01a81e61SBarry Smith ierr = MatRestoreRow(A,i,&nz,&cols,&vals);CHKERRQ(ierr); 804*01a81e61SBarry Smith } 805*01a81e61SBarry Smith 806*01a81e61SBarry Smith 807*01a81e61SBarry Smith /************************************************************/ 808*01a81e61SBarry Smith /************** Set up the column structure *****************/ 809*01a81e61SBarry Smith /*********************************************************/ 810*01a81e61SBarry Smith 811*01a81e61SBarry Smith if (AT) { 812*01a81e61SBarry Smith 813*01a81e61SBarry Smith /* count number of nonzeros in every column */ 814*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 815*01a81e61SBarry Smith ierr = MatGetRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 816*01a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&num_ptr[i-rstart],PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 817*01a81e61SBarry Smith } 818*01a81e61SBarry Smith 819*01a81e61SBarry Smith /* allocate buffers */ 820*01a81e61SBarry Smith len = 0; 821*01a81e61SBarry Smith for (i=0; i<mnl; i++) { 822*01a81e61SBarry Smith if (len < num_ptr[i]) len = num_ptr[i]; 823*01a81e61SBarry Smith } 824*01a81e61SBarry Smith 825*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 826*01a81e61SBarry Smith row_indx = i-rstart; 827*01a81e61SBarry Smith len = num_ptr[row_indx]; 828*01a81e61SBarry Smith rows->rptrs[row_indx] = (int*)malloc(len*sizeof(int)); 829*01a81e61SBarry Smith } 830*01a81e61SBarry Smith 831*01a81e61SBarry Smith /* copy the matrix (i.e., the structure) */ 832*01a81e61SBarry Smith for (i=rstart; i<rend; i++) { 833*01a81e61SBarry Smith row_indx = i - rstart; 834*01a81e61SBarry Smith ierr = MatGetRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 835*01a81e61SBarry Smith for (j=0; j<nz; j++) { 836*01a81e61SBarry Smith col = cols[j]; 837*01a81e61SBarry Smith len = rows->rlen[row_indx]++; 838*01a81e61SBarry Smith rows->rptrs[row_indx][len] = mapping[col]; 839*01a81e61SBarry Smith } 840*01a81e61SBarry Smith ierr = MatRestoreRow(AT,i,&nz,&cols,&vals);CHKERRQ(ierr); 841*01a81e61SBarry Smith } 842*01a81e61SBarry Smith } 843*01a81e61SBarry Smith 844*01a81e61SBarry Smith ierr = PetscFree(num_ptr);CHKERRQ(ierr); 845*01a81e61SBarry Smith ierr = PetscFree(mapping);CHKERRQ(ierr); 846*01a81e61SBarry Smith 847*01a81e61SBarry Smith order_pointers(M); 848*01a81e61SBarry Smith M->maxnz = calc_maxnz(M); 849*01a81e61SBarry Smith 850*01a81e61SBarry Smith *B = M; 851*01a81e61SBarry Smith 852*01a81e61SBarry Smith PetscFunctionReturn(0); 853*01a81e61SBarry Smith } 854*01a81e61SBarry Smith 855*01a81e61SBarry Smith /**********************************************************************/ 856*01a81e61SBarry Smith 857*01a81e61SBarry Smith /* 858*01a81e61SBarry Smith Converts from an SPAI matrix B to a PETSc matrix PB. 859*01a81e61SBarry Smith This assumes that the the SPAI matrix B is stored in 860*01a81e61SBarry Smith COMPRESSED-ROW format. 861*01a81e61SBarry Smith */ 862*01a81e61SBarry Smith #undef __FUNCT__ 863*01a81e61SBarry Smith #define __FUNCT__ "ConvertMatrixToMat" 864*01a81e61SBarry Smith PetscErrorCode ConvertMatrixToMat(MPI_Comm comm,matrix *B,Mat *PB) 865*01a81e61SBarry Smith { 866*01a81e61SBarry Smith int size,rank; 867*01a81e61SBarry Smith PetscErrorCode ierr; 868*01a81e61SBarry Smith int m,n,M,N; 869*01a81e61SBarry Smith int d_nz,o_nz; 870*01a81e61SBarry Smith int *d_nnz,*o_nnz; 871*01a81e61SBarry Smith int i,k,global_row,global_col,first_diag_col,last_diag_col; 872*01a81e61SBarry Smith PetscScalar val; 873*01a81e61SBarry Smith 874*01a81e61SBarry Smith PetscFunctionBegin; 875*01a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 876*01a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 877*01a81e61SBarry Smith 878*01a81e61SBarry Smith m = n = B->mnls[rank]; 879*01a81e61SBarry Smith d_nz = o_nz = 0; 880*01a81e61SBarry Smith 881*01a81e61SBarry Smith /* Determine preallocation for MatCreateMPIAIJ */ 882*01a81e61SBarry Smith ierr = PetscMalloc(m*sizeof(int),&d_nnz);CHKERRQ(ierr); 883*01a81e61SBarry Smith ierr = PetscMalloc(m*sizeof(int),&o_nnz);CHKERRQ(ierr); 884*01a81e61SBarry Smith for (i=0; i<m; i++) d_nnz[i] = o_nnz[i] = 0; 885*01a81e61SBarry Smith first_diag_col = B->start_indices[rank]; 886*01a81e61SBarry Smith last_diag_col = first_diag_col + B->mnls[rank]; 887*01a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 888*01a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 889*01a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 890*01a81e61SBarry Smith if ((global_col >= first_diag_col) && (global_col <= last_diag_col)) 891*01a81e61SBarry Smith d_nnz[i]++; 892*01a81e61SBarry Smith else 893*01a81e61SBarry Smith o_nnz[i]++; 894*01a81e61SBarry Smith } 895*01a81e61SBarry Smith } 896*01a81e61SBarry Smith 897*01a81e61SBarry Smith M = N = B->n; 898*01a81e61SBarry Smith /* Here we only know how to create AIJ format */ 899*01a81e61SBarry Smith ierr = MatCreate(comm,m,n,M,N,PB);CHKERRQ(ierr); 900*01a81e61SBarry Smith ierr = MatSetType(*PB,MATAIJ);CHKERRQ(ierr); 901*01a81e61SBarry Smith ierr = MatSeqAIJSetPreallocation(*PB,d_nz,d_nnz);CHKERRQ(ierr); 902*01a81e61SBarry Smith ierr = MatMPIAIJSetPreallocation(*PB,d_nz,d_nnz,o_nz,o_nnz);CHKERRQ(ierr); 903*01a81e61SBarry Smith 904*01a81e61SBarry Smith for (i=0; i<B->mnls[rank]; i++) { 905*01a81e61SBarry Smith global_row = B->start_indices[rank]+i; 906*01a81e61SBarry Smith for (k=0; k<B->lines->len[i]; k++) { 907*01a81e61SBarry Smith global_col = B->lines->ptrs[i][k]; 908*01a81e61SBarry Smith val = B->lines->A[i][k]; 909*01a81e61SBarry Smith ierr = MatSetValues(*PB,1,&global_row,1,&global_col,&val,ADD_VALUES);CHKERRQ(ierr); 910*01a81e61SBarry Smith } 911*01a81e61SBarry Smith } 912*01a81e61SBarry Smith 913*01a81e61SBarry Smith ierr = PetscFree(d_nnz);CHKERRQ(ierr); 914*01a81e61SBarry Smith ierr = PetscFree(o_nnz);CHKERRQ(ierr); 915*01a81e61SBarry Smith 916*01a81e61SBarry Smith ierr = MatAssemblyBegin(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 917*01a81e61SBarry Smith ierr = MatAssemblyEnd(*PB,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 918*01a81e61SBarry Smith 919*01a81e61SBarry Smith PetscFunctionReturn(0); 920*01a81e61SBarry Smith } 921*01a81e61SBarry Smith 922*01a81e61SBarry Smith /**********************************************************************/ 923*01a81e61SBarry Smith 924*01a81e61SBarry Smith /* 925*01a81e61SBarry Smith Converts from an SPAI vector v to a PETSc vec Pv. 926*01a81e61SBarry Smith */ 927*01a81e61SBarry Smith #undef __FUNCT__ 928*01a81e61SBarry Smith #define __FUNCT__ "ConvertVectorToVec" 929*01a81e61SBarry Smith PetscErrorCode ConvertVectorToVec(MPI_Comm comm,vector *v,Vec *Pv) 930*01a81e61SBarry Smith { 931*01a81e61SBarry Smith PetscErrorCode ierr; 932*01a81e61SBarry Smith int size,rank,m,M,i,*mnls,*start_indices,*global_indices; 933*01a81e61SBarry Smith 934*01a81e61SBarry Smith PetscFunctionBegin; 935*01a81e61SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 936*01a81e61SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 937*01a81e61SBarry Smith 938*01a81e61SBarry Smith m = v->mnl; 939*01a81e61SBarry Smith M = v->n; 940*01a81e61SBarry Smith 941*01a81e61SBarry Smith 942*01a81e61SBarry Smith ierr = VecCreateMPI(comm,m,M,Pv);CHKERRQ(ierr); 943*01a81e61SBarry Smith 944*01a81e61SBarry Smith ierr = PetscMalloc(size*sizeof(int),&mnls);CHKERRQ(ierr); 945*01a81e61SBarry Smith ierr = MPI_Allgather((void*)&v->mnl,1,MPI_INT,(void*)mnls,1,MPI_INT,comm);CHKERRQ(ierr); 946*01a81e61SBarry Smith 947*01a81e61SBarry Smith ierr = PetscMalloc(size*sizeof(int),&start_indices);CHKERRQ(ierr); 948*01a81e61SBarry Smith start_indices[0] = 0; 949*01a81e61SBarry Smith for (i=1; i<size; i++) 950*01a81e61SBarry Smith start_indices[i] = start_indices[i-1] +mnls[i-1]; 951*01a81e61SBarry Smith 952*01a81e61SBarry Smith ierr = PetscMalloc(v->mnl*sizeof(int),&global_indices);CHKERRQ(ierr); 953*01a81e61SBarry Smith for (i=0; i<v->mnl; i++) 954*01a81e61SBarry Smith global_indices[i] = start_indices[rank] + i; 955*01a81e61SBarry Smith 956*01a81e61SBarry Smith ierr = PetscFree(mnls);CHKERRQ(ierr); 957*01a81e61SBarry Smith ierr = PetscFree(start_indices);CHKERRQ(ierr); 958*01a81e61SBarry Smith 959*01a81e61SBarry Smith ierr = VecSetValues(*Pv,v->mnl,global_indices,v->v,INSERT_VALUES);CHKERRQ(ierr); 960*01a81e61SBarry Smith 961*01a81e61SBarry Smith ierr = PetscFree(global_indices);CHKERRQ(ierr); 962*01a81e61SBarry Smith 963*01a81e61SBarry Smith ierr = VecAssemblyBegin(*Pv);CHKERRQ(ierr); 964*01a81e61SBarry Smith ierr = VecAssemblyEnd(*Pv);CHKERRQ(ierr); 965*01a81e61SBarry Smith 966*01a81e61SBarry Smith PetscFunctionReturn(0); 967*01a81e61SBarry Smith } 968*01a81e61SBarry Smith 969*01a81e61SBarry Smith 970*01a81e61SBarry Smith 971*01a81e61SBarry Smith 972