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 8af0996ceSBarry Smith #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 22aaa7dc30SBarry Smith #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 static PetscErrorCode PCSetUp_PARMS(PC pc) 40f3161b27SJose E. Roman { 41f3161b27SJose E. Roman Mat pmat; 42f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 43f3161b27SJose E. Roman const PetscInt *mapptr0; 44f3161b27SJose E. Roman PetscInt n, lsize, low, high, i, pos, ncols, length; 45f3161b27SJose E. Roman int *maptmp, *mapptr, *ia, *ja, *ja1, *im; 46f3161b27SJose E. Roman PetscScalar *aa, *aa1; 47f3161b27SJose E. Roman const PetscInt *cols; 48f3161b27SJose E. Roman PetscInt meth[8]; 49f3161b27SJose E. Roman const PetscScalar *values; 50f3161b27SJose E. Roman MatInfo matinfo; 51f3161b27SJose E. Roman PetscMPIInt rank, npro; 52f3161b27SJose E. Roman 53f3161b27SJose E. Roman PetscFunctionBegin; 54f3161b27SJose E. Roman /* Get preconditioner matrix from PETSc and setup pARMS structs */ 559566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc,NULL,&pmat)); 56ce94432eSBarry Smith MPI_Comm_size(PetscObjectComm((PetscObject)pmat),&npro); 57ce94432eSBarry Smith MPI_Comm_rank(PetscObjectComm((PetscObject)pmat),&rank); 58f3161b27SJose E. Roman 599566063dSJacob Faibussowitsch PetscCall(MatGetSize(pmat,&n,NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npro+1,&mapptr)); 619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&maptmp)); 629566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pmat,&mapptr0)); 63f3161b27SJose E. Roman low = mapptr0[rank]; 64f3161b27SJose E. Roman high = mapptr0[rank+1]; 65f3161b27SJose E. Roman lsize = high - low; 66f3161b27SJose E. Roman 672fa5cd67SKarl Rupp for (i=0; i<npro+1; i++) mapptr[i] = mapptr0[i]+1; 682fa5cd67SKarl Rupp for (i = 0; i<n; i++) maptmp[i] = i+1; 69f3161b27SJose E. Roman 70f3161b27SJose E. Roman /* if created, destroy the previous map */ 71f3161b27SJose E. Roman if (parms->map) { 72f3161b27SJose E. Roman parms_MapFree(&parms->map); 730298fd71SBarry Smith parms->map = NULL; 74f3161b27SJose E. Roman } 75f3161b27SJose E. Roman 76f3161b27SJose E. Roman /* create pARMS map object */ 77ce94432eSBarry Smith parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,PetscObjectComm((PetscObject)pmat),1,NONINTERLACED); 78f3161b27SJose E. Roman 79f3161b27SJose E. Roman /* if created, destroy the previous pARMS matrix */ 80f3161b27SJose E. Roman if (parms->A) { 81f3161b27SJose E. Roman parms_MatFree(&parms->A); 820298fd71SBarry Smith parms->A = NULL; 83f3161b27SJose E. Roman } 84f3161b27SJose E. Roman 85f3161b27SJose E. Roman /* create pARMS mat object */ 86f3161b27SJose E. Roman parms_MatCreate(&parms->A,parms->map); 87f3161b27SJose E. Roman 88f3161b27SJose E. Roman /* setup and copy csr data structure for pARMS */ 899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize+1,&ia)); 90f3161b27SJose E. Roman ia[0] = 1; 919566063dSJacob Faibussowitsch PetscCall(MatGetInfo(pmat,MAT_LOCAL,&matinfo)); 92f3161b27SJose E. Roman length = matinfo.nz_used; 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length,&ja)); 949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length,&aa)); 95f3161b27SJose E. Roman 96f3161b27SJose E. Roman for (i = low; i<high; i++) { 97f3161b27SJose E. Roman pos = ia[i-low]-1; 989566063dSJacob Faibussowitsch PetscCall(MatGetRow(pmat,i,&ncols,&cols,&values)); 99f3161b27SJose E. Roman ia[i-low+1] = ia[i-low] + ncols; 100f3161b27SJose E. Roman 101f3161b27SJose E. Roman if (ia[i-low+1] >= length) { 102f3161b27SJose E. Roman length += ncols; 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length,&ja1)); 1049566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ja1,ja,ia[i-low]-1)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 106f3161b27SJose E. Roman ja = ja1; 1079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length,&aa1)); 1089566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(aa1,aa,ia[i-low]-1)); 1099566063dSJacob Faibussowitsch PetscCall(PetscFree(aa)); 110f3161b27SJose E. Roman aa = aa1; 111f3161b27SJose E. Roman } 1129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&ja[pos],cols,ncols)); 1139566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&aa[pos],values,ncols)); 1149566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pmat,i,&ncols,&cols,&values)); 115f3161b27SJose E. Roman } 116f3161b27SJose E. Roman 117f3161b27SJose E. Roman /* csr info is for local matrix so initialize im[] locally */ 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize,&im)); 1199566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(im,&maptmp[mapptr[rank]-1],lsize)); 120f3161b27SJose E. Roman 121f3161b27SJose E. Roman /* 1-based indexing */ 1222fa5cd67SKarl Rupp for (i=0; i<ia[lsize]-1; i++) ja[i] = ja[i]+1; 123f3161b27SJose E. Roman 124f3161b27SJose E. Roman /* Now copy csr matrix to parms_mat object */ 125f3161b27SJose E. Roman parms_MatSetValues(parms->A,(int)lsize,im,ia,ja,aa,INSERT); 126f3161b27SJose E. Roman 127f3161b27SJose E. Roman /* free memory */ 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(maptmp)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(mapptr)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(aa)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(ia)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(im)); 134f3161b27SJose E. Roman 135f3161b27SJose E. Roman /* setup parms matrix */ 136f3161b27SJose E. Roman parms_MatSetup(parms->A); 137f3161b27SJose E. Roman 138f3161b27SJose E. Roman /* if created, destroy the previous pARMS pc */ 139f3161b27SJose E. Roman if (parms->pc) { 140f3161b27SJose E. Roman parms_PCFree(&parms->pc); 1410298fd71SBarry Smith parms->pc = NULL; 142f3161b27SJose E. Roman } 143f3161b27SJose E. Roman 144f3161b27SJose E. Roman /* Now create pARMS preconditioner object based on A */ 145f3161b27SJose E. Roman parms_PCCreate(&parms->pc,parms->A); 146f3161b27SJose E. Roman 147f3161b27SJose E. Roman /* Transfer options from PC to pARMS */ 148f3161b27SJose E. Roman switch (parms->global) { 149f3161b27SJose E. Roman case 0: parms_PCSetType(parms->pc, PCRAS); break; 150f3161b27SJose E. Roman case 1: parms_PCSetType(parms->pc, PCSCHUR); break; 151f3161b27SJose E. Roman case 2: parms_PCSetType(parms->pc, PCBJ); break; 152f3161b27SJose E. Roman } 153f3161b27SJose E. Roman switch (parms->local) { 154f3161b27SJose E. Roman case 0: parms_PCSetILUType(parms->pc, PCILU0); break; 155f3161b27SJose E. Roman case 1: parms_PCSetILUType(parms->pc, PCILUK); break; 156f3161b27SJose E. Roman case 2: parms_PCSetILUType(parms->pc, PCILUT); break; 157f3161b27SJose E. Roman case 3: parms_PCSetILUType(parms->pc, PCARMS); break; 158f3161b27SJose E. Roman } 159f3161b27SJose E. Roman parms_PCSetInnerEps(parms->pc, parms->solvetol); 160f3161b27SJose E. Roman parms_PCSetNlevels(parms->pc, parms->levels); 161f3161b27SJose E. Roman parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0); 162f3161b27SJose E. Roman parms_PCSetBsize(parms->pc, parms->blocksize); 163f3161b27SJose E. Roman parms_PCSetTolInd(parms->pc, parms->indtol); 164f3161b27SJose E. Roman parms_PCSetInnerKSize(parms->pc, parms->maxdim); 165f3161b27SJose E. Roman parms_PCSetInnerMaxits(parms->pc, parms->maxits); 166f3161b27SJose E. Roman for (i=0; i<8; i++) meth[i] = parms->meth[i] ? 1 : 0; 167f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[0], 1); 168f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[4], 0); 169f3161b27SJose E. Roman parms_PCSetFill(parms->pc, parms->lfil); 170f3161b27SJose E. Roman parms_PCSetTol(parms->pc, parms->droptol); 171f3161b27SJose E. Roman 172f3161b27SJose E. Roman parms_PCSetup(parms->pc); 173f3161b27SJose E. Roman 174f3161b27SJose E. Roman /* Allocate two auxiliary vector of length lsize */ 1759566063dSJacob Faibussowitsch if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 1769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec0)); 1779566063dSJacob Faibussowitsch if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 1789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec1)); 179f3161b27SJose E. Roman PetscFunctionReturn(0); 180f3161b27SJose E. Roman } 181f3161b27SJose E. Roman 182f3161b27SJose E. Roman static PetscErrorCode PCView_PARMS(PC pc,PetscViewer viewer) 183f3161b27SJose E. Roman { 184f3161b27SJose E. Roman PetscBool iascii; 185f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 186f3161b27SJose E. Roman char *str; 187f3161b27SJose E. Roman double fill_fact; 188f3161b27SJose E. Roman 189f3161b27SJose E. Roman PetscFunctionBegin; 1909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 191f3161b27SJose E. Roman if (iascii) { 192f3161b27SJose E. Roman parms_PCGetName(parms->pc,&str); 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," global preconditioner: %s\n",str)); 194f3161b27SJose E. Roman parms_PCILUGetName(parms->pc,&str); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," local preconditioner: %s\n",str)); 196f3161b27SJose E. Roman parms_PCGetRatio(parms->pc,&fill_fact); 1979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," non-zero elements/original non-zero entries: %-4.2f\n",fill_fact)); 1989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Tolerance for local solve: %g\n",parms->solvetol)); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Number of levels: %d\n",parms->levels)); 200f3161b27SJose E. Roman if (parms->nonsymperm) { 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation\n")); 202f3161b27SJose E. Roman } 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Block size: %d\n",parms->blocksize)); 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Tolerance for independent sets: %g\n",parms->indtol)); 2059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Inner Krylov dimension: %d\n",parms->maxdim)); 2069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Maximum number of inner iterations: %d\n",parms->maxits)); 207f3161b27SJose E. Roman if (parms->meth[0]) { 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for interlevel blocks\n")); 209f3161b27SJose E. Roman } 210f3161b27SJose E. Roman if (parms->meth[1]) { 2119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using column permutation for interlevel blocks\n")); 212f3161b27SJose E. Roman } 213f3161b27SJose E. Roman if (parms->meth[2]) { 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using row scaling for interlevel blocks\n")); 215f3161b27SJose E. Roman } 216f3161b27SJose E. Roman if (parms->meth[3]) { 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using column scaling for interlevel blocks\n")); 218f3161b27SJose E. Roman } 219f3161b27SJose E. Roman if (parms->meth[4]) { 2209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using nonsymmetric permutation for last level blocks\n")); 221f3161b27SJose E. Roman } 222f3161b27SJose E. Roman if (parms->meth[5]) { 2239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using column permutation for last level blocks\n")); 224f3161b27SJose E. Roman } 225f3161b27SJose E. Roman if (parms->meth[6]) { 2269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using row scaling for last level blocks\n")); 227f3161b27SJose E. Roman } 228f3161b27SJose E. Roman if (parms->meth[7]) { 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using column scaling for last level blocks\n")); 230f3161b27SJose E. Roman } 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," amount of fill-in for ilut, iluk and arms: %d\n",parms->lfil[0])); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," amount of fill-in for schur: %d\n",parms->lfil[4])); 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," amount of fill-in for ILUT L and U: %d\n",parms->lfil[5])); 2349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n",parms->droptol[0])); 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," drop tolerance for schur complement at each level: %g\n",parms->droptol[4])); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," drop tolerance for ILUT in last level schur complement: %g\n",parms->droptol[5])); 237f3161b27SJose E. Roman } 238f3161b27SJose E. Roman PetscFunctionReturn(0); 239f3161b27SJose E. Roman } 240f3161b27SJose E. Roman 241f3161b27SJose E. Roman static PetscErrorCode PCDestroy_PARMS(PC pc) 242f3161b27SJose E. Roman { 243f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 244f3161b27SJose E. Roman 245f3161b27SJose E. Roman PetscFunctionBegin; 2462fa5cd67SKarl Rupp if (parms->map) parms_MapFree(&parms->map); 2472fa5cd67SKarl Rupp if (parms->A) parms_MatFree(&parms->A); 2482fa5cd67SKarl Rupp if (parms->pc) parms_PCFree(&parms->pc); 249*1baa6e33SBarry Smith if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 250*1baa6e33SBarry Smith if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 2519566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 252f3161b27SJose E. Roman 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc,0)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",NULL)); 2589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",NULL)); 260f3161b27SJose E. Roman PetscFunctionReturn(0); 261f3161b27SJose E. Roman } 262f3161b27SJose E. Roman 2634416b707SBarry Smith static PetscErrorCode PCSetFromOptions_PARMS(PetscOptionItems *PetscOptionsObject,PC pc) 264f3161b27SJose E. Roman { 265f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 266f3161b27SJose E. Roman PetscBool flag; 267f3161b27SJose E. Roman PCPARMSGlobalType global; 268f3161b27SJose E. Roman PCPARMSLocalType local; 269f3161b27SJose E. Roman 270f3161b27SJose E. Roman PetscFunctionBegin; 271d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"PARMS Options"); 2729566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_global","Global preconditioner","PCPARMSSetGlobal",PCPARMSGlobalTypes,(PetscEnum)parms->global,(PetscEnum*)&global,&flag)); 2739566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetGlobal(pc,global)); 2749566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_local","Local preconditioner","PCPARMSSetLocal",PCPARMSLocalTypes,(PetscEnum)parms->local,(PetscEnum*)&local,&flag)); 2759566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetLocal(pc,local)); 2769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_solve_tol","Tolerance for local solve","PCPARMSSetSolveTolerances",parms->solvetol,&parms->solvetol,NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,NULL)); 2789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,NULL)); 2799566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,NULL)); 2809566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,NULL)); 2819566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,NULL)); 2829566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,NULL)); 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],NULL)); 2849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],NULL)); 2859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],NULL)); 2869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],NULL)); 2879566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],NULL)); 2889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],NULL)); 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],NULL)); 2909566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],NULL)); 2919566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_lfil_ilu_arms","amount of fill-in for ilut, iluk and arms","PCPARMSSetFill",parms->lfil[0],&parms->lfil[0],&flag)); 292f3161b27SJose E. Roman if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0]; 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_lfil_schur","amount of fill-in for schur","PCPARMSSetFill",parms->lfil[4],&parms->lfil[4],NULL)); 2949566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_lfil_ilut_L_U","amount of fill-in for ILUT L and U","PCPARMSSetFill",parms->lfil[5],&parms->lfil[5],&flag)); 295f3161b27SJose E. Roman if (flag) parms->lfil[6] = parms->lfil[5]; 2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_droptol_factors","drop tolerance for L, U, L^{-1}F and EU^{-1}","None",parms->droptol[0],&parms->droptol[0],NULL)); 2979566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_droptol_schur_compl","drop tolerance for schur complement at each level","None",parms->droptol[4],&parms->droptol[4],&flag)); 298f3161b27SJose E. Roman if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0]; 2999566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_droptol_last_schur","drop tolerance for ILUT in last level schur complement","None",parms->droptol[5],&parms->droptol[5],&flag)); 300f3161b27SJose E. Roman if (flag) parms->droptol[6] = parms->droptol[5]; 301d0609cedSBarry Smith PetscOptionsHeadEnd(); 302f3161b27SJose E. Roman PetscFunctionReturn(0); 303f3161b27SJose E. Roman } 304f3161b27SJose E. Roman 305f3161b27SJose E. Roman static PetscErrorCode PCApply_PARMS(PC pc,Vec b,Vec x) 306f3161b27SJose E. Roman { 307f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 308f3161b27SJose E. Roman const PetscScalar *b1; 309f3161b27SJose E. Roman PetscScalar *x1; 310f3161b27SJose E. Roman 311f3161b27SJose E. Roman PetscFunctionBegin; 3129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b,&b1)); 3139566063dSJacob Faibussowitsch PetscCall(VecGetArray(x,&x1)); 314f3161b27SJose E. Roman parms_VecPermAux((PetscScalar*)b1,parms->lvec0,parms->map); 315f3161b27SJose E. Roman parms_PCApply(parms->pc,parms->lvec0,parms->lvec1); 316f3161b27SJose E. Roman parms_VecInvPermAux(parms->lvec1,x1,parms->map); 3179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b,&b1)); 3189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x,&x1)); 319f3161b27SJose E. Roman PetscFunctionReturn(0); 320f3161b27SJose E. Roman } 321f3161b27SJose E. Roman 322f7a08781SBarry Smith static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type) 323f3161b27SJose E. Roman { 324f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 325f3161b27SJose E. Roman 326f3161b27SJose E. Roman PetscFunctionBegin; 327f3161b27SJose E. Roman if (type != parms->global) { 328f3161b27SJose E. Roman parms->global = type; 329f3161b27SJose E. Roman pc->setupcalled = 0; 330f3161b27SJose E. Roman } 331f3161b27SJose E. Roman PetscFunctionReturn(0); 332f3161b27SJose E. Roman } 333f3161b27SJose E. Roman 3345de0dacdSJose E. Roman /*@ 335f3161b27SJose E. Roman PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS. 336f3161b27SJose E. Roman 337f3161b27SJose E. Roman Collective on PC 338f3161b27SJose E. Roman 339f3161b27SJose E. Roman Input Parameters: 340f3161b27SJose E. Roman + pc - the preconditioner context 341f3161b27SJose E. Roman - type - the global preconditioner type, one of 342f3161b27SJose E. Roman .vb 343f3161b27SJose E. Roman PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 344f3161b27SJose E. Roman PC_PARMS_GLOBAL_SCHUR - Schur complement 345f3161b27SJose E. Roman PC_PARMS_GLOBAL_BJ - Block Jacobi 346f3161b27SJose E. Roman .ve 347f3161b27SJose E. Roman 348f3161b27SJose E. Roman Options Database Keys: 349f3161b27SJose E. Roman -pc_parms_global [ras,schur,bj] - Sets global preconditioner 350f3161b27SJose E. Roman 351f3161b27SJose E. Roman Level: intermediate 352f3161b27SJose E. Roman 353f3161b27SJose E. Roman Notes: 354f3161b27SJose E. Roman See the pARMS function parms_PCSetType for more information. 355f3161b27SJose E. Roman 356db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetLocal()` 357f3161b27SJose E. Roman @*/ 358f3161b27SJose E. Roman PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type) 359f3161b27SJose E. Roman { 360f3161b27SJose E. Roman PetscFunctionBegin; 361f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 362f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 363cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type)); 364f3161b27SJose E. Roman PetscFunctionReturn(0); 365f3161b27SJose E. Roman } 366f3161b27SJose E. Roman 367f7a08781SBarry Smith static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type) 368f3161b27SJose E. Roman { 369f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 370f3161b27SJose E. Roman 371f3161b27SJose E. Roman PetscFunctionBegin; 372f3161b27SJose E. Roman if (type != parms->local) { 373f3161b27SJose E. Roman parms->local = type; 374f3161b27SJose E. Roman pc->setupcalled = 0; 375f3161b27SJose E. Roman } 376f3161b27SJose E. Roman PetscFunctionReturn(0); 377f3161b27SJose E. Roman } 378f3161b27SJose E. Roman 3795de0dacdSJose E. Roman /*@ 380f3161b27SJose E. Roman PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS. 381f3161b27SJose E. Roman 382f3161b27SJose E. Roman Collective on PC 383f3161b27SJose E. Roman 384f3161b27SJose E. Roman Input Parameters: 385f3161b27SJose E. Roman + pc - the preconditioner context 386f3161b27SJose E. Roman - type - the local preconditioner type, one of 387f3161b27SJose E. Roman .vb 388f3161b27SJose E. Roman PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 389f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 390f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUT - ILUT preconditioner 391f3161b27SJose E. Roman PC_PARMS_LOCAL_ARMS - ARMS preconditioner 392f3161b27SJose E. Roman .ve 393f3161b27SJose E. Roman 394f3161b27SJose E. Roman Options Database Keys: 395f3161b27SJose E. Roman -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 396f3161b27SJose E. Roman 397f3161b27SJose E. Roman Level: intermediate 398f3161b27SJose E. Roman 399f3161b27SJose E. Roman Notes: 400f3161b27SJose E. Roman For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 401f3161b27SJose E. Roman variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 402f3161b27SJose E. Roman 403f3161b27SJose E. Roman See the pARMS function parms_PCILUSetType for more information. 404f3161b27SJose E. Roman 405db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetGlobal()`, `PCPARMSSetNonsymPerm()` 406f3161b27SJose E. Roman 407f3161b27SJose E. Roman @*/ 408f3161b27SJose E. Roman PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type) 409f3161b27SJose E. Roman { 410f3161b27SJose E. Roman PetscFunctionBegin; 411f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 412f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc,type,2); 413cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type)); 414f3161b27SJose E. Roman PetscFunctionReturn(0); 415f3161b27SJose E. Roman } 416f3161b27SJose E. Roman 417f7a08781SBarry Smith static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits) 418f3161b27SJose E. Roman { 419f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 420f3161b27SJose E. Roman 421f3161b27SJose E. Roman PetscFunctionBegin; 422f3161b27SJose E. Roman if (tol != parms->solvetol) { 423f3161b27SJose E. Roman parms->solvetol = tol; 424f3161b27SJose E. Roman pc->setupcalled = 0; 425f3161b27SJose E. Roman } 426f3161b27SJose E. Roman if (maxits != parms->maxits) { 427f3161b27SJose E. Roman parms->maxits = maxits; 428f3161b27SJose E. Roman pc->setupcalled = 0; 429f3161b27SJose E. Roman } 430f3161b27SJose E. Roman PetscFunctionReturn(0); 431f3161b27SJose E. Roman } 432f3161b27SJose E. Roman 4335de0dacdSJose E. Roman /*@ 434f3161b27SJose E. Roman PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 435f3161b27SJose E. Roman inner GMRES solver, when the Schur global preconditioner is used. 436f3161b27SJose E. Roman 437f3161b27SJose E. Roman Collective on PC 438f3161b27SJose E. Roman 439f3161b27SJose E. Roman Input Parameters: 440f3161b27SJose E. Roman + pc - the preconditioner context 441f3161b27SJose E. Roman . tol - the convergence tolerance 442f3161b27SJose E. Roman - maxits - the maximum number of iterations to use 443f3161b27SJose E. Roman 444f3161b27SJose E. Roman Options Database Keys: 445f3161b27SJose E. Roman + -pc_parms_solve_tol - set the tolerance for local solve 446f3161b27SJose E. Roman - -pc_parms_max_it - set the maximum number of inner iterations 447f3161b27SJose E. Roman 448f3161b27SJose E. Roman Level: intermediate 449f3161b27SJose E. Roman 450f3161b27SJose E. Roman Notes: 451f3161b27SJose E. Roman See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information. 452f3161b27SJose E. Roman 453db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetSolveRestart()` 454f3161b27SJose E. Roman @*/ 455f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits) 456f3161b27SJose E. Roman { 457f3161b27SJose E. Roman PetscFunctionBegin; 458f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 459cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits)); 460f3161b27SJose E. Roman PetscFunctionReturn(0); 461f3161b27SJose E. Roman } 462f3161b27SJose E. Roman 463f7a08781SBarry Smith static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart) 464f3161b27SJose E. Roman { 465f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 466f3161b27SJose E. Roman 467f3161b27SJose E. Roman PetscFunctionBegin; 468f3161b27SJose E. Roman if (restart != parms->maxdim) { 469f3161b27SJose E. Roman parms->maxdim = restart; 470f3161b27SJose E. Roman pc->setupcalled = 0; 471f3161b27SJose E. Roman } 472f3161b27SJose E. Roman PetscFunctionReturn(0); 473f3161b27SJose E. Roman } 474f3161b27SJose E. Roman 4755de0dacdSJose E. Roman /*@ 476f3161b27SJose E. Roman PCPARMSSetSolveRestart - Sets the number of iterations at which the 477f3161b27SJose E. Roman inner GMRES solver restarts. 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 - restart - maximum dimension of the Krylov subspace 484f3161b27SJose E. Roman 485f3161b27SJose E. Roman Options Database Keys: 486f3161b27SJose E. Roman . -pc_parms_max_dim - sets the inner Krylov dimension 487f3161b27SJose E. Roman 488f3161b27SJose E. Roman Level: intermediate 489f3161b27SJose E. Roman 490f3161b27SJose E. Roman Notes: 491f3161b27SJose E. Roman See the pARMS function parms_PCSetInnerKSize for more information. 492f3161b27SJose E. Roman 493db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetSolveTolerances()` 494f3161b27SJose E. Roman @*/ 495f3161b27SJose E. Roman PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart) 496f3161b27SJose E. Roman { 497f3161b27SJose E. Roman PetscFunctionBegin; 498f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 499cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart)); 500f3161b27SJose E. Roman PetscFunctionReturn(0); 501f3161b27SJose E. Roman } 502f3161b27SJose E. Roman 503f7a08781SBarry Smith static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym) 504f3161b27SJose E. Roman { 505f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 506f3161b27SJose E. Roman 507f3161b27SJose E. Roman PetscFunctionBegin; 5085de0dacdSJose E. Roman if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 509f3161b27SJose E. Roman parms->nonsymperm = nonsym; 510f3161b27SJose E. Roman pc->setupcalled = 0; 511f3161b27SJose E. Roman } 512f3161b27SJose E. Roman PetscFunctionReturn(0); 513f3161b27SJose E. Roman } 514f3161b27SJose E. Roman 5155de0dacdSJose E. Roman /*@ 516f3161b27SJose E. Roman PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 517f3161b27SJose E. Roman symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 518f3161b27SJose E. Roman 519f3161b27SJose E. Roman Collective on PC 520f3161b27SJose E. Roman 521f3161b27SJose E. Roman Input Parameters: 522f3161b27SJose E. Roman + pc - the preconditioner context 523f3161b27SJose E. Roman - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used; 524f3161b27SJose E. Roman PETSC_FALSE indicates the symmetric ARMS is used 525f3161b27SJose E. Roman 526f3161b27SJose E. Roman Options Database Keys: 527f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 528f3161b27SJose E. Roman 529f3161b27SJose E. Roman Level: intermediate 530f3161b27SJose E. Roman 531f3161b27SJose E. Roman Notes: 532f3161b27SJose E. Roman See the pARMS function parms_PCSetPermType for more information. 533f3161b27SJose E. Roman 534db781477SPatrick Sanan .seealso: `PCPARMS` 535f3161b27SJose E. Roman @*/ 536f3161b27SJose E. Roman PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym) 537f3161b27SJose E. Roman { 538f3161b27SJose E. Roman PetscFunctionBegin; 539f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 540cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym)); 541f3161b27SJose E. Roman PetscFunctionReturn(0); 542f3161b27SJose E. Roman } 543f3161b27SJose E. Roman 544f7a08781SBarry Smith static PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 545f3161b27SJose E. Roman { 546f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS*)pc->data; 547f3161b27SJose E. Roman 548f3161b27SJose E. Roman PetscFunctionBegin; 549f3161b27SJose E. Roman if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 550f3161b27SJose E. Roman parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 551f3161b27SJose E. Roman pc->setupcalled = 0; 552f3161b27SJose E. Roman } 553f3161b27SJose E. Roman if (lfil1 != parms->lfil[4]) { 554f3161b27SJose E. Roman parms->lfil[4] = lfil1; 555f3161b27SJose E. Roman pc->setupcalled = 0; 556f3161b27SJose E. Roman } 557f3161b27SJose E. Roman if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 558f3161b27SJose E. Roman parms->lfil[5] = parms->lfil[6] = lfil2; 559f3161b27SJose E. Roman pc->setupcalled = 0; 560f3161b27SJose E. Roman } 561f3161b27SJose E. Roman PetscFunctionReturn(0); 562f3161b27SJose E. Roman } 563f3161b27SJose E. Roman 5645de0dacdSJose E. Roman /*@ 565f3161b27SJose E. Roman PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 566f3161b27SJose E. Roman Consider the original matrix A = [B F; E C] and the approximate version 567f3161b27SJose E. Roman M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 568f3161b27SJose E. Roman 569f3161b27SJose E. Roman Collective on PC 570f3161b27SJose E. Roman 571f3161b27SJose E. Roman Input Parameters: 572f3161b27SJose E. Roman + pc - the preconditioner context 573f3161b27SJose E. Roman . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 574f3161b27SJose E. Roman . fil1 - the level of fill-in kept in S 575f3161b27SJose E. Roman - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 576f3161b27SJose E. Roman 577f3161b27SJose E. Roman Options Database Keys: 578f3161b27SJose E. Roman + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 579f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 580f3161b27SJose E. Roman - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 581f3161b27SJose E. Roman 582f3161b27SJose E. Roman Level: intermediate 583f3161b27SJose E. Roman 584f3161b27SJose E. Roman Notes: 585f3161b27SJose E. Roman See the pARMS function parms_PCSetFill for more information. 586f3161b27SJose E. Roman 587db781477SPatrick Sanan .seealso: `PCPARMS` 588f3161b27SJose E. Roman @*/ 589f3161b27SJose E. Roman PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 590f3161b27SJose E. Roman { 591f3161b27SJose E. Roman PetscFunctionBegin; 592f3161b27SJose E. Roman PetscValidHeaderSpecific(pc,PC_CLASSID,1); 593cac4c232SBarry Smith PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2)); 594f3161b27SJose E. Roman PetscFunctionReturn(0); 595f3161b27SJose E. Roman } 596f3161b27SJose E. Roman 597f3161b27SJose E. Roman /*MC 598f3161b27SJose E. Roman PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 599f3161b27SJose E. Roman available in the package pARMS 600f3161b27SJose E. Roman 601f3161b27SJose E. Roman Options Database Keys: 602f3161b27SJose E. Roman + -pc_parms_global - one of ras, schur, bj 603f3161b27SJose E. Roman . -pc_parms_local - one of ilu0, iluk, ilut, arms 604f3161b27SJose E. Roman . -pc_parms_solve_tol - set the tolerance for local solve 605f3161b27SJose E. Roman . -pc_parms_levels - set the number of levels 606f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 607f3161b27SJose E. Roman . -pc_parms_blocksize - set the block size 608f3161b27SJose E. Roman . -pc_parms_ind_tol - set the tolerance for independent sets 609f3161b27SJose E. Roman . -pc_parms_max_dim - set the inner krylov dimension 610f3161b27SJose E. Roman . -pc_parms_max_it - set the maximum number of inner iterations 6114b62eef3SJose E. Roman . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 612f3161b27SJose E. Roman . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 613f3161b27SJose E. Roman . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 614f3161b27SJose E. Roman . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 6154b62eef3SJose E. Roman . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 616f3161b27SJose E. Roman . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 617f3161b27SJose E. Roman . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 618f3161b27SJose E. Roman . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 619f3161b27SJose E. Roman . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 620f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 621f3161b27SJose E. Roman . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 622f3161b27SJose E. Roman . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 623f3161b27SJose E. Roman . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 624f3161b27SJose E. Roman - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 625f3161b27SJose E. Roman 626f3161b27SJose E. Roman IMPORTANT: 627f3161b27SJose E. Roman Unless configured appropriately, this preconditioner performs an inexact solve 628f3161b27SJose E. Roman as part of the preconditioner application. Therefore, it must be used in combination 629f6a5184aSRichard Tran Mills with flexible variants of iterative solvers, such as KSPFGMRES or KSPGCR. 630f3161b27SJose E. Roman 631f3161b27SJose E. Roman Level: intermediate 632f3161b27SJose E. Roman 633db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC` 634f3161b27SJose E. Roman M*/ 635f3161b27SJose E. Roman 6368cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc) 637f3161b27SJose E. Roman { 638f3161b27SJose E. Roman PC_PARMS *parms; 639f3161b27SJose E. Roman 640f3161b27SJose E. Roman PetscFunctionBegin; 6419566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&parms)); 6422fa5cd67SKarl Rupp 643f3161b27SJose E. Roman parms->map = 0; 644f3161b27SJose E. Roman parms->A = 0; 645f3161b27SJose E. Roman parms->pc = 0; 646abd6c125SJose E. Roman parms->global = PC_PARMS_GLOBAL_RAS; 647f3161b27SJose E. Roman parms->local = PC_PARMS_LOCAL_ARMS; 648f3161b27SJose E. Roman parms->levels = 10; 649f3161b27SJose E. Roman parms->nonsymperm = PETSC_TRUE; 650f3161b27SJose E. Roman parms->blocksize = 250; 651abd6c125SJose E. Roman parms->maxdim = 0; 652abd6c125SJose E. Roman parms->maxits = 0; 6534b62eef3SJose E. Roman parms->meth[0] = PETSC_FALSE; 6544b62eef3SJose E. Roman parms->meth[1] = PETSC_FALSE; 6554b62eef3SJose E. Roman parms->meth[2] = PETSC_FALSE; 6564b62eef3SJose E. Roman parms->meth[3] = PETSC_FALSE; 6574b62eef3SJose E. Roman parms->meth[4] = PETSC_FALSE; 6584b62eef3SJose E. Roman parms->meth[5] = PETSC_FALSE; 6594b62eef3SJose E. Roman parms->meth[6] = PETSC_FALSE; 6604b62eef3SJose E. Roman parms->meth[7] = PETSC_FALSE; 661f3161b27SJose E. Roman parms->solvetol = 0.01; 662f3161b27SJose E. Roman parms->indtol = 0.4; 663f3161b27SJose E. Roman parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 664f3161b27SJose E. Roman parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 665f3161b27SJose E. Roman parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 666f3161b27SJose E. Roman parms->droptol[4] = 0.001; 667f3161b27SJose E. Roman parms->droptol[5] = parms->droptol[6] = 0.001; 6682fa5cd67SKarl Rupp 669f3161b27SJose E. Roman pc->data = parms; 670f3161b27SJose E. Roman pc->ops->destroy = PCDestroy_PARMS; 671f3161b27SJose E. Roman pc->ops->setfromoptions = PCSetFromOptions_PARMS; 672f3161b27SJose E. Roman pc->ops->setup = PCSetUp_PARMS; 673f3161b27SJose E. Roman pc->ops->apply = PCApply_PARMS; 674f3161b27SJose E. Roman pc->ops->view = PCView_PARMS; 6752fa5cd67SKarl Rupp 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetGlobal_C",PCPARMSSetGlobal_PARMS)); 6779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetLocal_C",PCPARMSSetLocal_PARMS)); 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveTolerances_C",PCPARMSSetSolveTolerances_PARMS)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetSolveRestart_C",PCPARMSSetSolveRestart_PARMS)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetNonsymPerm_C",PCPARMSSetNonsymPerm_PARMS)); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPARMSSetFill_C",PCPARMSSetFill_PARMS)); 682f3161b27SJose E. Roman PetscFunctionReturn(0); 683f3161b27SJose E. Roman } 684