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 39d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_PARMS(PC pc) 40d71ae5a4SJacob Faibussowitsch { 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)); 56f1580f4eSBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pmat), &npro)); 57f1580f4eSBarry Smith PetscCallMPI(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) { 149d71ae5a4SJacob Faibussowitsch case 0: 150d71ae5a4SJacob Faibussowitsch parms_PCSetType(parms->pc, PCRAS); 151d71ae5a4SJacob Faibussowitsch break; 152d71ae5a4SJacob Faibussowitsch case 1: 153d71ae5a4SJacob Faibussowitsch parms_PCSetType(parms->pc, PCSCHUR); 154d71ae5a4SJacob Faibussowitsch break; 155d71ae5a4SJacob Faibussowitsch case 2: 156d71ae5a4SJacob Faibussowitsch parms_PCSetType(parms->pc, PCBJ); 157d71ae5a4SJacob Faibussowitsch break; 158f3161b27SJose E. Roman } 159f3161b27SJose E. Roman switch (parms->local) { 160d71ae5a4SJacob Faibussowitsch case 0: 161d71ae5a4SJacob Faibussowitsch parms_PCSetILUType(parms->pc, PCILU0); 162d71ae5a4SJacob Faibussowitsch break; 163d71ae5a4SJacob Faibussowitsch case 1: 164d71ae5a4SJacob Faibussowitsch parms_PCSetILUType(parms->pc, PCILUK); 165d71ae5a4SJacob Faibussowitsch break; 166d71ae5a4SJacob Faibussowitsch case 2: 167d71ae5a4SJacob Faibussowitsch parms_PCSetILUType(parms->pc, PCILUT); 168d71ae5a4SJacob Faibussowitsch break; 169d71ae5a4SJacob Faibussowitsch case 3: 170d71ae5a4SJacob Faibussowitsch parms_PCSetILUType(parms->pc, PCARMS); 171d71ae5a4SJacob Faibussowitsch break; 172f3161b27SJose E. Roman } 173f3161b27SJose E. Roman parms_PCSetInnerEps(parms->pc, parms->solvetol); 174f3161b27SJose E. Roman parms_PCSetNlevels(parms->pc, parms->levels); 175f3161b27SJose E. Roman parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0); 176f3161b27SJose E. Roman parms_PCSetBsize(parms->pc, parms->blocksize); 177f3161b27SJose E. Roman parms_PCSetTolInd(parms->pc, parms->indtol); 178f3161b27SJose E. Roman parms_PCSetInnerKSize(parms->pc, parms->maxdim); 179f3161b27SJose E. Roman parms_PCSetInnerMaxits(parms->pc, parms->maxits); 180f3161b27SJose E. Roman for (i = 0; i < 8; i++) meth[i] = parms->meth[i] ? 1 : 0; 181f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[0], 1); 182f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[4], 0); 183f3161b27SJose E. Roman parms_PCSetFill(parms->pc, parms->lfil); 184f3161b27SJose E. Roman parms_PCSetTol(parms->pc, parms->droptol); 185f3161b27SJose E. Roman 186f3161b27SJose E. Roman parms_PCSetup(parms->pc); 187f3161b27SJose E. Roman 188f3161b27SJose E. Roman /* Allocate two auxiliary vector of length lsize */ 1899566063dSJacob Faibussowitsch if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 1909566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec0)); 1919566063dSJacob Faibussowitsch if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 1929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec1)); 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 194f3161b27SJose E. Roman } 195f3161b27SJose E. Roman 196d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_PARMS(PC pc, PetscViewer viewer) 197d71ae5a4SJacob Faibussowitsch { 198f3161b27SJose E. Roman PetscBool iascii; 199f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 200f3161b27SJose E. Roman char *str; 201f3161b27SJose E. Roman double fill_fact; 202f3161b27SJose E. Roman 203f3161b27SJose E. Roman PetscFunctionBegin; 2049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 205f3161b27SJose E. Roman if (iascii) { 206f3161b27SJose E. Roman parms_PCGetName(parms->pc, &str); 2079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " global preconditioner: %s\n", str)); 208f3161b27SJose E. Roman parms_PCILUGetName(parms->pc, &str); 2099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " local preconditioner: %s\n", str)); 210f3161b27SJose E. Roman parms_PCGetRatio(parms->pc, &fill_fact); 2119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " non-zero elements/original non-zero entries: %-4.2f\n", fill_fact)); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Tolerance for local solve: %g\n", parms->solvetol)); 2139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels: %d\n", parms->levels)); 21448a46eb9SPierre Jolivet if (parms->nonsymperm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation\n")); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Block size: %d\n", parms->blocksize)); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Tolerance for independent sets: %g\n", parms->indtol)); 2179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Inner Krylov dimension: %d\n", parms->maxdim)); 2189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum number of inner iterations: %d\n", parms->maxits)); 21948a46eb9SPierre Jolivet if (parms->meth[0]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation for interlevel blocks\n")); 22048a46eb9SPierre Jolivet if (parms->meth[1]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column permutation for interlevel blocks\n")); 22148a46eb9SPierre Jolivet if (parms->meth[2]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using row scaling for interlevel blocks\n")); 22248a46eb9SPierre Jolivet if (parms->meth[3]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column scaling for interlevel blocks\n")); 22348a46eb9SPierre Jolivet if (parms->meth[4]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation for last level blocks\n")); 22448a46eb9SPierre Jolivet if (parms->meth[5]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column permutation for last level blocks\n")); 22548a46eb9SPierre Jolivet if (parms->meth[6]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using row scaling for last level blocks\n")); 22648a46eb9SPierre Jolivet if (parms->meth[7]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column scaling for last level blocks\n")); 2279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for ilut, iluk and arms: %d\n", parms->lfil[0])); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for schur: %d\n", parms->lfil[4])); 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for ILUT L and U: %d\n", parms->lfil[5])); 2309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n", parms->droptol[0])); 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for schur complement at each level: %g\n", parms->droptol[4])); 2329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for ILUT in last level schur complement: %g\n", parms->droptol[5])); 233f3161b27SJose E. Roman } 2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 235f3161b27SJose E. Roman } 236f3161b27SJose E. Roman 237d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_PARMS(PC pc) 238d71ae5a4SJacob Faibussowitsch { 239f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 240f3161b27SJose E. Roman 241f3161b27SJose E. Roman PetscFunctionBegin; 2422fa5cd67SKarl Rupp if (parms->map) parms_MapFree(&parms->map); 2432fa5cd67SKarl Rupp if (parms->A) parms_MatFree(&parms->A); 2442fa5cd67SKarl Rupp if (parms->pc) parms_PCFree(&parms->pc); 2451baa6e33SBarry Smith if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 2461baa6e33SBarry Smith if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 2479566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 248f3161b27SJose E. Roman 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, 0)); 2509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetGlobal_C", NULL)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetLocal_C", NULL)); 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveTolerances_C", NULL)); 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveRestart_C", NULL)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetNonsymPerm_C", NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetFill_C", NULL)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257f3161b27SJose E. Roman } 258f3161b27SJose E. Roman 259d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_PARMS(PC pc, PetscOptionItems *PetscOptionsObject) 260d71ae5a4SJacob Faibussowitsch { 261f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 262f3161b27SJose E. Roman PetscBool flag; 263f3161b27SJose E. Roman PCPARMSGlobalType global; 264f3161b27SJose E. Roman PCPARMSLocalType local; 265f3161b27SJose E. Roman 266f3161b27SJose E. Roman PetscFunctionBegin; 267d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PARMS Options"); 2689566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_global", "Global preconditioner", "PCPARMSSetGlobal", PCPARMSGlobalTypes, (PetscEnum)parms->global, (PetscEnum *)&global, &flag)); 2699566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetGlobal(pc, global)); 2709566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_local", "Local preconditioner", "PCPARMSSetLocal", PCPARMSLocalTypes, (PetscEnum)parms->local, (PetscEnum *)&local, &flag)); 2719566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetLocal(pc, local)); 2729566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_solve_tol", "Tolerance for local solve", "PCPARMSSetSolveTolerances", parms->solvetol, &parms->solvetol, NULL)); 2739566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_levels", "Number of levels", "None", parms->levels, &parms->levels, NULL)); 2749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_nonsymmetric_perm", "Use nonsymmetric permutation", "PCPARMSSetNonsymPerm", parms->nonsymperm, &parms->nonsymperm, NULL)); 2759566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_blocksize", "Block size", "None", parms->blocksize, &parms->blocksize, NULL)); 2769566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_ind_tol", "Tolerance for independent sets", "None", parms->indtol, &parms->indtol, NULL)); 2779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_dim", "Inner Krylov dimension", "PCPARMSSetSolveRestart", parms->maxdim, &parms->maxdim, NULL)); 2789566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_it", "Maximum number of inner iterations", "PCPARMSSetSolveTolerances", parms->maxits, &parms->maxits, NULL)); 2799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm", "nonsymmetric permutation for interlevel blocks", "None", parms->meth[0], &parms->meth[0], NULL)); 2809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_perm", "column permutation for interlevel blocks", "None", parms->meth[1], &parms->meth[1], NULL)); 2819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_row_scaling", "row scaling for interlevel blocks", "None", parms->meth[2], &parms->meth[2], NULL)); 2829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_scaling", "column scaling for interlevel blocks", "None", parms->meth[3], &parms->meth[3], NULL)); 2839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_nonsymmetric_perm", "nonsymmetric permutation for last level blocks", "None", parms->meth[4], &parms->meth[4], NULL)); 2849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_perm", "column permutation for last level blocks", "None", parms->meth[5], &parms->meth[5], NULL)); 2859566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_row_scaling", "row scaling for last level blocks", "None", parms->meth[6], &parms->meth[6], NULL)); 2869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_scaling", "column scaling for last level blocks", "None", parms->meth[7], &parms->meth[7], NULL)); 2879566063dSJacob 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)); 288f3161b27SJose E. Roman if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0]; 2899566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_lfil_schur", "amount of fill-in for schur", "PCPARMSSetFill", parms->lfil[4], &parms->lfil[4], NULL)); 2909566063dSJacob 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)); 291f3161b27SJose E. Roman if (flag) parms->lfil[6] = parms->lfil[5]; 2929566063dSJacob 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)); 2939566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_droptol_schur_compl", "drop tolerance for schur complement at each level", "None", parms->droptol[4], &parms->droptol[4], &flag)); 294f3161b27SJose E. Roman if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0]; 2959566063dSJacob 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)); 296f3161b27SJose E. Roman if (flag) parms->droptol[6] = parms->droptol[5]; 297d0609cedSBarry Smith PetscOptionsHeadEnd(); 2983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 299f3161b27SJose E. Roman } 300f3161b27SJose E. Roman 301d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_PARMS(PC pc, Vec b, Vec x) 302d71ae5a4SJacob Faibussowitsch { 303f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 304f3161b27SJose E. Roman const PetscScalar *b1; 305f3161b27SJose E. Roman PetscScalar *x1; 306f3161b27SJose E. Roman 307f3161b27SJose E. Roman PetscFunctionBegin; 3089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &b1)); 3099566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x1)); 310f3161b27SJose E. Roman parms_VecPermAux((PetscScalar *)b1, parms->lvec0, parms->map); 311f3161b27SJose E. Roman parms_PCApply(parms->pc, parms->lvec0, parms->lvec1); 312f3161b27SJose E. Roman parms_VecInvPermAux(parms->lvec1, x1, parms->map); 3139566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &b1)); 3149566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x1)); 3153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 316f3161b27SJose E. Roman } 317f3161b27SJose E. Roman 318d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc, PCPARMSGlobalType type) 319d71ae5a4SJacob Faibussowitsch { 320f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 321f3161b27SJose E. Roman 322f3161b27SJose E. Roman PetscFunctionBegin; 323f3161b27SJose E. Roman if (type != parms->global) { 324f3161b27SJose E. Roman parms->global = type; 325f3161b27SJose E. Roman pc->setupcalled = 0; 326f3161b27SJose E. Roman } 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328f3161b27SJose E. Roman } 329f3161b27SJose E. Roman 3305de0dacdSJose E. Roman /*@ 331f1580f4eSBarry Smith PCPARMSSetGlobal - Sets the global preconditioner to be used in `PCPARMS`. 332f3161b27SJose E. Roman 333c3339decSBarry Smith Collective 334f3161b27SJose E. Roman 335f3161b27SJose E. Roman Input Parameters: 336f3161b27SJose E. Roman + pc - the preconditioner context 337f3161b27SJose E. Roman - type - the global preconditioner type, one of 338f3161b27SJose E. Roman .vb 339f3161b27SJose E. Roman PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 340f3161b27SJose E. Roman PC_PARMS_GLOBAL_SCHUR - Schur complement 341f3161b27SJose E. Roman PC_PARMS_GLOBAL_BJ - Block Jacobi 342f3161b27SJose E. Roman .ve 343f3161b27SJose E. Roman 344f1580f4eSBarry Smith Options Database Key: 345f1580f4eSBarry Smith . -pc_parms_global [ras,schur,bj] - Sets global preconditioner 346f3161b27SJose E. Roman 347f3161b27SJose E. Roman Level: intermediate 348f3161b27SJose E. Roman 349f1580f4eSBarry Smith Note: 350f1580f4eSBarry Smith See the pARMS function `parms_PCSetType()` for more information. 351f3161b27SJose E. Roman 352*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS`, `PCPARMSSetLocal()` 353f3161b27SJose E. Roman @*/ 354d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetGlobal(PC pc, PCPARMSGlobalType type) 355d71ae5a4SJacob Faibussowitsch { 356f3161b27SJose E. Roman PetscFunctionBegin; 357f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 358f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc, type, 2); 359cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetGlobal_C", (PC, PCPARMSGlobalType), (pc, type)); 3603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 361f3161b27SJose E. Roman } 362f3161b27SJose E. Roman 363d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc, PCPARMSLocalType type) 364d71ae5a4SJacob Faibussowitsch { 365f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 366f3161b27SJose E. Roman 367f3161b27SJose E. Roman PetscFunctionBegin; 368f3161b27SJose E. Roman if (type != parms->local) { 369f3161b27SJose E. Roman parms->local = type; 370f3161b27SJose E. Roman pc->setupcalled = 0; 371f3161b27SJose E. Roman } 3723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 373f3161b27SJose E. Roman } 374f3161b27SJose E. Roman 3755de0dacdSJose E. Roman /*@ 376f1580f4eSBarry Smith PCPARMSSetLocal - Sets the local preconditioner to be used in `PCPARMS`. 377f3161b27SJose E. Roman 378c3339decSBarry Smith Collective 379f3161b27SJose E. Roman 380f3161b27SJose E. Roman Input Parameters: 381f3161b27SJose E. Roman + pc - the preconditioner context 382f3161b27SJose E. Roman - type - the local preconditioner type, one of 383f3161b27SJose E. Roman .vb 384f3161b27SJose E. Roman PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 385f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 386f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUT - ILUT preconditioner 387f3161b27SJose E. Roman PC_PARMS_LOCAL_ARMS - ARMS preconditioner 388f3161b27SJose E. Roman .ve 389f3161b27SJose E. Roman 390f3161b27SJose E. Roman Options Database Keys: 391feefa0e1SJacob Faibussowitsch . pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 392f3161b27SJose E. Roman 393f3161b27SJose E. Roman Level: intermediate 394f3161b27SJose E. Roman 395f3161b27SJose E. Roman Notes: 396f3161b27SJose E. Roman For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 397f3161b27SJose E. Roman variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 398f3161b27SJose E. Roman 399f1580f4eSBarry Smith See the pARMS function `parms_PCILUSetType()` for more information. 400f3161b27SJose E. Roman 401*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS`, `PCPARMSSetGlobal()`, `PCPARMSSetNonsymPerm()` 402f3161b27SJose E. Roman 403f3161b27SJose E. Roman @*/ 404d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetLocal(PC pc, PCPARMSLocalType type) 405d71ae5a4SJacob Faibussowitsch { 406f3161b27SJose E. Roman PetscFunctionBegin; 407f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 408f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc, type, 2); 409cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetLocal_C", (PC, PCPARMSLocalType), (pc, type)); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 411f3161b27SJose E. Roman } 412f3161b27SJose E. Roman 413d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc, PetscReal tol, PetscInt maxits) 414d71ae5a4SJacob Faibussowitsch { 415f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 416f3161b27SJose E. Roman 417f3161b27SJose E. Roman PetscFunctionBegin; 418f3161b27SJose E. Roman if (tol != parms->solvetol) { 419f3161b27SJose E. Roman parms->solvetol = tol; 420f3161b27SJose E. Roman pc->setupcalled = 0; 421f3161b27SJose E. Roman } 422f3161b27SJose E. Roman if (maxits != parms->maxits) { 423f3161b27SJose E. Roman parms->maxits = maxits; 424f3161b27SJose E. Roman pc->setupcalled = 0; 425f3161b27SJose E. Roman } 4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 427f3161b27SJose E. Roman } 428f3161b27SJose E. Roman 4295de0dacdSJose E. Roman /*@ 430f3161b27SJose E. Roman PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 431f3161b27SJose E. Roman inner GMRES solver, when the Schur global preconditioner is used. 432f3161b27SJose E. Roman 433c3339decSBarry Smith Collective 434f3161b27SJose E. Roman 435f3161b27SJose E. Roman Input Parameters: 436f3161b27SJose E. Roman + pc - the preconditioner context 437f3161b27SJose E. Roman . tol - the convergence tolerance 438f3161b27SJose E. Roman - maxits - the maximum number of iterations to use 439f3161b27SJose E. Roman 440f3161b27SJose E. Roman Options Database Keys: 441f3161b27SJose E. Roman + -pc_parms_solve_tol - set the tolerance for local solve 442f3161b27SJose E. Roman - -pc_parms_max_it - set the maximum number of inner iterations 443f3161b27SJose E. Roman 444f3161b27SJose E. Roman Level: intermediate 445f3161b27SJose E. Roman 446f1580f4eSBarry Smith Note: 447f1580f4eSBarry Smith See the pARMS functions `parms_PCSetInnerEps()` and `parms_PCSetInnerMaxits()` for more information. 448f3161b27SJose E. Roman 449*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS`, `PCPARMSSetSolveRestart()` 450f3161b27SJose E. Roman @*/ 451d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetSolveTolerances(PC pc, PetscReal tol, PetscInt maxits) 452d71ae5a4SJacob Faibussowitsch { 453f3161b27SJose E. Roman PetscFunctionBegin; 454f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 455cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetSolveTolerances_C", (PC, PetscReal, PetscInt), (pc, tol, maxits)); 4563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 457f3161b27SJose E. Roman } 458f3161b27SJose E. Roman 459d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc, PetscInt restart) 460d71ae5a4SJacob Faibussowitsch { 461f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 462f3161b27SJose E. Roman 463f3161b27SJose E. Roman PetscFunctionBegin; 464f3161b27SJose E. Roman if (restart != parms->maxdim) { 465f3161b27SJose E. Roman parms->maxdim = restart; 466f3161b27SJose E. Roman pc->setupcalled = 0; 467f3161b27SJose E. Roman } 4683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 469f3161b27SJose E. Roman } 470f3161b27SJose E. Roman 4715de0dacdSJose E. Roman /*@ 472f3161b27SJose E. Roman PCPARMSSetSolveRestart - Sets the number of iterations at which the 473f3161b27SJose E. Roman inner GMRES solver restarts. 474f3161b27SJose E. Roman 475c3339decSBarry Smith Collective 476f3161b27SJose E. Roman 477f3161b27SJose E. Roman Input Parameters: 478f3161b27SJose E. Roman + pc - the preconditioner context 479f3161b27SJose E. Roman - restart - maximum dimension of the Krylov subspace 480f3161b27SJose E. Roman 481f1580f4eSBarry Smith Options Database Key: 482f3161b27SJose E. Roman . -pc_parms_max_dim - sets the inner Krylov dimension 483f3161b27SJose E. Roman 484f3161b27SJose E. Roman Level: intermediate 485f3161b27SJose E. Roman 486f1580f4eSBarry Smith Note: 487f3161b27SJose E. Roman See the pARMS function parms_PCSetInnerKSize for more information. 488f3161b27SJose E. Roman 489*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS`, `PCPARMSSetSolveTolerances()` 490f3161b27SJose E. Roman @*/ 491d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetSolveRestart(PC pc, PetscInt restart) 492d71ae5a4SJacob Faibussowitsch { 493f3161b27SJose E. Roman PetscFunctionBegin; 494f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 495cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetSolveRestart_C", (PC, PetscInt), (pc, restart)); 4963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 497f3161b27SJose E. Roman } 498f3161b27SJose E. Roman 499d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc, PetscBool nonsym) 500d71ae5a4SJacob Faibussowitsch { 501f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 502f3161b27SJose E. Roman 503f3161b27SJose E. Roman PetscFunctionBegin; 5045de0dacdSJose E. Roman if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 505f3161b27SJose E. Roman parms->nonsymperm = nonsym; 506f3161b27SJose E. Roman pc->setupcalled = 0; 507f3161b27SJose E. Roman } 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 509f3161b27SJose E. Roman } 510f3161b27SJose E. Roman 5115de0dacdSJose E. Roman /*@ 512f3161b27SJose E. Roman PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 513f3161b27SJose E. Roman symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 514f3161b27SJose E. Roman 515c3339decSBarry Smith Collective 516f3161b27SJose E. Roman 517f3161b27SJose E. Roman Input Parameters: 518f3161b27SJose E. Roman + pc - the preconditioner context 519f1580f4eSBarry Smith - nonsym - `PETSC_TRUE` indicates the non-symmetric ARMS is used; 520f1580f4eSBarry Smith `PETSC_FALSE` indicates the symmetric ARMS is used 521f3161b27SJose E. Roman 522f1580f4eSBarry Smith Options Database Key: 523f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 524f3161b27SJose E. Roman 525f3161b27SJose E. Roman Level: intermediate 526f3161b27SJose E. Roman 527f1580f4eSBarry Smith Note: 528f1580f4eSBarry Smith See the pARMS function `parms_PCSetPermType()` for more information. 529f3161b27SJose E. Roman 530*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS` 531f3161b27SJose E. Roman @*/ 532d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetNonsymPerm(PC pc, PetscBool nonsym) 533d71ae5a4SJacob Faibussowitsch { 534f3161b27SJose E. Roman PetscFunctionBegin; 535f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 536cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetNonsymPerm_C", (PC, PetscBool), (pc, nonsym)); 5373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 538f3161b27SJose E. Roman } 539f3161b27SJose E. Roman 540d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCPARMSSetFill_PARMS(PC pc, PetscInt lfil0, PetscInt lfil1, PetscInt lfil2) 541d71ae5a4SJacob Faibussowitsch { 542f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 543f3161b27SJose E. Roman 544f3161b27SJose E. Roman PetscFunctionBegin; 545f3161b27SJose E. Roman if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 546f3161b27SJose E. Roman parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 547f3161b27SJose E. Roman pc->setupcalled = 0; 548f3161b27SJose E. Roman } 549f3161b27SJose E. Roman if (lfil1 != parms->lfil[4]) { 550f3161b27SJose E. Roman parms->lfil[4] = lfil1; 551f3161b27SJose E. Roman pc->setupcalled = 0; 552f3161b27SJose E. Roman } 553f3161b27SJose E. Roman if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 554f3161b27SJose E. Roman parms->lfil[5] = parms->lfil[6] = lfil2; 555f3161b27SJose E. Roman pc->setupcalled = 0; 556f3161b27SJose E. Roman } 5573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 558f3161b27SJose E. Roman } 559f3161b27SJose E. Roman 5605de0dacdSJose E. Roman /*@ 561f3161b27SJose E. Roman PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 562f3161b27SJose E. Roman Consider the original matrix A = [B F; E C] and the approximate version 563f3161b27SJose E. Roman M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 564f3161b27SJose E. Roman 565c3339decSBarry Smith Collective 566f3161b27SJose E. Roman 567f3161b27SJose E. Roman Input Parameters: 568f3161b27SJose E. Roman + pc - the preconditioner context 569feefa0e1SJacob Faibussowitsch . lfil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 570feefa0e1SJacob Faibussowitsch . lfil1 - the level of fill-in kept in S 571feefa0e1SJacob Faibussowitsch - lfil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 572f3161b27SJose E. Roman 573f3161b27SJose E. Roman Options Database Keys: 574f3161b27SJose E. Roman + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 575f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 576f3161b27SJose E. Roman - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 577f3161b27SJose E. Roman 578f3161b27SJose E. Roman Level: intermediate 579f3161b27SJose E. Roman 580f1580f4eSBarry Smith Note: 581f1580f4eSBarry Smith See the pARMS function `parms_PCSetFill()` for more information. 582f3161b27SJose E. Roman 583*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCPARMS` 584f3161b27SJose E. Roman @*/ 585d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPARMSSetFill(PC pc, PetscInt lfil0, PetscInt lfil1, PetscInt lfil2) 586d71ae5a4SJacob Faibussowitsch { 587f3161b27SJose E. Roman PetscFunctionBegin; 588f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 589cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetFill_C", (PC, PetscInt, PetscInt, PetscInt), (pc, lfil0, lfil1, lfil2)); 5903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 591f3161b27SJose E. Roman } 592f3161b27SJose E. Roman 593f3161b27SJose E. Roman /*MC 594f3161b27SJose E. Roman PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 595f3161b27SJose E. Roman available in the package pARMS 596f3161b27SJose E. Roman 597f3161b27SJose E. Roman Options Database Keys: 598f3161b27SJose E. Roman + -pc_parms_global - one of ras, schur, bj 599f3161b27SJose E. Roman . -pc_parms_local - one of ilu0, iluk, ilut, arms 600f3161b27SJose E. Roman . -pc_parms_solve_tol - set the tolerance for local solve 601f3161b27SJose E. Roman . -pc_parms_levels - set the number of levels 602f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 603f3161b27SJose E. Roman . -pc_parms_blocksize - set the block size 604f3161b27SJose E. Roman . -pc_parms_ind_tol - set the tolerance for independent sets 605f3161b27SJose E. Roman . -pc_parms_max_dim - set the inner krylov dimension 606f3161b27SJose E. Roman . -pc_parms_max_it - set the maximum number of inner iterations 6074b62eef3SJose E. Roman . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 608f3161b27SJose E. Roman . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 609f3161b27SJose E. Roman . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 610f3161b27SJose E. Roman . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 6114b62eef3SJose E. Roman . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 612f3161b27SJose E. Roman . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 613f3161b27SJose E. Roman . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 614f3161b27SJose E. Roman . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 615f3161b27SJose E. Roman . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 616f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 617f3161b27SJose E. Roman . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 618f3161b27SJose E. Roman . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 619f3161b27SJose E. Roman . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 620f3161b27SJose E. Roman - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 621f3161b27SJose E. Roman 622f1580f4eSBarry Smith Note: 623f3161b27SJose E. Roman Unless configured appropriately, this preconditioner performs an inexact solve 624f3161b27SJose E. Roman as part of the preconditioner application. Therefore, it must be used in combination 625f1580f4eSBarry Smith with flexible variants of iterative solvers, such as `KSPFGMRES` or `KSPGCR`. 626f3161b27SJose E. Roman 627f3161b27SJose E. Roman Level: intermediate 628f3161b27SJose E. Roman 629*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMG`, `PCGAMG`, `PCHYPRE`, `PCPARMSSetGlobal()`, 630f1580f4eSBarry Smith `PCPARMSSetLocal()`, `PCPARMSSetSolveTolerances()`, `PCPARMSSetSolveRestart()`, `PCPARMSSetNonsymPerm()`, 631f1580f4eSBarry Smith `PCPARMSSetFill()` 632f3161b27SJose E. Roman M*/ 633f3161b27SJose E. Roman 634d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc) 635d71ae5a4SJacob Faibussowitsch { 636f3161b27SJose E. Roman PC_PARMS *parms; 637f3161b27SJose E. Roman 638f3161b27SJose E. Roman PetscFunctionBegin; 6394dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&parms)); 6402fa5cd67SKarl Rupp 641f3161b27SJose E. Roman parms->map = 0; 642f3161b27SJose E. Roman parms->A = 0; 643f3161b27SJose E. Roman parms->pc = 0; 644abd6c125SJose E. Roman parms->global = PC_PARMS_GLOBAL_RAS; 645f3161b27SJose E. Roman parms->local = PC_PARMS_LOCAL_ARMS; 646f3161b27SJose E. Roman parms->levels = 10; 647f3161b27SJose E. Roman parms->nonsymperm = PETSC_TRUE; 648f3161b27SJose E. Roman parms->blocksize = 250; 649abd6c125SJose E. Roman parms->maxdim = 0; 650abd6c125SJose E. Roman parms->maxits = 0; 6514b62eef3SJose E. Roman parms->meth[0] = PETSC_FALSE; 6524b62eef3SJose E. Roman parms->meth[1] = PETSC_FALSE; 6534b62eef3SJose E. Roman parms->meth[2] = PETSC_FALSE; 6544b62eef3SJose E. Roman parms->meth[3] = PETSC_FALSE; 6554b62eef3SJose E. Roman parms->meth[4] = PETSC_FALSE; 6564b62eef3SJose E. Roman parms->meth[5] = PETSC_FALSE; 6574b62eef3SJose E. Roman parms->meth[6] = PETSC_FALSE; 6584b62eef3SJose E. Roman parms->meth[7] = PETSC_FALSE; 659f3161b27SJose E. Roman parms->solvetol = 0.01; 660f3161b27SJose E. Roman parms->indtol = 0.4; 661f3161b27SJose E. Roman parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 662f3161b27SJose E. Roman parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 663f3161b27SJose E. Roman parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 664f3161b27SJose E. Roman parms->droptol[4] = 0.001; 665f3161b27SJose E. Roman parms->droptol[5] = parms->droptol[6] = 0.001; 6662fa5cd67SKarl Rupp 667f3161b27SJose E. Roman pc->data = parms; 668f3161b27SJose E. Roman pc->ops->destroy = PCDestroy_PARMS; 669f3161b27SJose E. Roman pc->ops->setfromoptions = PCSetFromOptions_PARMS; 670f3161b27SJose E. Roman pc->ops->setup = PCSetUp_PARMS; 671f3161b27SJose E. Roman pc->ops->apply = PCApply_PARMS; 672f3161b27SJose E. Roman pc->ops->view = PCView_PARMS; 6732fa5cd67SKarl Rupp 6749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetGlobal_C", PCPARMSSetGlobal_PARMS)); 6759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetLocal_C", PCPARMSSetLocal_PARMS)); 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveTolerances_C", PCPARMSSetSolveTolerances_PARMS)); 6779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveRestart_C", PCPARMSSetSolveRestart_PARMS)); 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetNonsymPerm_C", PCPARMSSetNonsymPerm_PARMS)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetFill_C", PCPARMSSetFill_PARMS)); 6803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 681f3161b27SJose E. Roman } 682