1f3161b27SJose E. Roman #define PETSCKSP_DLL 2f3161b27SJose E. Roman 3f3161b27SJose E. Roman /* 4f3161b27SJose E. Roman Provides an interface to pARMS. 5f3161b27SJose E. Roman Requires pARMS 3.2 or later. 6f3161b27SJose E. Roman */ 7f3161b27SJose E. Roman 8b45d2f2cSJed Brown #include "petsc-private/pcimpl.h" /*I "petscpc.h" I*/ 9f3161b27SJose E. Roman 10519f805aSKarl Rupp #if defined(PETSC_USE_COMPLEX) 11f3161b27SJose E. Roman #define DBL_CMPLX 12f3161b27SJose E. Roman #else 13f3161b27SJose E. Roman #define DBL 14f3161b27SJose E. Roman #endif 15f3161b27SJose E. Roman #define USE_MPI 16f3161b27SJose E. Roman #define REAL double 17f3161b27SJose E. Roman #define HAS_BLAS 18f3161b27SJose E. Roman #define FORTRAN_UNDERSCORE 19f3161b27SJose E. Roman #include "parms_sys.h" 20f3161b27SJose E. Roman #undef FLOAT 21f3161b27SJose E. Roman #define FLOAT PetscScalar 22f3161b27SJose E. Roman #include "parms.h" 23f3161b27SJose E. Roman 24f3161b27SJose E. Roman /* 25f3161b27SJose E. Roman Private context (data structure) for the preconditioner. 26f3161b27SJose E. Roman */ 27f3161b27SJose E. Roman typedef struct { 28f3161b27SJose E. Roman parms_Map map; 29f3161b27SJose E. Roman parms_Mat A; 30f3161b27SJose E. Roman parms_PC pc; 31f3161b27SJose E. Roman PCPARMSGlobalType global; 32f3161b27SJose E. Roman PCPARMSLocalType local; 33f3161b27SJose E. Roman PetscInt levels, blocksize, maxdim, maxits, lfil[7]; 34f3161b27SJose E. Roman PetscBool nonsymperm, meth[8]; 35f3161b27SJose E. Roman PetscReal solvetol, indtol, droptol[7]; 36f3161b27SJose E. Roman PetscScalar *lvec0, *lvec1; 37f3161b27SJose E. Roman } PC_PARMS; 38f3161b27SJose E. Roman 39f3161b27SJose E. Roman 40f3161b27SJose E. Roman #undef __FUNCT__ 41f3161b27SJose E. Roman #define __FUNCT__ "PCSetUp_PARMS" 42f3161b27SJose E. Roman static PetscErrorCode PCSetUp_PARMS(PC pc) 43f3161b27SJose E. Roman { 44f3161b27SJose E. Roman Mat pmat; 45f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 46f3161b27SJose E. Roman const PetscInt *mapptr0; 47f3161b27SJose E. Roman PetscInt n, lsize, low, high, i, pos, ncols, length; 48f3161b27SJose E. Roman int *maptmp, *mapptr, *ia, *ja, *ja1, *im; 49f3161b27SJose E. Roman PetscScalar *aa, *aa1; 50f3161b27SJose E. Roman const PetscInt *cols; 51f3161b27SJose E. Roman PetscInt meth[8]; 52f3161b27SJose E. Roman const PetscScalar *values; 53f3161b27SJose E. Roman PetscErrorCode ierr; 54f3161b27SJose E. Roman MatInfo matinfo; 55f3161b27SJose E. Roman PetscMPIInt rank, npro; 56f3161b27SJose E. Roman 57f3161b27SJose E. Roman PetscFunctionBegin; 58f3161b27SJose E. Roman /* Get preconditioner matrix from PETSc and setup pARMS structs */ 590298fd71SBarry Smith ierr = PCGetOperators(pc,NULL,&pmat,NULL);CHKERRQ(ierr); 60*ce94432eSBarry Smith MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro); 61*ce94432eSBarry Smith MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank); 62f3161b27SJose E. Roman 630298fd71SBarry Smith ierr = MatGetSize(pmat,&n,NULL);CHKERRQ(ierr); 64f3161b27SJose E. Roman ierr = PetscMalloc((npro+1)*sizeof(int),&mapptr);CHKERRQ(ierr); 65f3161b27SJose E. Roman ierr = PetscMalloc(n*sizeof(int),&maptmp);CHKERRQ(ierr); 66f3161b27SJose E. Roman ierr = MatGetOwnershipRanges(pmat,&mapptr0);CHKERRQ(ierr); 67f3161b27SJose E. Roman low = mapptr0[rank]; 68f3161b27SJose E. Roman high = mapptr0[rank+1]; 69f3161b27SJose E. Roman lsize = high - low; 70f3161b27SJose E. Roman 712fa5cd67SKarl Rupp for (i=0; i<npro+1; i++) mapptr[i] = mapptr0[i]+1; 722fa5cd67SKarl Rupp for (i = 0; i<n; i++) maptmp[i] = i+1; 73f3161b27SJose E. Roman 74f3161b27SJose E. Roman /* if created, destroy the previous map */ 75f3161b27SJose E. Roman if (parms->map) { 76f3161b27SJose E. Roman parms_MapFree(&parms->map); 770298fd71SBarry Smith parms->map = NULL; 78f3161b27SJose E. Roman } 79f3161b27SJose E. Roman 80f3161b27SJose E. Roman /* create pARMS map object */ 81*ce94432eSBarry Smith parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,PetscObjectComm((PetscObject)pmat),1,NONINTERLACED); 82f3161b27SJose E. Roman 83f3161b27SJose E. Roman /* if created, destroy the previous pARMS matrix */ 84f3161b27SJose E. Roman if (parms->A) { 85f3161b27SJose E. Roman parms_MatFree(&parms->A); 860298fd71SBarry Smith parms->A = NULL; 87f3161b27SJose E. Roman } 88f3161b27SJose E. Roman 89f3161b27SJose E. Roman /* create pARMS mat object */ 90f3161b27SJose E. Roman parms_MatCreate(&parms->A,parms->map); 91f3161b27SJose E. Roman 92f3161b27SJose E. Roman /* setup and copy csr data structure for pARMS */ 93f3161b27SJose E. Roman ierr = PetscMalloc((lsize+1)*sizeof(int),&ia);CHKERRQ(ierr); 94f3161b27SJose E. Roman ia[0] = 1; 95f3161b27SJose E. Roman ierr = MatGetInfo(pmat,MAT_LOCAL,&matinfo);CHKERRQ(ierr); 96f3161b27SJose E. Roman length = matinfo.nz_used; 97f3161b27SJose E. Roman ierr = PetscMalloc(length*sizeof(int),&ja);CHKERRQ(ierr); 98f3161b27SJose E. Roman ierr = PetscMalloc(length*sizeof(PetscScalar),&aa);CHKERRQ(ierr); 99f3161b27SJose E. Roman 100f3161b27SJose E. Roman for (i = low; i<high; i++) { 101f3161b27SJose E. Roman pos = ia[i-low]-1; 102f3161b27SJose E. Roman ierr = MatGetRow(pmat,i,&ncols,&cols,&values);CHKERRQ(ierr); 103f3161b27SJose E. Roman ia[i-low+1] = ia[i-low] + ncols; 104f3161b27SJose E. Roman 105f3161b27SJose E. Roman if (ia[i-low+1] >= length) { 106f3161b27SJose E. Roman length += ncols; 107f3161b27SJose E. Roman ierr = PetscMalloc(length*sizeof(int),&ja1);CHKERRQ(ierr); 108f3161b27SJose E. Roman ierr = PetscMemcpy(ja1,ja,(ia[i-low]-1)*sizeof(int));CHKERRQ(ierr); 109f3161b27SJose E. Roman ierr = PetscFree(ja);CHKERRQ(ierr); 110f3161b27SJose E. Roman ja = ja1; 111f3161b27SJose E. Roman ierr = PetscMalloc(length*sizeof(PetscScalar),&aa1);CHKERRQ(ierr); 112f3161b27SJose E. Roman ierr = PetscMemcpy(aa1,aa,(ia[i-low]-1)*sizeof(PetscScalar));CHKERRQ(ierr); 113f3161b27SJose E. Roman ierr = PetscFree(aa);CHKERRQ(ierr); 114f3161b27SJose E. Roman aa = aa1; 115f3161b27SJose E. Roman } 116f3161b27SJose E. Roman ierr = PetscMemcpy(&ja[pos],cols,ncols*sizeof(int));CHKERRQ(ierr); 117f3161b27SJose E. Roman ierr = PetscMemcpy(&aa[pos],values,ncols*sizeof(PetscScalar));CHKERRQ(ierr); 118f3161b27SJose E. Roman ierr = MatRestoreRow(pmat,i,&ncols,&cols,&values);CHKERRQ(ierr); 119f3161b27SJose E. Roman } 120f3161b27SJose E. Roman 121f3161b27SJose E. Roman /* csr info is for local matrix so initialize im[] locally */ 122f3161b27SJose E. Roman ierr = PetscMalloc(lsize*sizeof(int),&im);CHKERRQ(ierr); 123f3161b27SJose E. Roman ierr = PetscMemcpy(im,&maptmp[mapptr[rank]-1],lsize*sizeof(int));CHKERRQ(ierr); 124f3161b27SJose E. Roman 125f3161b27SJose E. Roman /* 1-based indexing */ 1262fa5cd67SKarl Rupp for (i=0; i<ia[lsize]-1; i++) ja[i] = ja[i]+1; 127f3161b27SJose E. Roman 128f3161b27SJose E. Roman /* Now copy csr matrix to parms_mat object */ 129f3161b27SJose E. Roman parms_MatSetValues(parms->A,(int)lsize,im,ia,ja,aa,INSERT); 130f3161b27SJose E. Roman 131f3161b27SJose E. Roman /* free memory */ 132f3161b27SJose E. Roman ierr = PetscFree(maptmp);CHKERRQ(ierr); 133f3161b27SJose E. Roman ierr = PetscFree(mapptr);CHKERRQ(ierr); 134f3161b27SJose E. Roman ierr = PetscFree(aa);CHKERRQ(ierr); 135f3161b27SJose E. Roman ierr = PetscFree(ja);CHKERRQ(ierr); 136f3161b27SJose E. Roman ierr = PetscFree(ia);CHKERRQ(ierr); 137f3161b27SJose E. Roman ierr = PetscFree(im);CHKERRQ(ierr); 138f3161b27SJose E. Roman 139f3161b27SJose E. Roman /* setup parms matrix */ 140f3161b27SJose E. Roman parms_MatSetup(parms->A); 141f3161b27SJose E. Roman 142f3161b27SJose E. Roman /* if created, destroy the previous pARMS pc */ 143f3161b27SJose E. Roman if (parms->pc) { 144f3161b27SJose E. Roman parms_PCFree(&parms->pc); 1450298fd71SBarry Smith parms->pc = NULL; 146f3161b27SJose E. Roman } 147f3161b27SJose E. Roman 148f3161b27SJose E. Roman /* Now create pARMS preconditioner object based on A */ 149f3161b27SJose E. Roman parms_PCCreate(&parms->pc,parms->A); 150f3161b27SJose E. Roman 151f3161b27SJose E. Roman /* Transfer options from PC to pARMS */ 152f3161b27SJose E. Roman switch (parms->global) { 153f3161b27SJose E. Roman case 0: parms_PCSetType(parms->pc, PCRAS); break; 154f3161b27SJose E. Roman case 1: parms_PCSetType(parms->pc, PCSCHUR); break; 155f3161b27SJose E. Roman case 2: parms_PCSetType(parms->pc, PCBJ); break; 156f3161b27SJose E. Roman } 157f3161b27SJose E. Roman switch (parms->local) { 158f3161b27SJose E. Roman case 0: parms_PCSetILUType(parms->pc, PCILU0); break; 159f3161b27SJose E. Roman case 1: parms_PCSetILUType(parms->pc, PCILUK); break; 160f3161b27SJose E. Roman case 2: parms_PCSetILUType(parms->pc, PCILUT); break; 161f3161b27SJose E. Roman case 3: parms_PCSetILUType(parms->pc, PCARMS); break; 162f3161b27SJose E. Roman } 163f3161b27SJose E. Roman parms_PCSetInnerEps(parms->pc, parms->solvetol); 164f3161b27SJose E. Roman parms_PCSetNlevels(parms->pc, parms->levels); 165f3161b27SJose E. Roman parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0); 166f3161b27SJose E. Roman parms_PCSetBsize(parms->pc, parms->blocksize); 167f3161b27SJose E. Roman parms_PCSetTolInd(parms->pc, parms->indtol); 168f3161b27SJose E. Roman parms_PCSetInnerKSize(parms->pc, parms->maxdim); 169f3161b27SJose E. Roman parms_PCSetInnerMaxits(parms->pc, parms->maxits); 170f3161b27SJose E. Roman for (i=0; i<8; i++) meth[i] = parms->meth[i] ? 1 : 0; 171f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[0], 1); 172f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[4], 0); 173f3161b27SJose E. Roman parms_PCSetFill(parms->pc, parms->lfil); 174f3161b27SJose E. Roman parms_PCSetTol(parms->pc, parms->droptol); 175f3161b27SJose E. Roman 176f3161b27SJose E. Roman parms_PCSetup(parms->pc); 177f3161b27SJose E. Roman 178f3161b27SJose E. Roman /* Allocate two auxiliary vector of length lsize */ 179f3161b27SJose E. Roman if (parms->lvec0) { ierr = PetscFree(parms->lvec0);CHKERRQ(ierr); } 180f3161b27SJose E. Roman ierr = PetscMalloc(lsize*sizeof(PetscScalar), &parms->lvec0);CHKERRQ(ierr); 181f3161b27SJose E. Roman if (parms->lvec1) { ierr = PetscFree(parms->lvec1);CHKERRQ(ierr); } 182f3161b27SJose E. Roman ierr = PetscMalloc(lsize*sizeof(PetscScalar), &parms->lvec1);CHKERRQ(ierr); 183f3161b27SJose E. Roman PetscFunctionReturn(0); 184f3161b27SJose E. Roman } 185f3161b27SJose E. Roman 186f3161b27SJose E. Roman #undef __FUNCT__ 187f3161b27SJose E. Roman #define __FUNCT__ "PCView_PARMS" 188f3161b27SJose E. Roman static PetscErrorCode PCView_PARMS(PC pc,PetscViewer viewer) 189f3161b27SJose E. Roman { 190f3161b27SJose E. Roman PetscErrorCode ierr; 191f3161b27SJose E. Roman PetscBool iascii; 192f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 193f3161b27SJose E. Roman char *str; 194f3161b27SJose E. Roman double fill_fact; 195f3161b27SJose E. Roman 196f3161b27SJose E. Roman PetscFunctionBegin; 197251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 198f3161b27SJose E. Roman if (iascii) { 199f3161b27SJose E. Roman parms_PCGetName(parms->pc,&str); 200f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," global preconditioner: %s\n",str);CHKERRQ(ierr); 201f3161b27SJose E. Roman parms_PCILUGetName(parms->pc,&str); 202f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," local preconditioner: %s\n",str);CHKERRQ(ierr); 203f3161b27SJose E. Roman parms_PCGetRatio(parms->pc,&fill_fact); 204f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," non-zero elements/original non-zero entries: %-4.2f\n",fill_fact);CHKERRQ(ierr); 205f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Tolerance for local solve: %g\n",parms->solvetol);CHKERRQ(ierr); 206f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Number of levels: %d\n",parms->levels);CHKERRQ(ierr); 207f3161b27SJose E. Roman if (parms->nonsymperm) { 208f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation\n");CHKERRQ(ierr); 209f3161b27SJose E. Roman } 210f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Block size: %d\n",parms->blocksize);CHKERRQ(ierr); 211f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Tolerance for independent sets: %g\n",parms->indtol);CHKERRQ(ierr); 212f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Inner Krylov dimension: %d\n",parms->maxdim);CHKERRQ(ierr); 213f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Maximum number of inner iterations: %d\n",parms->maxits);CHKERRQ(ierr); 214f3161b27SJose E. Roman if (parms->meth[0]) { 215f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for interlevel blocks\n");CHKERRQ(ierr); 216f3161b27SJose E. Roman } 217f3161b27SJose E. Roman if (parms->meth[1]) { 218f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using column permutation for interlevel blocks\n");CHKERRQ(ierr); 219f3161b27SJose E. Roman } 220f3161b27SJose E. Roman if (parms->meth[2]) { 221f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using row scaling for interlevel blocks\n");CHKERRQ(ierr); 222f3161b27SJose E. Roman } 223f3161b27SJose E. Roman if (parms->meth[3]) { 224f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using column scaling for interlevel blocks\n");CHKERRQ(ierr); 225f3161b27SJose E. Roman } 226f3161b27SJose E. Roman if (parms->meth[4]) { 227f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for last level blocks\n");CHKERRQ(ierr); 228f3161b27SJose E. Roman } 229f3161b27SJose E. Roman if (parms->meth[5]) { 230f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using column permutation for last level blocks\n");CHKERRQ(ierr); 231f3161b27SJose E. Roman } 232f3161b27SJose E. Roman if (parms->meth[6]) { 233f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using row scaling for last level blocks\n");CHKERRQ(ierr); 234f3161b27SJose E. Roman } 235f3161b27SJose E. Roman if (parms->meth[7]) { 236f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," Using column scaling for last level blocks\n");CHKERRQ(ierr); 237f3161b27SJose E. Roman } 238f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," amount of fill-in for ilut, iluk and arms: %d\n",parms->lfil[0]);CHKERRQ(ierr); 239f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," amount of fill-in for schur: %d\n",parms->lfil[4]);CHKERRQ(ierr); 240f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," amount of fill-in for ILUT L and U: %d\n",parms->lfil[5]);CHKERRQ(ierr); 241f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n",parms->droptol[0]);CHKERRQ(ierr); 242f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," drop tolerance for schur complement at each level: %g\n",parms->droptol[4]);CHKERRQ(ierr); 243f3161b27SJose E. Roman ierr = PetscViewerASCIIPrintf(viewer," drop tolerance for ILUT in last level schur complement: %g\n",parms->droptol[5]);CHKERRQ(ierr); 244f3161b27SJose E. Roman } 245f3161b27SJose E. Roman PetscFunctionReturn(0); 246f3161b27SJose E. Roman } 247f3161b27SJose E. Roman 248f3161b27SJose E. Roman #undef __FUNCT__ 249f3161b27SJose E. Roman #define __FUNCT__ "PCDestroy_PARMS" 250f3161b27SJose E. Roman static PetscErrorCode PCDestroy_PARMS(PC pc) 251f3161b27SJose E. Roman { 252f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 253f3161b27SJose E. Roman PetscErrorCode ierr; 254f3161b27SJose E. Roman 255f3161b27SJose E. Roman PetscFunctionBegin; 2562fa5cd67SKarl Rupp if (parms->map) parms_MapFree(&parms->map); 2572fa5cd67SKarl Rupp if (parms->A) parms_MatFree(&parms->A); 2582fa5cd67SKarl Rupp if (parms->pc) parms_PCFree(&parms->pc); 259f3161b27SJose E. Roman if (parms->lvec0) { 260f3161b27SJose E. Roman ierr = PetscFree(parms->lvec0);CHKERRQ(ierr); 261f3161b27SJose E. Roman } 262f3161b27SJose E. Roman if (parms->lvec1) { 263f3161b27SJose E. Roman ierr = PetscFree(parms->lvec1);CHKERRQ(ierr); 264f3161b27SJose E. Roman } 265f3161b27SJose E. Roman ierr = PetscFree(pc->data);CHKERRQ(ierr); 266f3161b27SJose E. Roman 267f3161b27SJose E. Roman ierr = PetscObjectChangeTypeName((PetscObject)pc,0);CHKERRQ(ierr); 2680298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetGlobal_C","",NULL);CHKERRQ(ierr); 2690298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetLocal_C","",NULL);CHKERRQ(ierr); 2700298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveTolerances_C","",NULL);CHKERRQ(ierr); 2710298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveRestart_C","",NULL);CHKERRQ(ierr); 2720298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetNonsymPerm_C","",NULL);CHKERRQ(ierr); 2730298fd71SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetFill_C","",NULL);CHKERRQ(ierr); 274f3161b27SJose E. Roman PetscFunctionReturn(0); 275f3161b27SJose E. Roman } 276f3161b27SJose E. Roman 277f3161b27SJose E. Roman #undef __FUNCT__ 278f3161b27SJose E. Roman #define __FUNCT__ "PCSetFromOptions_PARMS" 279f3161b27SJose E. Roman static PetscErrorCode PCSetFromOptions_PARMS(PC pc) 280f3161b27SJose E. Roman { 281f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 282f3161b27SJose E. Roman PetscBool flag; 283f3161b27SJose E. Roman PCPARMSGlobalType global; 284f3161b27SJose E. Roman PCPARMSLocalType local; 285f3161b27SJose E. Roman PetscErrorCode ierr; 286f3161b27SJose E. Roman 287f3161b27SJose E. Roman PetscFunctionBegin; 288f3161b27SJose E. Roman ierr = PetscOptionsHead("PARMS Options");CHKERRQ(ierr); 289f3161b27SJose E. Roman ierr = PetscOptionsEnum("-pc_parms_global","Global preconditioner","PCPARMSSetGlobal",PCPARMSGlobalTypes,(PetscEnum)parms->global,(PetscEnum*)&global,&flag);CHKERRQ(ierr); 290f3161b27SJose E. Roman if (flag) { ierr = PCPARMSSetGlobal(pc,global);CHKERRQ(ierr); } 291f3161b27SJose E. Roman ierr = PetscOptionsEnum("-pc_parms_local","Local preconditioner","PCPARMSSetLocal",PCPARMSLocalTypes,(PetscEnum)parms->local,(PetscEnum*)&local,&flag);CHKERRQ(ierr); 292f3161b27SJose E. Roman if (flag) { ierr = PCPARMSSetLocal(pc,local);CHKERRQ(ierr); } 293f3161b27SJose E. Roman ierr = PetscOptionsReal("-pc_parms_solve_tol","Tolerance for local solve","PCPARMSSetSolveTolerances",parms->solvetol,&parms->solvetol,&flag);CHKERRQ(ierr); 294f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,&flag);CHKERRQ(ierr); 295f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,&flag);CHKERRQ(ierr); 296f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,&flag);CHKERRQ(ierr); 297f3161b27SJose E. Roman ierr = PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,&flag);CHKERRQ(ierr); 298f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,&flag);CHKERRQ(ierr); 299f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,&flag);CHKERRQ(ierr); 300f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],&flag);CHKERRQ(ierr); 301f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],&flag);CHKERRQ(ierr); 302f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],&flag);CHKERRQ(ierr); 303f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],&flag);CHKERRQ(ierr); 304f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],&flag);CHKERRQ(ierr); 305f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],&flag);CHKERRQ(ierr); 306f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],&flag);CHKERRQ(ierr); 307f3161b27SJose E. Roman ierr = PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],&flag);CHKERRQ(ierr); 308f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_lfil_ilu_arms","amount of fill-in for ilut, iluk and arms","PCPARMSSetFill",parms->lfil[0],&parms->lfil[0],&flag);CHKERRQ(ierr); 309f3161b27SJose E. Roman if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0]; 310f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_lfil_schur","amount of fill-in for schur","PCPARMSSetFill",parms->lfil[4],&parms->lfil[4],&flag);CHKERRQ(ierr); 311f3161b27SJose E. Roman ierr = PetscOptionsInt("-pc_parms_lfil_ilut_L_U","amount of fill-in for ILUT L and U","PCPARMSSetFill",parms->lfil[5],&parms->lfil[5],&flag);CHKERRQ(ierr); 312f3161b27SJose E. Roman if (flag) parms->lfil[6] = parms->lfil[5]; 313f3161b27SJose E. Roman ierr = PetscOptionsReal("-pc_parms_droptol_factors","drop tolerance for L, U, L^{-1}F and EU^{-1}","None",parms->droptol[0],&parms->droptol[0],&flag);CHKERRQ(ierr); 314f3161b27SJose E. Roman ierr = PetscOptionsReal("-pc_parms_droptol_schur_compl","drop tolerance for schur complement at each level","None",parms->droptol[4],&parms->droptol[4],&flag);CHKERRQ(ierr); 315f3161b27SJose E. Roman if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0]; 316f3161b27SJose E. Roman ierr = PetscOptionsReal("-pc_parms_droptol_last_schur","drop tolerance for ILUT in last level schur complement","None",parms->droptol[5],&parms->droptol[5],&flag);CHKERRQ(ierr); 317f3161b27SJose E. Roman if (flag) parms->droptol[6] = parms->droptol[5]; 318f3161b27SJose E. Roman ierr = PetscOptionsTail();CHKERRQ(ierr); 319f3161b27SJose E. Roman PetscFunctionReturn(0); 320f3161b27SJose E. Roman } 321f3161b27SJose E. Roman 322f3161b27SJose E. Roman #undef __FUNCT__ 323f3161b27SJose E. Roman #define __FUNCT__ "PCApply_PARMS" 324f3161b27SJose E. Roman static PetscErrorCode PCApply_PARMS(PC pc,Vec b,Vec x) 325f3161b27SJose E. Roman { 326f3161b27SJose E. Roman PetscErrorCode ierr; 327f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 328f3161b27SJose E. Roman const PetscScalar *b1; 329f3161b27SJose E. Roman PetscScalar *x1; 330f3161b27SJose E. Roman 331f3161b27SJose E. Roman PetscFunctionBegin; 332f3161b27SJose E. Roman ierr = VecGetArrayRead(b,&b1);CHKERRQ(ierr); 333f3161b27SJose E. Roman ierr = VecGetArray(x,&x1);CHKERRQ(ierr); 334f3161b27SJose E. Roman parms_VecPermAux((PetscScalar*)b1,parms->lvec0,parms->map); 335f3161b27SJose E. Roman parms_PCApply(parms->pc,parms->lvec0,parms->lvec1); 336f3161b27SJose E. Roman parms_VecInvPermAux(parms->lvec1,x1,parms->map); 337f3161b27SJose E. Roman ierr = VecRestoreArrayRead(b,&b1);CHKERRQ(ierr); 338f3161b27SJose E. Roman ierr = VecRestoreArray(x,&x1);CHKERRQ(ierr); 339f3161b27SJose E. Roman PetscFunctionReturn(0); 340f3161b27SJose E. Roman } 341f3161b27SJose E. Roman 342f3161b27SJose E. Roman EXTERN_C_BEGIN 343f3161b27SJose E. Roman #undef __FUNCT__ 344f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetGlobal_PARMS" 345f3161b27SJose E. Roman PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type) 346f3161b27SJose E. Roman { 347f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 348f3161b27SJose E. Roman 349f3161b27SJose E. Roman PetscFunctionBegin; 350f3161b27SJose E. Roman if (type != parms->global) { 351f3161b27SJose E. Roman parms->global = type; 352f3161b27SJose E. Roman pc->setupcalled = 0; 353f3161b27SJose E. Roman } 354f3161b27SJose E. Roman PetscFunctionReturn(0); 355f3161b27SJose E. Roman } 356f3161b27SJose E. Roman EXTERN_C_END 357f3161b27SJose E. Roman 358f3161b27SJose E. Roman #undef __FUNCT__ 359f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetGlobal" 3605de0dacdSJose E. Roman /*@ 361f3161b27SJose E. Roman PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS. 362f3161b27SJose E. Roman 363f3161b27SJose E. Roman Collective on PC 364f3161b27SJose E. Roman 365f3161b27SJose E. Roman Input Parameters: 366f3161b27SJose E. Roman + pc - the preconditioner context 367f3161b27SJose E. Roman - type - the global preconditioner type, one of 368f3161b27SJose E. Roman .vb 369f3161b27SJose E. Roman PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 370f3161b27SJose E. Roman PC_PARMS_GLOBAL_SCHUR - Schur complement 371f3161b27SJose E. Roman PC_PARMS_GLOBAL_BJ - Block Jacobi 372f3161b27SJose E. Roman .ve 373f3161b27SJose E. Roman 374f3161b27SJose E. Roman Options Database Keys: 375f3161b27SJose E. Roman -pc_parms_global [ras,schur,bj] - Sets global preconditioner 376f3161b27SJose E. Roman 377f3161b27SJose E. Roman Level: intermediate 378f3161b27SJose E. Roman 379f3161b27SJose E. Roman Notes: 380f3161b27SJose E. Roman See the pARMS function parms_PCSetType for more information. 381f3161b27SJose E. Roman 382f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetLocal() 383f3161b27SJose E. Roman @*/ 384f3161b27SJose E. Roman PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type) 385f3161b27SJose E. Roman { 386f3161b27SJose E. Roman PetscErrorCode ierr; 387f3161b27SJose E. Roman 388f3161b27SJose E. Roman PetscFunctionBegin; 389f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 391f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));CHKERRQ(ierr); 392f3161b27SJose E. Roman PetscFunctionReturn(0); 393f3161b27SJose E. Roman } 394f3161b27SJose E. Roman 395f3161b27SJose E. Roman EXTERN_C_BEGIN 396f3161b27SJose E. Roman #undef __FUNCT__ 397f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetLocal_PARMS" 398f3161b27SJose E. Roman PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type) 399f3161b27SJose E. Roman { 400f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 401f3161b27SJose E. Roman 402f3161b27SJose E. Roman PetscFunctionBegin; 403f3161b27SJose E. Roman if (type != parms->local) { 404f3161b27SJose E. Roman parms->local = type; 405f3161b27SJose E. Roman pc->setupcalled = 0; 406f3161b27SJose E. Roman } 407f3161b27SJose E. Roman PetscFunctionReturn(0); 408f3161b27SJose E. Roman } 409f3161b27SJose E. Roman EXTERN_C_END 410f3161b27SJose E. Roman 411f3161b27SJose E. Roman #undef __FUNCT__ 412f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetLocal" 4135de0dacdSJose E. Roman /*@ 414f3161b27SJose E. Roman PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS. 415f3161b27SJose E. Roman 416f3161b27SJose E. Roman Collective on PC 417f3161b27SJose E. Roman 418f3161b27SJose E. Roman Input Parameters: 419f3161b27SJose E. Roman + pc - the preconditioner context 420f3161b27SJose E. Roman - type - the local preconditioner type, one of 421f3161b27SJose E. Roman .vb 422f3161b27SJose E. Roman PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 423f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 424f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUT - ILUT preconditioner 425f3161b27SJose E. Roman PC_PARMS_LOCAL_ARMS - ARMS preconditioner 426f3161b27SJose E. Roman .ve 427f3161b27SJose E. Roman 428f3161b27SJose E. Roman Options Database Keys: 429f3161b27SJose E. Roman -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 430f3161b27SJose E. Roman 431f3161b27SJose E. Roman Level: intermediate 432f3161b27SJose E. Roman 433f3161b27SJose E. Roman Notes: 434f3161b27SJose E. Roman For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 435f3161b27SJose E. Roman variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 436f3161b27SJose E. Roman 437f3161b27SJose E. Roman See the pARMS function parms_PCILUSetType for more information. 438f3161b27SJose E. Roman 439f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm() 440f3161b27SJose E. Roman 441f3161b27SJose E. Roman @*/ 442f3161b27SJose E. Roman PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type) 443f3161b27SJose E. Roman { 444f3161b27SJose E. Roman PetscErrorCode ierr; 445f3161b27SJose E. Roman 446f3161b27SJose E. Roman PetscFunctionBegin; 447f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 448f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 449f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));CHKERRQ(ierr); 450f3161b27SJose E. Roman PetscFunctionReturn(0); 451f3161b27SJose E. Roman } 452f3161b27SJose E. Roman 453f3161b27SJose E. Roman EXTERN_C_BEGIN 454f3161b27SJose E. Roman #undef __FUNCT__ 455f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveTolerances_PARMS" 456f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits) 457f3161b27SJose E. Roman { 458f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 459f3161b27SJose E. Roman 460f3161b27SJose E. Roman PetscFunctionBegin; 461f3161b27SJose E. Roman if (tol != parms->solvetol) { 462f3161b27SJose E. Roman parms->solvetol = tol; 463f3161b27SJose E. Roman pc->setupcalled = 0; 464f3161b27SJose E. Roman } 465f3161b27SJose E. Roman if (maxits != parms->maxits) { 466f3161b27SJose E. Roman parms->maxits = maxits; 467f3161b27SJose E. Roman pc->setupcalled = 0; 468f3161b27SJose E. Roman } 469f3161b27SJose E. Roman PetscFunctionReturn(0); 470f3161b27SJose E. Roman } 471f3161b27SJose E. Roman EXTERN_C_END 472f3161b27SJose E. Roman 473f3161b27SJose E. Roman #undef __FUNCT__ 474f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveTolerances" 4755de0dacdSJose E. Roman /*@ 476f3161b27SJose E. Roman PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 477f3161b27SJose E. Roman inner GMRES solver, when the Schur global preconditioner is used. 478f3161b27SJose E. Roman 479f3161b27SJose E. Roman Collective on PC 480f3161b27SJose E. Roman 481f3161b27SJose E. Roman Input Parameters: 482f3161b27SJose E. Roman + pc - the preconditioner context 483f3161b27SJose E. Roman . tol - the convergence tolerance 484f3161b27SJose E. Roman - maxits - the maximum number of iterations to use 485f3161b27SJose E. Roman 486f3161b27SJose E. Roman Options Database Keys: 487f3161b27SJose E. Roman + -pc_parms_solve_tol - set the tolerance for local solve 488f3161b27SJose E. Roman - -pc_parms_max_it - set the maximum number of inner iterations 489f3161b27SJose E. Roman 490f3161b27SJose E. Roman Level: intermediate 491f3161b27SJose E. Roman 492f3161b27SJose E. Roman Notes: 493f3161b27SJose E. Roman See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information. 494f3161b27SJose E. Roman 495f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetSolveRestart() 496f3161b27SJose E. Roman @*/ 497f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits) 498f3161b27SJose E. Roman { 499f3161b27SJose E. Roman PetscErrorCode ierr; 500f3161b27SJose E. Roman 501f3161b27SJose E. Roman PetscFunctionBegin; 502f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 503f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));CHKERRQ(ierr); 504f3161b27SJose E. Roman PetscFunctionReturn(0); 505f3161b27SJose E. Roman } 506f3161b27SJose E. Roman 507f3161b27SJose E. Roman EXTERN_C_BEGIN 508f3161b27SJose E. Roman #undef __FUNCT__ 509f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveRestart_PARMS" 510f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart) 511f3161b27SJose E. Roman { 512f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 513f3161b27SJose E. Roman 514f3161b27SJose E. Roman PetscFunctionBegin; 515f3161b27SJose E. Roman if (restart != parms->maxdim) { 516f3161b27SJose E. Roman parms->maxdim = restart; 517f3161b27SJose E. Roman pc->setupcalled = 0; 518f3161b27SJose E. Roman } 519f3161b27SJose E. Roman PetscFunctionReturn(0); 520f3161b27SJose E. Roman } 521f3161b27SJose E. Roman EXTERN_C_END 522f3161b27SJose E. Roman 523f3161b27SJose E. Roman #undef __FUNCT__ 524f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveRestart" 5255de0dacdSJose E. Roman /*@ 526f3161b27SJose E. Roman PCPARMSSetSolveRestart - Sets the number of iterations at which the 527f3161b27SJose E. Roman inner GMRES solver restarts. 528f3161b27SJose E. Roman 529f3161b27SJose E. Roman Collective on PC 530f3161b27SJose E. Roman 531f3161b27SJose E. Roman Input Parameters: 532f3161b27SJose E. Roman + pc - the preconditioner context 533f3161b27SJose E. Roman - restart - maximum dimension of the Krylov subspace 534f3161b27SJose E. Roman 535f3161b27SJose E. Roman Options Database Keys: 536f3161b27SJose E. Roman . -pc_parms_max_dim - sets the inner Krylov dimension 537f3161b27SJose E. Roman 538f3161b27SJose E. Roman Level: intermediate 539f3161b27SJose E. Roman 540f3161b27SJose E. Roman Notes: 541f3161b27SJose E. Roman See the pARMS function parms_PCSetInnerKSize for more information. 542f3161b27SJose E. Roman 543f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetSolveTolerances() 544f3161b27SJose E. Roman @*/ 545f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart) 546f3161b27SJose E. Roman { 547f3161b27SJose E. Roman PetscErrorCode ierr; 548f3161b27SJose E. Roman 549f3161b27SJose E. Roman PetscFunctionBegin; 550f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 551f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));CHKERRQ(ierr); 552f3161b27SJose E. Roman PetscFunctionReturn(0); 553f3161b27SJose E. Roman } 554f3161b27SJose E. Roman 555f3161b27SJose E. Roman EXTERN_C_BEGIN 556f3161b27SJose E. Roman #undef __FUNCT__ 557f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetNonsymPerm_PARMS" 558f3161b27SJose E. Roman PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym) 559f3161b27SJose E. Roman { 560f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 561f3161b27SJose E. Roman 562f3161b27SJose E. Roman PetscFunctionBegin; 5635de0dacdSJose E. Roman if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 564f3161b27SJose E. Roman parms->nonsymperm = nonsym; 565f3161b27SJose E. Roman pc->setupcalled = 0; 566f3161b27SJose E. Roman } 567f3161b27SJose E. Roman PetscFunctionReturn(0); 568f3161b27SJose E. Roman } 569f3161b27SJose E. Roman EXTERN_C_END 570f3161b27SJose E. Roman 571f3161b27SJose E. Roman #undef __FUNCT__ 572f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetNonsymPerm" 5735de0dacdSJose E. Roman /*@ 574f3161b27SJose E. Roman PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 575f3161b27SJose E. Roman symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 576f3161b27SJose E. Roman 577f3161b27SJose E. Roman Collective on PC 578f3161b27SJose E. Roman 579f3161b27SJose E. Roman Input Parameters: 580f3161b27SJose E. Roman + pc - the preconditioner context 581f3161b27SJose E. Roman - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used; 582f3161b27SJose E. Roman PETSC_FALSE indicates the symmetric ARMS is used 583f3161b27SJose E. Roman 584f3161b27SJose E. Roman Options Database Keys: 585f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 586f3161b27SJose E. Roman 587f3161b27SJose E. Roman Level: intermediate 588f3161b27SJose E. Roman 589f3161b27SJose E. Roman Notes: 590f3161b27SJose E. Roman See the pARMS function parms_PCSetPermType for more information. 591f3161b27SJose E. Roman 592f3161b27SJose E. Roman .seealso: PCPARMS 593f3161b27SJose E. Roman @*/ 594f3161b27SJose E. Roman PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym) 595f3161b27SJose E. Roman { 596f3161b27SJose E. Roman PetscErrorCode ierr; 597f3161b27SJose E. Roman 598f3161b27SJose E. Roman PetscFunctionBegin; 599f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 600f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));CHKERRQ(ierr); 601f3161b27SJose E. Roman PetscFunctionReturn(0); 602f3161b27SJose E. Roman } 603f3161b27SJose E. Roman 604f3161b27SJose E. Roman EXTERN_C_BEGIN 605f3161b27SJose E. Roman #undef __FUNCT__ 606f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetFill_PARMS" 607f3161b27SJose E. Roman PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 608f3161b27SJose E. Roman { 609f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 610f3161b27SJose E. Roman 611f3161b27SJose E. Roman PetscFunctionBegin; 612f3161b27SJose E. Roman if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 613f3161b27SJose E. Roman parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 614f3161b27SJose E. Roman pc->setupcalled = 0; 615f3161b27SJose E. Roman } 616f3161b27SJose E. Roman if (lfil1 != parms->lfil[4]) { 617f3161b27SJose E. Roman parms->lfil[4] = lfil1; 618f3161b27SJose E. Roman pc->setupcalled = 0; 619f3161b27SJose E. Roman } 620f3161b27SJose E. Roman if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 621f3161b27SJose E. Roman parms->lfil[5] = parms->lfil[6] = lfil2; 622f3161b27SJose E. Roman pc->setupcalled = 0; 623f3161b27SJose E. Roman } 624f3161b27SJose E. Roman PetscFunctionReturn(0); 625f3161b27SJose E. Roman } 626f3161b27SJose E. Roman EXTERN_C_END 627f3161b27SJose E. Roman 628f3161b27SJose E. Roman #undef __FUNCT__ 629f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetFill" 6305de0dacdSJose E. Roman /*@ 631f3161b27SJose E. Roman PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 632f3161b27SJose E. Roman Consider the original matrix A = [B F; E C] and the approximate version 633f3161b27SJose E. Roman M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 634f3161b27SJose E. Roman 635f3161b27SJose E. Roman Collective on PC 636f3161b27SJose E. Roman 637f3161b27SJose E. Roman Input Parameters: 638f3161b27SJose E. Roman + pc - the preconditioner context 639f3161b27SJose E. Roman . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 640f3161b27SJose E. Roman . fil1 - the level of fill-in kept in S 641f3161b27SJose E. Roman - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 642f3161b27SJose E. Roman 643f3161b27SJose E. Roman Options Database Keys: 644f3161b27SJose E. Roman + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 645f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 646f3161b27SJose E. Roman - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 647f3161b27SJose E. Roman 648f3161b27SJose E. Roman Level: intermediate 649f3161b27SJose E. Roman 650f3161b27SJose E. Roman Notes: 651f3161b27SJose E. Roman See the pARMS function parms_PCSetFill for more information. 652f3161b27SJose E. Roman 653f3161b27SJose E. Roman .seealso: PCPARMS 654f3161b27SJose E. Roman @*/ 655f3161b27SJose E. Roman PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 656f3161b27SJose E. Roman { 657f3161b27SJose E. Roman PetscErrorCode ierr; 658f3161b27SJose E. Roman 659f3161b27SJose E. Roman PetscFunctionBegin; 660f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 661f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));CHKERRQ(ierr); 662f3161b27SJose E. Roman PetscFunctionReturn(0); 663f3161b27SJose E. Roman } 664f3161b27SJose E. Roman 665f3161b27SJose E. Roman /*MC 666f3161b27SJose E. Roman PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 667f3161b27SJose E. Roman available in the package pARMS 668f3161b27SJose E. Roman 669f3161b27SJose E. Roman Options Database Keys: 670f3161b27SJose E. Roman + -pc_parms_global - one of ras, schur, bj 671f3161b27SJose E. Roman . -pc_parms_local - one of ilu0, iluk, ilut, arms 672f3161b27SJose E. Roman . -pc_parms_solve_tol - set the tolerance for local solve 673f3161b27SJose E. Roman . -pc_parms_levels - set the number of levels 674f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 675f3161b27SJose E. Roman . -pc_parms_blocksize - set the block size 676f3161b27SJose E. Roman . -pc_parms_ind_tol - set the tolerance for independent sets 677f3161b27SJose E. Roman . -pc_parms_max_dim - set the inner krylov dimension 678f3161b27SJose E. Roman . -pc_parms_max_it - set the maximum number of inner iterations 6794b62eef3SJose E. Roman . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 680f3161b27SJose E. Roman . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 681f3161b27SJose E. Roman . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 682f3161b27SJose E. Roman . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 6834b62eef3SJose E. Roman . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 684f3161b27SJose E. Roman . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 685f3161b27SJose E. Roman . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 686f3161b27SJose E. Roman . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 687f3161b27SJose E. Roman . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 688f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 689f3161b27SJose E. Roman . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 690f3161b27SJose E. Roman . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 691f3161b27SJose E. Roman . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 692f3161b27SJose E. Roman - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 693f3161b27SJose E. Roman 694f3161b27SJose E. Roman IMPORTANT: 695f3161b27SJose E. Roman Unless configured appropriately, this preconditioner performs an inexact solve 696f3161b27SJose E. Roman as part of the preconditioner application. Therefore, it must be used in combination 697f3161b27SJose E. Roman with flexible variants of iterative solvers, such as KSPFGMRES or KSPCGR. 698f3161b27SJose E. Roman 699f3161b27SJose E. Roman Level: intermediate 700f3161b27SJose E. Roman 701f3161b27SJose E. Roman .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 702f3161b27SJose E. Roman M*/ 703f3161b27SJose E. Roman 704f3161b27SJose E. Roman EXTERN_C_BEGIN 705f3161b27SJose E. Roman #undef __FUNCT__ 706f3161b27SJose E. Roman #define __FUNCT__ "PCCreate_PARMS" 707f3161b27SJose E. Roman PetscErrorCode PCCreate_PARMS(PC pc) 708f3161b27SJose E. Roman { 709f3161b27SJose E. Roman PC_PARMS *parms; 710f3161b27SJose E. Roman PetscErrorCode ierr; 711f3161b27SJose E. Roman 712f3161b27SJose E. Roman PetscFunctionBegin; 713f3161b27SJose E. Roman ierr = PetscNewLog(pc,PC_PARMS,&parms);CHKERRQ(ierr); 7142fa5cd67SKarl Rupp 715f3161b27SJose E. Roman parms->map = 0; 716f3161b27SJose E. Roman parms->A = 0; 717f3161b27SJose E. Roman parms->pc = 0; 718abd6c125SJose E. Roman parms->global = PC_PARMS_GLOBAL_RAS; 719f3161b27SJose E. Roman parms->local = PC_PARMS_LOCAL_ARMS; 720f3161b27SJose E. Roman parms->levels = 10; 721f3161b27SJose E. Roman parms->nonsymperm = PETSC_TRUE; 722f3161b27SJose E. Roman parms->blocksize = 250; 723abd6c125SJose E. Roman parms->maxdim = 0; 724abd6c125SJose E. Roman parms->maxits = 0; 7254b62eef3SJose E. Roman parms->meth[0] = PETSC_FALSE; 7264b62eef3SJose E. Roman parms->meth[1] = PETSC_FALSE; 7274b62eef3SJose E. Roman parms->meth[2] = PETSC_FALSE; 7284b62eef3SJose E. Roman parms->meth[3] = PETSC_FALSE; 7294b62eef3SJose E. Roman parms->meth[4] = PETSC_FALSE; 7304b62eef3SJose E. Roman parms->meth[5] = PETSC_FALSE; 7314b62eef3SJose E. Roman parms->meth[6] = PETSC_FALSE; 7324b62eef3SJose E. Roman parms->meth[7] = PETSC_FALSE; 733f3161b27SJose E. Roman parms->solvetol = 0.01; 734f3161b27SJose E. Roman parms->indtol = 0.4; 735f3161b27SJose E. Roman parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 736f3161b27SJose E. Roman parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 737f3161b27SJose E. Roman parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 738f3161b27SJose E. Roman parms->droptol[4] = 0.001; 739f3161b27SJose E. Roman parms->droptol[5] = parms->droptol[6] = 0.001; 7402fa5cd67SKarl Rupp 741f3161b27SJose E. Roman pc->data = parms; 742f3161b27SJose E. Roman pc->ops->destroy = PCDestroy_PARMS; 743f3161b27SJose E. Roman pc->ops->setfromoptions = PCSetFromOptions_PARMS; 744f3161b27SJose E. Roman pc->ops->setup = PCSetUp_PARMS; 745f3161b27SJose E. Roman pc->ops->apply = PCApply_PARMS; 746f3161b27SJose E. Roman pc->ops->view = PCView_PARMS; 7472fa5cd67SKarl Rupp 748f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetGlobal_C","PCPARMSSetGlobal_PARMS",PCPARMSSetGlobal_PARMS);CHKERRQ(ierr); 749f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetLocal_C","PCPARMSSetLocal_PARMS",PCPARMSSetLocal_PARMS);CHKERRQ(ierr); 750f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveTolerances_C","PCPARMSSetSolveTolerances_PARMS",PCPARMSSetSolveTolerances_PARMS);CHKERRQ(ierr); 751f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveRestart_C","PCPARMSSetSolveRestart_PARMS",PCPARMSSetSolveRestart_PARMS);CHKERRQ(ierr); 752f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetNonsymPerm_C","PCPARMSSetNonsymPerm_PARMS",PCPARMSSetNonsymPerm_PARMS);CHKERRQ(ierr); 753f3161b27SJose E. Roman ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetFill_C","PCPARMSSetFill_PARMS",PCPARMSSetFill_PARMS);CHKERRQ(ierr); 754f3161b27SJose E. Roman PetscFunctionReturn(0); 755f3161b27SJose E. Roman } 756f3161b27SJose E. Roman EXTERN_C_END 757