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); 60ce94432eSBarry Smith MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro); 61ce94432eSBarry Smith MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank); 62f3161b27SJose E. Roman 630298fd71SBarry Smith ierr = MatGetSize(pmat,&n,NULL);CHKERRQ(ierr); 64785e854fSJed Brown ierr = PetscMalloc1((npro+1),&mapptr);CHKERRQ(ierr); 65785e854fSJed Brown ierr = PetscMalloc1(n,&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 */ 81ce94432eSBarry 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 */ 93785e854fSJed Brown ierr = PetscMalloc1((lsize+1),&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; 97785e854fSJed Brown ierr = PetscMalloc1(length,&ja);CHKERRQ(ierr); 98785e854fSJed Brown ierr = PetscMalloc1(length,&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; 107785e854fSJed Brown ierr = PetscMalloc1(length,&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; 111785e854fSJed Brown ierr = PetscMalloc1(length,&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 */ 122785e854fSJed Brown ierr = PetscMalloc1(lsize,&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); } 180785e854fSJed Brown ierr = PetscMalloc1(lsize, &parms->lvec0);CHKERRQ(ierr); 181f3161b27SJose E. Roman if (parms->lvec1) { ierr = PetscFree(parms->lvec1);CHKERRQ(ierr); } 182785e854fSJed Brown ierr = PetscMalloc1(lsize, &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); 268bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",NULL);CHKERRQ(ierr); 269bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",NULL);CHKERRQ(ierr); 270bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",NULL);CHKERRQ(ierr); 271bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",NULL);CHKERRQ(ierr); 272bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",NULL);CHKERRQ(ierr); 273bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((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 #undef __FUNCT__ 343f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetGlobal_PARMS" 344f7a08781SBarry Smith static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type) 345f3161b27SJose E. Roman { 346f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 347f3161b27SJose E. Roman 348f3161b27SJose E. Roman PetscFunctionBegin; 349f3161b27SJose E. Roman if (type != parms->global) { 350f3161b27SJose E. Roman parms->global = type; 351f3161b27SJose E. Roman pc->setupcalled = 0; 352f3161b27SJose E. Roman } 353f3161b27SJose E. Roman PetscFunctionReturn(0); 354f3161b27SJose E. Roman } 355f3161b27SJose E. Roman 356f3161b27SJose E. Roman #undef __FUNCT__ 357f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetGlobal" 3585de0dacdSJose E. Roman /*@ 359f3161b27SJose E. Roman PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS. 360f3161b27SJose E. Roman 361f3161b27SJose E. Roman Collective on PC 362f3161b27SJose E. Roman 363f3161b27SJose E. Roman Input Parameters: 364f3161b27SJose E. Roman + pc - the preconditioner context 365f3161b27SJose E. Roman - type - the global preconditioner type, one of 366f3161b27SJose E. Roman .vb 367f3161b27SJose E. Roman PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 368f3161b27SJose E. Roman PC_PARMS_GLOBAL_SCHUR - Schur complement 369f3161b27SJose E. Roman PC_PARMS_GLOBAL_BJ - Block Jacobi 370f3161b27SJose E. Roman .ve 371f3161b27SJose E. Roman 372f3161b27SJose E. Roman Options Database Keys: 373f3161b27SJose E. Roman -pc_parms_global [ras,schur,bj] - Sets global preconditioner 374f3161b27SJose E. Roman 375f3161b27SJose E. Roman Level: intermediate 376f3161b27SJose E. Roman 377f3161b27SJose E. Roman Notes: 378f3161b27SJose E. Roman See the pARMS function parms_PCSetType for more information. 379f3161b27SJose E. Roman 380f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetLocal() 381f3161b27SJose E. Roman @*/ 382f3161b27SJose E. Roman PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type) 383f3161b27SJose E. Roman { 384f3161b27SJose E. Roman PetscErrorCode ierr; 385f3161b27SJose E. Roman 386f3161b27SJose E. Roman PetscFunctionBegin; 387f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 388f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 389f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));CHKERRQ(ierr); 390f3161b27SJose E. Roman PetscFunctionReturn(0); 391f3161b27SJose E. Roman } 392f3161b27SJose E. Roman 393f3161b27SJose E. Roman #undef __FUNCT__ 394f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetLocal_PARMS" 395f7a08781SBarry Smith static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type) 396f3161b27SJose E. Roman { 397f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 398f3161b27SJose E. Roman 399f3161b27SJose E. Roman PetscFunctionBegin; 400f3161b27SJose E. Roman if (type != parms->local) { 401f3161b27SJose E. Roman parms->local = type; 402f3161b27SJose E. Roman pc->setupcalled = 0; 403f3161b27SJose E. Roman } 404f3161b27SJose E. Roman PetscFunctionReturn(0); 405f3161b27SJose E. Roman } 406f3161b27SJose E. Roman 407f3161b27SJose E. Roman #undef __FUNCT__ 408f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetLocal" 4095de0dacdSJose E. Roman /*@ 410f3161b27SJose E. Roman PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS. 411f3161b27SJose E. Roman 412f3161b27SJose E. Roman Collective on PC 413f3161b27SJose E. Roman 414f3161b27SJose E. Roman Input Parameters: 415f3161b27SJose E. Roman + pc - the preconditioner context 416f3161b27SJose E. Roman - type - the local preconditioner type, one of 417f3161b27SJose E. Roman .vb 418f3161b27SJose E. Roman PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 419f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 420f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUT - ILUT preconditioner 421f3161b27SJose E. Roman PC_PARMS_LOCAL_ARMS - ARMS preconditioner 422f3161b27SJose E. Roman .ve 423f3161b27SJose E. Roman 424f3161b27SJose E. Roman Options Database Keys: 425f3161b27SJose E. Roman -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 426f3161b27SJose E. Roman 427f3161b27SJose E. Roman Level: intermediate 428f3161b27SJose E. Roman 429f3161b27SJose E. Roman Notes: 430f3161b27SJose E. Roman For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 431f3161b27SJose E. Roman variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 432f3161b27SJose E. Roman 433f3161b27SJose E. Roman See the pARMS function parms_PCILUSetType for more information. 434f3161b27SJose E. Roman 435f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm() 436f3161b27SJose E. Roman 437f3161b27SJose E. Roman @*/ 438f3161b27SJose E. Roman PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type) 439f3161b27SJose E. Roman { 440f3161b27SJose E. Roman PetscErrorCode ierr; 441f3161b27SJose E. Roman 442f3161b27SJose E. Roman PetscFunctionBegin; 443f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 444f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 445f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));CHKERRQ(ierr); 446f3161b27SJose E. Roman PetscFunctionReturn(0); 447f3161b27SJose E. Roman } 448f3161b27SJose E. Roman 449f3161b27SJose E. Roman #undef __FUNCT__ 450f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveTolerances_PARMS" 451f7a08781SBarry Smith static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits) 452f3161b27SJose E. Roman { 453f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 454f3161b27SJose E. Roman 455f3161b27SJose E. Roman PetscFunctionBegin; 456f3161b27SJose E. Roman if (tol != parms->solvetol) { 457f3161b27SJose E. Roman parms->solvetol = tol; 458f3161b27SJose E. Roman pc->setupcalled = 0; 459f3161b27SJose E. Roman } 460f3161b27SJose E. Roman if (maxits != parms->maxits) { 461f3161b27SJose E. Roman parms->maxits = maxits; 462f3161b27SJose E. Roman pc->setupcalled = 0; 463f3161b27SJose E. Roman } 464f3161b27SJose E. Roman PetscFunctionReturn(0); 465f3161b27SJose E. Roman } 466f3161b27SJose E. Roman 467f3161b27SJose E. Roman #undef __FUNCT__ 468f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveTolerances" 4695de0dacdSJose E. Roman /*@ 470f3161b27SJose E. Roman PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 471f3161b27SJose E. Roman inner GMRES solver, when the Schur global preconditioner is used. 472f3161b27SJose E. Roman 473f3161b27SJose E. Roman Collective on PC 474f3161b27SJose E. Roman 475f3161b27SJose E. Roman Input Parameters: 476f3161b27SJose E. Roman + pc - the preconditioner context 477f3161b27SJose E. Roman . tol - the convergence tolerance 478f3161b27SJose E. Roman - maxits - the maximum number of iterations to use 479f3161b27SJose E. Roman 480f3161b27SJose E. Roman Options Database Keys: 481f3161b27SJose E. Roman + -pc_parms_solve_tol - set the tolerance for local solve 482f3161b27SJose E. Roman - -pc_parms_max_it - set the maximum number of inner iterations 483f3161b27SJose E. Roman 484f3161b27SJose E. Roman Level: intermediate 485f3161b27SJose E. Roman 486f3161b27SJose E. Roman Notes: 487f3161b27SJose E. Roman See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information. 488f3161b27SJose E. Roman 489f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetSolveRestart() 490f3161b27SJose E. Roman @*/ 491f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits) 492f3161b27SJose E. Roman { 493f3161b27SJose E. Roman PetscErrorCode ierr; 494f3161b27SJose E. Roman 495f3161b27SJose E. Roman PetscFunctionBegin; 496f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 497f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));CHKERRQ(ierr); 498f3161b27SJose E. Roman PetscFunctionReturn(0); 499f3161b27SJose E. Roman } 500f3161b27SJose E. Roman 501f3161b27SJose E. Roman #undef __FUNCT__ 502f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveRestart_PARMS" 503f7a08781SBarry Smith static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart) 504f3161b27SJose E. Roman { 505f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 506f3161b27SJose E. Roman 507f3161b27SJose E. Roman PetscFunctionBegin; 508f3161b27SJose E. Roman if (restart != parms->maxdim) { 509f3161b27SJose E. Roman parms->maxdim = restart; 510f3161b27SJose E. Roman pc->setupcalled = 0; 511f3161b27SJose E. Roman } 512f3161b27SJose E. Roman PetscFunctionReturn(0); 513f3161b27SJose E. Roman } 514f3161b27SJose E. Roman 515f3161b27SJose E. Roman #undef __FUNCT__ 516f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetSolveRestart" 5175de0dacdSJose E. Roman /*@ 518f3161b27SJose E. Roman PCPARMSSetSolveRestart - Sets the number of iterations at which the 519f3161b27SJose E. Roman inner GMRES solver restarts. 520f3161b27SJose E. Roman 521f3161b27SJose E. Roman Collective on PC 522f3161b27SJose E. Roman 523f3161b27SJose E. Roman Input Parameters: 524f3161b27SJose E. Roman + pc - the preconditioner context 525f3161b27SJose E. Roman - restart - maximum dimension of the Krylov subspace 526f3161b27SJose E. Roman 527f3161b27SJose E. Roman Options Database Keys: 528f3161b27SJose E. Roman . -pc_parms_max_dim - sets the inner Krylov dimension 529f3161b27SJose E. Roman 530f3161b27SJose E. Roman Level: intermediate 531f3161b27SJose E. Roman 532f3161b27SJose E. Roman Notes: 533f3161b27SJose E. Roman See the pARMS function parms_PCSetInnerKSize for more information. 534f3161b27SJose E. Roman 535f3161b27SJose E. Roman .seealso: PCPARMS, PCPARMSSetSolveTolerances() 536f3161b27SJose E. Roman @*/ 537f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart) 538f3161b27SJose E. Roman { 539f3161b27SJose E. Roman PetscErrorCode ierr; 540f3161b27SJose E. Roman 541f3161b27SJose E. Roman PetscFunctionBegin; 542f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 543f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));CHKERRQ(ierr); 544f3161b27SJose E. Roman PetscFunctionReturn(0); 545f3161b27SJose E. Roman } 546f3161b27SJose E. Roman 547f3161b27SJose E. Roman #undef __FUNCT__ 548f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetNonsymPerm_PARMS" 549f7a08781SBarry Smith static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym) 550f3161b27SJose E. Roman { 551f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 552f3161b27SJose E. Roman 553f3161b27SJose E. Roman PetscFunctionBegin; 5545de0dacdSJose E. Roman if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 555f3161b27SJose E. Roman parms->nonsymperm = nonsym; 556f3161b27SJose E. Roman pc->setupcalled = 0; 557f3161b27SJose E. Roman } 558f3161b27SJose E. Roman PetscFunctionReturn(0); 559f3161b27SJose E. Roman } 560f3161b27SJose E. Roman 561f3161b27SJose E. Roman #undef __FUNCT__ 562f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetNonsymPerm" 5635de0dacdSJose E. Roman /*@ 564f3161b27SJose E. Roman PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 565f3161b27SJose E. Roman symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 566f3161b27SJose E. Roman 567f3161b27SJose E. Roman Collective on PC 568f3161b27SJose E. Roman 569f3161b27SJose E. Roman Input Parameters: 570f3161b27SJose E. Roman + pc - the preconditioner context 571f3161b27SJose E. Roman - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used; 572f3161b27SJose E. Roman PETSC_FALSE indicates the symmetric ARMS is used 573f3161b27SJose E. Roman 574f3161b27SJose E. Roman Options Database Keys: 575f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 576f3161b27SJose E. Roman 577f3161b27SJose E. Roman Level: intermediate 578f3161b27SJose E. Roman 579f3161b27SJose E. Roman Notes: 580f3161b27SJose E. Roman See the pARMS function parms_PCSetPermType for more information. 581f3161b27SJose E. Roman 582f3161b27SJose E. Roman .seealso: PCPARMS 583f3161b27SJose E. Roman @*/ 584f3161b27SJose E. Roman PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym) 585f3161b27SJose E. Roman { 586f3161b27SJose E. Roman PetscErrorCode ierr; 587f3161b27SJose E. Roman 588f3161b27SJose E. Roman PetscFunctionBegin; 589f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 590f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));CHKERRQ(ierr); 591f3161b27SJose E. Roman PetscFunctionReturn(0); 592f3161b27SJose E. Roman } 593f3161b27SJose E. Roman 594f3161b27SJose E. Roman #undef __FUNCT__ 595f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetFill_PARMS" 596f7a08781SBarry Smith static PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 597f3161b27SJose E. Roman { 598f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 599f3161b27SJose E. Roman 600f3161b27SJose E. Roman PetscFunctionBegin; 601f3161b27SJose E. Roman if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 602f3161b27SJose E. Roman parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 603f3161b27SJose E. Roman pc->setupcalled = 0; 604f3161b27SJose E. Roman } 605f3161b27SJose E. Roman if (lfil1 != parms->lfil[4]) { 606f3161b27SJose E. Roman parms->lfil[4] = lfil1; 607f3161b27SJose E. Roman pc->setupcalled = 0; 608f3161b27SJose E. Roman } 609f3161b27SJose E. Roman if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 610f3161b27SJose E. Roman parms->lfil[5] = parms->lfil[6] = lfil2; 611f3161b27SJose E. Roman pc->setupcalled = 0; 612f3161b27SJose E. Roman } 613f3161b27SJose E. Roman PetscFunctionReturn(0); 614f3161b27SJose E. Roman } 615f3161b27SJose E. Roman 616f3161b27SJose E. Roman #undef __FUNCT__ 617f3161b27SJose E. Roman #define __FUNCT__ "PCPARMSSetFill" 6185de0dacdSJose E. Roman /*@ 619f3161b27SJose E. Roman PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 620f3161b27SJose E. Roman Consider the original matrix A = [B F; E C] and the approximate version 621f3161b27SJose E. Roman M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 622f3161b27SJose E. Roman 623f3161b27SJose E. Roman Collective on PC 624f3161b27SJose E. Roman 625f3161b27SJose E. Roman Input Parameters: 626f3161b27SJose E. Roman + pc - the preconditioner context 627f3161b27SJose E. Roman . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 628f3161b27SJose E. Roman . fil1 - the level of fill-in kept in S 629f3161b27SJose E. Roman - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 630f3161b27SJose E. Roman 631f3161b27SJose E. Roman Options Database Keys: 632f3161b27SJose E. Roman + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 633f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 634f3161b27SJose E. Roman - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 635f3161b27SJose E. Roman 636f3161b27SJose E. Roman Level: intermediate 637f3161b27SJose E. Roman 638f3161b27SJose E. Roman Notes: 639f3161b27SJose E. Roman See the pARMS function parms_PCSetFill for more information. 640f3161b27SJose E. Roman 641f3161b27SJose E. Roman .seealso: PCPARMS 642f3161b27SJose E. Roman @*/ 643f3161b27SJose E. Roman PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 644f3161b27SJose E. Roman { 645f3161b27SJose E. Roman PetscErrorCode ierr; 646f3161b27SJose E. Roman 647f3161b27SJose E. Roman PetscFunctionBegin; 648f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 649f3161b27SJose E. Roman ierr = PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));CHKERRQ(ierr); 650f3161b27SJose E. Roman PetscFunctionReturn(0); 651f3161b27SJose E. Roman } 652f3161b27SJose E. Roman 653f3161b27SJose E. Roman /*MC 654f3161b27SJose E. Roman PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 655f3161b27SJose E. Roman available in the package pARMS 656f3161b27SJose E. Roman 657f3161b27SJose E. Roman Options Database Keys: 658f3161b27SJose E. Roman + -pc_parms_global - one of ras, schur, bj 659f3161b27SJose E. Roman . -pc_parms_local - one of ilu0, iluk, ilut, arms 660f3161b27SJose E. Roman . -pc_parms_solve_tol - set the tolerance for local solve 661f3161b27SJose E. Roman . -pc_parms_levels - set the number of levels 662f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 663f3161b27SJose E. Roman . -pc_parms_blocksize - set the block size 664f3161b27SJose E. Roman . -pc_parms_ind_tol - set the tolerance for independent sets 665f3161b27SJose E. Roman . -pc_parms_max_dim - set the inner krylov dimension 666f3161b27SJose E. Roman . -pc_parms_max_it - set the maximum number of inner iterations 6674b62eef3SJose E. Roman . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 668f3161b27SJose E. Roman . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 669f3161b27SJose E. Roman . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 670f3161b27SJose E. Roman . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 6714b62eef3SJose E. Roman . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 672f3161b27SJose E. Roman . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 673f3161b27SJose E. Roman . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 674f3161b27SJose E. Roman . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 675f3161b27SJose E. Roman . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 676f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 677f3161b27SJose E. Roman . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 678f3161b27SJose E. Roman . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 679f3161b27SJose E. Roman . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 680f3161b27SJose E. Roman - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 681f3161b27SJose E. Roman 682f3161b27SJose E. Roman IMPORTANT: 683f3161b27SJose E. Roman Unless configured appropriately, this preconditioner performs an inexact solve 684f3161b27SJose E. Roman as part of the preconditioner application. Therefore, it must be used in combination 685f3161b27SJose E. Roman with flexible variants of iterative solvers, such as KSPFGMRES or KSPCGR. 686f3161b27SJose E. Roman 687f3161b27SJose E. Roman Level: intermediate 688f3161b27SJose E. Roman 689f3161b27SJose E. Roman .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 690f3161b27SJose E. Roman M*/ 691f3161b27SJose E. Roman 692f3161b27SJose E. Roman #undef __FUNCT__ 693f3161b27SJose E. Roman #define __FUNCT__ "PCCreate_PARMS" 6948cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc) 695f3161b27SJose E. Roman { 696f3161b27SJose E. Roman PC_PARMS *parms; 697f3161b27SJose E. Roman PetscErrorCode ierr; 698f3161b27SJose E. Roman 699f3161b27SJose E. Roman PetscFunctionBegin; 700*b00a9115SJed Brown ierr = PetscNewLog(pc,&parms);CHKERRQ(ierr); 7012fa5cd67SKarl Rupp 702f3161b27SJose E. Roman parms->map = 0; 703f3161b27SJose E. Roman parms->A = 0; 704f3161b27SJose E. Roman parms->pc = 0; 705abd6c125SJose E. Roman parms->global = PC_PARMS_GLOBAL_RAS; 706f3161b27SJose E. Roman parms->local = PC_PARMS_LOCAL_ARMS; 707f3161b27SJose E. Roman parms->levels = 10; 708f3161b27SJose E. Roman parms->nonsymperm = PETSC_TRUE; 709f3161b27SJose E. Roman parms->blocksize = 250; 710abd6c125SJose E. Roman parms->maxdim = 0; 711abd6c125SJose E. Roman parms->maxits = 0; 7124b62eef3SJose E. Roman parms->meth[0] = PETSC_FALSE; 7134b62eef3SJose E. Roman parms->meth[1] = PETSC_FALSE; 7144b62eef3SJose E. Roman parms->meth[2] = PETSC_FALSE; 7154b62eef3SJose E. Roman parms->meth[3] = PETSC_FALSE; 7164b62eef3SJose E. Roman parms->meth[4] = PETSC_FALSE; 7174b62eef3SJose E. Roman parms->meth[5] = PETSC_FALSE; 7184b62eef3SJose E. Roman parms->meth[6] = PETSC_FALSE; 7194b62eef3SJose E. Roman parms->meth[7] = PETSC_FALSE; 720f3161b27SJose E. Roman parms->solvetol = 0.01; 721f3161b27SJose E. Roman parms->indtol = 0.4; 722f3161b27SJose E. Roman parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 723f3161b27SJose E. Roman parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 724f3161b27SJose E. Roman parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 725f3161b27SJose E. Roman parms->droptol[4] = 0.001; 726f3161b27SJose E. Roman parms->droptol[5] = parms->droptol[6] = 0.001; 7272fa5cd67SKarl Rupp 728f3161b27SJose E. Roman pc->data = parms; 729f3161b27SJose E. Roman pc->ops->destroy = PCDestroy_PARMS; 730f3161b27SJose E. Roman pc->ops->setfromoptions = PCSetFromOptions_PARMS; 731f3161b27SJose E. Roman pc->ops->setup = PCSetUp_PARMS; 732f3161b27SJose E. Roman pc->ops->apply = PCApply_PARMS; 733f3161b27SJose E. Roman pc->ops->view = PCView_PARMS; 7342fa5cd67SKarl Rupp 735bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",PCPARMSSetGlobal_PARMS);CHKERRQ(ierr); 736bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",PCPARMSSetLocal_PARMS);CHKERRQ(ierr); 737bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",PCPARMSSetSolveTolerances_PARMS);CHKERRQ(ierr); 738bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",PCPARMSSetSolveRestart_PARMS);CHKERRQ(ierr); 739bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",PCPARMSSetNonsymPerm_PARMS);CHKERRQ(ierr); 740bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",PCPARMSSetFill_PARMS);CHKERRQ(ierr); 741f3161b27SJose E. Roman PetscFunctionReturn(0); 742f3161b27SJose E. Roman } 743