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 399371c9d4SSatish Balay static PetscErrorCode PCSetUp_PARMS(PC pc) { 40f3161b27SJose E. Roman Mat pmat; 41f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 42f3161b27SJose E. Roman const PetscInt *mapptr0; 43f3161b27SJose E. Roman PetscInt n, lsize, low, high, i, pos, ncols, length; 44f3161b27SJose E. Roman int *maptmp, *mapptr, *ia, *ja, *ja1, *im; 45f3161b27SJose E. Roman PetscScalar *aa, *aa1; 46f3161b27SJose E. Roman const PetscInt *cols; 47f3161b27SJose E. Roman PetscInt meth[8]; 48f3161b27SJose E. Roman const PetscScalar *values; 49f3161b27SJose E. Roman MatInfo matinfo; 50f3161b27SJose E. Roman PetscMPIInt rank, npro; 51f3161b27SJose E. Roman 52f3161b27SJose E. Roman PetscFunctionBegin; 53f3161b27SJose E. Roman /* Get preconditioner matrix from PETSc and setup pARMS structs */ 549566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &pmat)); 55*f1580f4eSBarry Smith PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pmat), &npro)); 56*f1580f4eSBarry Smith PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pmat), &rank)); 57f3161b27SJose E. Roman 589566063dSJacob Faibussowitsch PetscCall(MatGetSize(pmat, &n, NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npro + 1, &mapptr)); 609566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &maptmp)); 619566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(pmat, &mapptr0)); 62f3161b27SJose E. Roman low = mapptr0[rank]; 63f3161b27SJose E. Roman high = mapptr0[rank + 1]; 64f3161b27SJose E. Roman lsize = high - low; 65f3161b27SJose E. Roman 662fa5cd67SKarl Rupp for (i = 0; i < npro + 1; i++) mapptr[i] = mapptr0[i] + 1; 672fa5cd67SKarl Rupp for (i = 0; i < n; i++) maptmp[i] = i + 1; 68f3161b27SJose E. Roman 69f3161b27SJose E. Roman /* if created, destroy the previous map */ 70f3161b27SJose E. Roman if (parms->map) { 71f3161b27SJose E. Roman parms_MapFree(&parms->map); 720298fd71SBarry Smith parms->map = NULL; 73f3161b27SJose E. Roman } 74f3161b27SJose E. Roman 75f3161b27SJose E. Roman /* create pARMS map object */ 76ce94432eSBarry Smith parms_MapCreateFromPtr(&parms->map, (int)n, maptmp, mapptr, PetscObjectComm((PetscObject)pmat), 1, NONINTERLACED); 77f3161b27SJose E. Roman 78f3161b27SJose E. Roman /* if created, destroy the previous pARMS matrix */ 79f3161b27SJose E. Roman if (parms->A) { 80f3161b27SJose E. Roman parms_MatFree(&parms->A); 810298fd71SBarry Smith parms->A = NULL; 82f3161b27SJose E. Roman } 83f3161b27SJose E. Roman 84f3161b27SJose E. Roman /* create pARMS mat object */ 85f3161b27SJose E. Roman parms_MatCreate(&parms->A, parms->map); 86f3161b27SJose E. Roman 87f3161b27SJose E. Roman /* setup and copy csr data structure for pARMS */ 889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize + 1, &ia)); 89f3161b27SJose E. Roman ia[0] = 1; 909566063dSJacob Faibussowitsch PetscCall(MatGetInfo(pmat, MAT_LOCAL, &matinfo)); 91f3161b27SJose E. Roman length = matinfo.nz_used; 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length, &ja)); 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length, &aa)); 94f3161b27SJose E. Roman 95f3161b27SJose E. Roman for (i = low; i < high; i++) { 96f3161b27SJose E. Roman pos = ia[i - low] - 1; 979566063dSJacob Faibussowitsch PetscCall(MatGetRow(pmat, i, &ncols, &cols, &values)); 98f3161b27SJose E. Roman ia[i - low + 1] = ia[i - low] + ncols; 99f3161b27SJose E. Roman 100f3161b27SJose E. Roman if (ia[i - low + 1] >= length) { 101f3161b27SJose E. Roman length += ncols; 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length, &ja1)); 1039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ja1, ja, ia[i - low] - 1)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 105f3161b27SJose E. Roman ja = ja1; 1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(length, &aa1)); 1079566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(aa1, aa, ia[i - low] - 1)); 1089566063dSJacob Faibussowitsch PetscCall(PetscFree(aa)); 109f3161b27SJose E. Roman aa = aa1; 110f3161b27SJose E. Roman } 1119566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&ja[pos], cols, ncols)); 1129566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&aa[pos], values, ncols)); 1139566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(pmat, i, &ncols, &cols, &values)); 114f3161b27SJose E. Roman } 115f3161b27SJose E. Roman 116f3161b27SJose E. Roman /* csr info is for local matrix so initialize im[] locally */ 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &im)); 1189566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(im, &maptmp[mapptr[rank] - 1], lsize)); 119f3161b27SJose E. Roman 120f3161b27SJose E. Roman /* 1-based indexing */ 1212fa5cd67SKarl Rupp for (i = 0; i < ia[lsize] - 1; i++) ja[i] = ja[i] + 1; 122f3161b27SJose E. Roman 123f3161b27SJose E. Roman /* Now copy csr matrix to parms_mat object */ 124f3161b27SJose E. Roman parms_MatSetValues(parms->A, (int)lsize, im, ia, ja, aa, INSERT); 125f3161b27SJose E. Roman 126f3161b27SJose E. Roman /* free memory */ 1279566063dSJacob Faibussowitsch PetscCall(PetscFree(maptmp)); 1289566063dSJacob Faibussowitsch PetscCall(PetscFree(mapptr)); 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(aa)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(ja)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(ia)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(im)); 133f3161b27SJose E. Roman 134f3161b27SJose E. Roman /* setup parms matrix */ 135f3161b27SJose E. Roman parms_MatSetup(parms->A); 136f3161b27SJose E. Roman 137f3161b27SJose E. Roman /* if created, destroy the previous pARMS pc */ 138f3161b27SJose E. Roman if (parms->pc) { 139f3161b27SJose E. Roman parms_PCFree(&parms->pc); 1400298fd71SBarry Smith parms->pc = NULL; 141f3161b27SJose E. Roman } 142f3161b27SJose E. Roman 143f3161b27SJose E. Roman /* Now create pARMS preconditioner object based on A */ 144f3161b27SJose E. Roman parms_PCCreate(&parms->pc, parms->A); 145f3161b27SJose E. Roman 146f3161b27SJose E. Roman /* Transfer options from PC to pARMS */ 147f3161b27SJose E. Roman switch (parms->global) { 148f3161b27SJose E. Roman case 0: parms_PCSetType(parms->pc, PCRAS); break; 149f3161b27SJose E. Roman case 1: parms_PCSetType(parms->pc, PCSCHUR); break; 150f3161b27SJose E. Roman case 2: parms_PCSetType(parms->pc, PCBJ); break; 151f3161b27SJose E. Roman } 152f3161b27SJose E. Roman switch (parms->local) { 153f3161b27SJose E. Roman case 0: parms_PCSetILUType(parms->pc, PCILU0); break; 154f3161b27SJose E. Roman case 1: parms_PCSetILUType(parms->pc, PCILUK); break; 155f3161b27SJose E. Roman case 2: parms_PCSetILUType(parms->pc, PCILUT); break; 156f3161b27SJose E. Roman case 3: parms_PCSetILUType(parms->pc, PCARMS); break; 157f3161b27SJose E. Roman } 158f3161b27SJose E. Roman parms_PCSetInnerEps(parms->pc, parms->solvetol); 159f3161b27SJose E. Roman parms_PCSetNlevels(parms->pc, parms->levels); 160f3161b27SJose E. Roman parms_PCSetPermType(parms->pc, parms->nonsymperm ? 1 : 0); 161f3161b27SJose E. Roman parms_PCSetBsize(parms->pc, parms->blocksize); 162f3161b27SJose E. Roman parms_PCSetTolInd(parms->pc, parms->indtol); 163f3161b27SJose E. Roman parms_PCSetInnerKSize(parms->pc, parms->maxdim); 164f3161b27SJose E. Roman parms_PCSetInnerMaxits(parms->pc, parms->maxits); 165f3161b27SJose E. Roman for (i = 0; i < 8; i++) meth[i] = parms->meth[i] ? 1 : 0; 166f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[0], 1); 167f3161b27SJose E. Roman parms_PCSetPermScalOptions(parms->pc, &meth[4], 0); 168f3161b27SJose E. Roman parms_PCSetFill(parms->pc, parms->lfil); 169f3161b27SJose E. Roman parms_PCSetTol(parms->pc, parms->droptol); 170f3161b27SJose E. Roman 171f3161b27SJose E. Roman parms_PCSetup(parms->pc); 172f3161b27SJose E. Roman 173f3161b27SJose E. Roman /* Allocate two auxiliary vector of length lsize */ 1749566063dSJacob Faibussowitsch if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 1759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec0)); 1769566063dSJacob Faibussowitsch if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 1779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(lsize, &parms->lvec1)); 178f3161b27SJose E. Roman PetscFunctionReturn(0); 179f3161b27SJose E. Roman } 180f3161b27SJose E. Roman 1819371c9d4SSatish Balay static PetscErrorCode PCView_PARMS(PC pc, PetscViewer viewer) { 182f3161b27SJose E. Roman PetscBool iascii; 183f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 184f3161b27SJose E. Roman char *str; 185f3161b27SJose E. Roman double fill_fact; 186f3161b27SJose E. Roman 187f3161b27SJose E. Roman PetscFunctionBegin; 1889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 189f3161b27SJose E. Roman if (iascii) { 190f3161b27SJose E. Roman parms_PCGetName(parms->pc, &str); 1919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " global preconditioner: %s\n", str)); 192f3161b27SJose E. Roman parms_PCILUGetName(parms->pc, &str); 1939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " local preconditioner: %s\n", str)); 194f3161b27SJose E. Roman parms_PCGetRatio(parms->pc, &fill_fact); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " non-zero elements/original non-zero entries: %-4.2f\n", fill_fact)); 1969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Tolerance for local solve: %g\n", parms->solvetol)); 1979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Number of levels: %d\n", parms->levels)); 19848a46eb9SPierre Jolivet if (parms->nonsymperm) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation\n")); 1999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Block size: %d\n", parms->blocksize)); 2009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Tolerance for independent sets: %g\n", parms->indtol)); 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Inner Krylov dimension: %d\n", parms->maxdim)); 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Maximum number of inner iterations: %d\n", parms->maxits)); 20348a46eb9SPierre Jolivet if (parms->meth[0]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation for interlevel blocks\n")); 20448a46eb9SPierre Jolivet if (parms->meth[1]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column permutation for interlevel blocks\n")); 20548a46eb9SPierre Jolivet if (parms->meth[2]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using row scaling for interlevel blocks\n")); 20648a46eb9SPierre Jolivet if (parms->meth[3]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column scaling for interlevel blocks\n")); 20748a46eb9SPierre Jolivet if (parms->meth[4]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using nonsymmetric permutation for last level blocks\n")); 20848a46eb9SPierre Jolivet if (parms->meth[5]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column permutation for last level blocks\n")); 20948a46eb9SPierre Jolivet if (parms->meth[6]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using row scaling for last level blocks\n")); 21048a46eb9SPierre Jolivet if (parms->meth[7]) PetscCall(PetscViewerASCIIPrintf(viewer, " Using column scaling for last level blocks\n")); 2119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for ilut, iluk and arms: %d\n", parms->lfil[0])); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for schur: %d\n", parms->lfil[4])); 2139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " amount of fill-in for ILUT L and U: %d\n", parms->lfil[5])); 2149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for L, U, L^{-1}F and EU^{-1}: %g\n", parms->droptol[0])); 2159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for schur complement at each level: %g\n", parms->droptol[4])); 2169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " drop tolerance for ILUT in last level schur complement: %g\n", parms->droptol[5])); 217f3161b27SJose E. Roman } 218f3161b27SJose E. Roman PetscFunctionReturn(0); 219f3161b27SJose E. Roman } 220f3161b27SJose E. Roman 2219371c9d4SSatish Balay static PetscErrorCode PCDestroy_PARMS(PC pc) { 222f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 223f3161b27SJose E. Roman 224f3161b27SJose E. Roman PetscFunctionBegin; 2252fa5cd67SKarl Rupp if (parms->map) parms_MapFree(&parms->map); 2262fa5cd67SKarl Rupp if (parms->A) parms_MatFree(&parms->A); 2272fa5cd67SKarl Rupp if (parms->pc) parms_PCFree(&parms->pc); 2281baa6e33SBarry Smith if (parms->lvec0) PetscCall(PetscFree(parms->lvec0)); 2291baa6e33SBarry Smith if (parms->lvec1) PetscCall(PetscFree(parms->lvec1)); 2309566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 231f3161b27SJose E. Roman 2329566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)pc, 0)); 2339566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetGlobal_C", NULL)); 2349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetLocal_C", NULL)); 2359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveTolerances_C", NULL)); 2369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveRestart_C", NULL)); 2379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetNonsymPerm_C", NULL)); 2389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetFill_C", NULL)); 239f3161b27SJose E. Roman PetscFunctionReturn(0); 240f3161b27SJose E. Roman } 241f3161b27SJose E. Roman 2429371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_PARMS(PC pc, PetscOptionItems *PetscOptionsObject) { 243f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 244f3161b27SJose E. Roman PetscBool flag; 245f3161b27SJose E. Roman PCPARMSGlobalType global; 246f3161b27SJose E. Roman PCPARMSLocalType local; 247f3161b27SJose E. Roman 248f3161b27SJose E. Roman PetscFunctionBegin; 249d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PARMS Options"); 2509566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_global", "Global preconditioner", "PCPARMSSetGlobal", PCPARMSGlobalTypes, (PetscEnum)parms->global, (PetscEnum *)&global, &flag)); 2519566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetGlobal(pc, global)); 2529566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_parms_local", "Local preconditioner", "PCPARMSSetLocal", PCPARMSLocalTypes, (PetscEnum)parms->local, (PetscEnum *)&local, &flag)); 2539566063dSJacob Faibussowitsch if (flag) PetscCall(PCPARMSSetLocal(pc, local)); 2549566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_solve_tol", "Tolerance for local solve", "PCPARMSSetSolveTolerances", parms->solvetol, &parms->solvetol, NULL)); 2559566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_levels", "Number of levels", "None", parms->levels, &parms->levels, NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_nonsymmetric_perm", "Use nonsymmetric permutation", "PCPARMSSetNonsymPerm", parms->nonsymperm, &parms->nonsymperm, NULL)); 2579566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_blocksize", "Block size", "None", parms->blocksize, &parms->blocksize, NULL)); 2589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_ind_tol", "Tolerance for independent sets", "None", parms->indtol, &parms->indtol, NULL)); 2599566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_dim", "Inner Krylov dimension", "PCPARMSSetSolveRestart", parms->maxdim, &parms->maxdim, NULL)); 2609566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_max_it", "Maximum number of inner iterations", "PCPARMSSetSolveTolerances", parms->maxits, &parms->maxits, NULL)); 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm", "nonsymmetric permutation for interlevel blocks", "None", parms->meth[0], &parms->meth[0], NULL)); 2629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_perm", "column permutation for interlevel blocks", "None", parms->meth[1], &parms->meth[1], NULL)); 2639566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_row_scaling", "row scaling for interlevel blocks", "None", parms->meth[2], &parms->meth[2], NULL)); 2649566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_inter_column_scaling", "column scaling for interlevel blocks", "None", parms->meth[3], &parms->meth[3], NULL)); 2659566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_nonsymmetric_perm", "nonsymmetric permutation for last level blocks", "None", parms->meth[4], &parms->meth[4], NULL)); 2669566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_perm", "column permutation for last level blocks", "None", parms->meth[5], &parms->meth[5], NULL)); 2679566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_row_scaling", "row scaling for last level blocks", "None", parms->meth[6], &parms->meth[6], NULL)); 2689566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_parms_last_column_scaling", "column scaling for last level blocks", "None", parms->meth[7], &parms->meth[7], NULL)); 2699566063dSJacob 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)); 270f3161b27SJose E. Roman if (flag) parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0]; 2719566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_parms_lfil_schur", "amount of fill-in for schur", "PCPARMSSetFill", parms->lfil[4], &parms->lfil[4], NULL)); 2729566063dSJacob 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)); 273f3161b27SJose E. Roman if (flag) parms->lfil[6] = parms->lfil[5]; 2749566063dSJacob 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)); 2759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_parms_droptol_schur_compl", "drop tolerance for schur complement at each level", "None", parms->droptol[4], &parms->droptol[4], &flag)); 276f3161b27SJose E. Roman if (flag) parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = parms->droptol[0]; 2779566063dSJacob 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)); 278f3161b27SJose E. Roman if (flag) parms->droptol[6] = parms->droptol[5]; 279d0609cedSBarry Smith PetscOptionsHeadEnd(); 280f3161b27SJose E. Roman PetscFunctionReturn(0); 281f3161b27SJose E. Roman } 282f3161b27SJose E. Roman 2839371c9d4SSatish Balay static PetscErrorCode PCApply_PARMS(PC pc, Vec b, Vec x) { 284f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 285f3161b27SJose E. Roman const PetscScalar *b1; 286f3161b27SJose E. Roman PetscScalar *x1; 287f3161b27SJose E. Roman 288f3161b27SJose E. Roman PetscFunctionBegin; 2899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(b, &b1)); 2909566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x1)); 291f3161b27SJose E. Roman parms_VecPermAux((PetscScalar *)b1, parms->lvec0, parms->map); 292f3161b27SJose E. Roman parms_PCApply(parms->pc, parms->lvec0, parms->lvec1); 293f3161b27SJose E. Roman parms_VecInvPermAux(parms->lvec1, x1, parms->map); 2949566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(b, &b1)); 2959566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x1)); 296f3161b27SJose E. Roman PetscFunctionReturn(0); 297f3161b27SJose E. Roman } 298f3161b27SJose E. Roman 2999371c9d4SSatish Balay static PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc, PCPARMSGlobalType type) { 300f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 301f3161b27SJose E. Roman 302f3161b27SJose E. Roman PetscFunctionBegin; 303f3161b27SJose E. Roman if (type != parms->global) { 304f3161b27SJose E. Roman parms->global = type; 305f3161b27SJose E. Roman pc->setupcalled = 0; 306f3161b27SJose E. Roman } 307f3161b27SJose E. Roman PetscFunctionReturn(0); 308f3161b27SJose E. Roman } 309f3161b27SJose E. Roman 3105de0dacdSJose E. Roman /*@ 311*f1580f4eSBarry Smith PCPARMSSetGlobal - Sets the global preconditioner to be used in `PCPARMS`. 312f3161b27SJose E. Roman 313*f1580f4eSBarry Smith Collective on pc 314f3161b27SJose E. Roman 315f3161b27SJose E. Roman Input Parameters: 316f3161b27SJose E. Roman + pc - the preconditioner context 317f3161b27SJose E. Roman - type - the global preconditioner type, one of 318f3161b27SJose E. Roman .vb 319f3161b27SJose E. Roman PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 320f3161b27SJose E. Roman PC_PARMS_GLOBAL_SCHUR - Schur complement 321f3161b27SJose E. Roman PC_PARMS_GLOBAL_BJ - Block Jacobi 322f3161b27SJose E. Roman .ve 323f3161b27SJose E. Roman 324*f1580f4eSBarry Smith Options Database Key: 325*f1580f4eSBarry Smith . -pc_parms_global [ras,schur,bj] - Sets global preconditioner 326f3161b27SJose E. Roman 327f3161b27SJose E. Roman Level: intermediate 328f3161b27SJose E. Roman 329*f1580f4eSBarry Smith Note: 330*f1580f4eSBarry Smith See the pARMS function `parms_PCSetType()` for more information. 331f3161b27SJose E. Roman 332db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetLocal()` 333f3161b27SJose E. Roman @*/ 3349371c9d4SSatish Balay PetscErrorCode PCPARMSSetGlobal(PC pc, PCPARMSGlobalType type) { 335f3161b27SJose E. Roman PetscFunctionBegin; 336f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 337f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc, type, 2); 338cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetGlobal_C", (PC, PCPARMSGlobalType), (pc, type)); 339f3161b27SJose E. Roman PetscFunctionReturn(0); 340f3161b27SJose E. Roman } 341f3161b27SJose E. Roman 3429371c9d4SSatish Balay static PetscErrorCode PCPARMSSetLocal_PARMS(PC pc, PCPARMSLocalType type) { 343f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 344f3161b27SJose E. Roman 345f3161b27SJose E. Roman PetscFunctionBegin; 346f3161b27SJose E. Roman if (type != parms->local) { 347f3161b27SJose E. Roman parms->local = type; 348f3161b27SJose E. Roman pc->setupcalled = 0; 349f3161b27SJose E. Roman } 350f3161b27SJose E. Roman PetscFunctionReturn(0); 351f3161b27SJose E. Roman } 352f3161b27SJose E. Roman 3535de0dacdSJose E. Roman /*@ 354*f1580f4eSBarry Smith PCPARMSSetLocal - Sets the local preconditioner to be used in `PCPARMS`. 355f3161b27SJose E. Roman 356*f1580f4eSBarry Smith Collective on pc 357f3161b27SJose E. Roman 358f3161b27SJose E. Roman Input Parameters: 359f3161b27SJose E. Roman + pc - the preconditioner context 360f3161b27SJose E. Roman - type - the local preconditioner type, one of 361f3161b27SJose E. Roman .vb 362f3161b27SJose E. Roman PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 363f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 364f3161b27SJose E. Roman PC_PARMS_LOCAL_ILUT - ILUT preconditioner 365f3161b27SJose E. Roman PC_PARMS_LOCAL_ARMS - ARMS preconditioner 366f3161b27SJose E. Roman .ve 367f3161b27SJose E. Roman 368f3161b27SJose E. Roman Options Database Keys: 369f3161b27SJose E. Roman -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 370f3161b27SJose E. Roman 371f3161b27SJose E. Roman Level: intermediate 372f3161b27SJose E. Roman 373f3161b27SJose E. Roman Notes: 374f3161b27SJose E. Roman For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 375f3161b27SJose E. Roman variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 376f3161b27SJose E. Roman 377*f1580f4eSBarry Smith See the pARMS function `parms_PCILUSetType()` for more information. 378f3161b27SJose E. Roman 379db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetGlobal()`, `PCPARMSSetNonsymPerm()` 380f3161b27SJose E. Roman 381f3161b27SJose E. Roman @*/ 3829371c9d4SSatish Balay PetscErrorCode PCPARMSSetLocal(PC pc, PCPARMSLocalType type) { 383f3161b27SJose E. Roman PetscFunctionBegin; 384f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 385f3161b27SJose E. Roman PetscValidLogicalCollectiveEnum(pc, type, 2); 386cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetLocal_C", (PC, PCPARMSLocalType), (pc, type)); 387f3161b27SJose E. Roman PetscFunctionReturn(0); 388f3161b27SJose E. Roman } 389f3161b27SJose E. Roman 3909371c9d4SSatish Balay static PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc, PetscReal tol, PetscInt maxits) { 391f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 392f3161b27SJose E. Roman 393f3161b27SJose E. Roman PetscFunctionBegin; 394f3161b27SJose E. Roman if (tol != parms->solvetol) { 395f3161b27SJose E. Roman parms->solvetol = tol; 396f3161b27SJose E. Roman pc->setupcalled = 0; 397f3161b27SJose E. Roman } 398f3161b27SJose E. Roman if (maxits != parms->maxits) { 399f3161b27SJose E. Roman parms->maxits = maxits; 400f3161b27SJose E. Roman pc->setupcalled = 0; 401f3161b27SJose E. Roman } 402f3161b27SJose E. Roman PetscFunctionReturn(0); 403f3161b27SJose E. Roman } 404f3161b27SJose E. Roman 4055de0dacdSJose E. Roman /*@ 406f3161b27SJose E. Roman PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 407f3161b27SJose E. Roman inner GMRES solver, when the Schur global preconditioner is used. 408f3161b27SJose E. Roman 409*f1580f4eSBarry Smith Collective on pc 410f3161b27SJose E. Roman 411f3161b27SJose E. Roman Input Parameters: 412f3161b27SJose E. Roman + pc - the preconditioner context 413f3161b27SJose E. Roman . tol - the convergence tolerance 414f3161b27SJose E. Roman - maxits - the maximum number of iterations to use 415f3161b27SJose E. Roman 416f3161b27SJose E. Roman Options Database Keys: 417f3161b27SJose E. Roman + -pc_parms_solve_tol - set the tolerance for local solve 418f3161b27SJose E. Roman - -pc_parms_max_it - set the maximum number of inner iterations 419f3161b27SJose E. Roman 420f3161b27SJose E. Roman Level: intermediate 421f3161b27SJose E. Roman 422*f1580f4eSBarry Smith Note: 423*f1580f4eSBarry Smith See the pARMS functions `parms_PCSetInnerEps()` and `parms_PCSetInnerMaxits()` for more information. 424f3161b27SJose E. Roman 425db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetSolveRestart()` 426f3161b27SJose E. Roman @*/ 4279371c9d4SSatish Balay PetscErrorCode PCPARMSSetSolveTolerances(PC pc, PetscReal tol, PetscInt maxits) { 428f3161b27SJose E. Roman PetscFunctionBegin; 429f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 430cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetSolveTolerances_C", (PC, PetscReal, PetscInt), (pc, tol, maxits)); 431f3161b27SJose E. Roman PetscFunctionReturn(0); 432f3161b27SJose E. Roman } 433f3161b27SJose E. Roman 4349371c9d4SSatish Balay static PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc, PetscInt restart) { 435f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 436f3161b27SJose E. Roman 437f3161b27SJose E. Roman PetscFunctionBegin; 438f3161b27SJose E. Roman if (restart != parms->maxdim) { 439f3161b27SJose E. Roman parms->maxdim = restart; 440f3161b27SJose E. Roman pc->setupcalled = 0; 441f3161b27SJose E. Roman } 442f3161b27SJose E. Roman PetscFunctionReturn(0); 443f3161b27SJose E. Roman } 444f3161b27SJose E. Roman 4455de0dacdSJose E. Roman /*@ 446f3161b27SJose E. Roman PCPARMSSetSolveRestart - Sets the number of iterations at which the 447f3161b27SJose E. Roman inner GMRES solver restarts. 448f3161b27SJose E. Roman 449*f1580f4eSBarry Smith Collective on pc 450f3161b27SJose E. Roman 451f3161b27SJose E. Roman Input Parameters: 452f3161b27SJose E. Roman + pc - the preconditioner context 453f3161b27SJose E. Roman - restart - maximum dimension of the Krylov subspace 454f3161b27SJose E. Roman 455*f1580f4eSBarry Smith Options Database Key: 456f3161b27SJose E. Roman . -pc_parms_max_dim - sets the inner Krylov dimension 457f3161b27SJose E. Roman 458f3161b27SJose E. Roman Level: intermediate 459f3161b27SJose E. Roman 460*f1580f4eSBarry Smith Note: 461f3161b27SJose E. Roman See the pARMS function parms_PCSetInnerKSize for more information. 462f3161b27SJose E. Roman 463db781477SPatrick Sanan .seealso: `PCPARMS`, `PCPARMSSetSolveTolerances()` 464f3161b27SJose E. Roman @*/ 4659371c9d4SSatish Balay PetscErrorCode PCPARMSSetSolveRestart(PC pc, PetscInt restart) { 466f3161b27SJose E. Roman PetscFunctionBegin; 467f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 468cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetSolveRestart_C", (PC, PetscInt), (pc, restart)); 469f3161b27SJose E. Roman PetscFunctionReturn(0); 470f3161b27SJose E. Roman } 471f3161b27SJose E. Roman 4729371c9d4SSatish Balay static PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc, PetscBool nonsym) { 473f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 474f3161b27SJose E. Roman 475f3161b27SJose E. Roman PetscFunctionBegin; 4765de0dacdSJose E. Roman if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 477f3161b27SJose E. Roman parms->nonsymperm = nonsym; 478f3161b27SJose E. Roman pc->setupcalled = 0; 479f3161b27SJose E. Roman } 480f3161b27SJose E. Roman PetscFunctionReturn(0); 481f3161b27SJose E. Roman } 482f3161b27SJose E. Roman 4835de0dacdSJose E. Roman /*@ 484f3161b27SJose E. Roman PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 485f3161b27SJose E. Roman symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 486f3161b27SJose E. Roman 487*f1580f4eSBarry Smith Collective on pc 488f3161b27SJose E. Roman 489f3161b27SJose E. Roman Input Parameters: 490f3161b27SJose E. Roman + pc - the preconditioner context 491*f1580f4eSBarry Smith - nonsym - `PETSC_TRUE` indicates the non-symmetric ARMS is used; 492*f1580f4eSBarry Smith `PETSC_FALSE` indicates the symmetric ARMS is used 493f3161b27SJose E. Roman 494*f1580f4eSBarry Smith Options Database Key: 495f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 496f3161b27SJose E. Roman 497f3161b27SJose E. Roman Level: intermediate 498f3161b27SJose E. Roman 499*f1580f4eSBarry Smith Note: 500*f1580f4eSBarry Smith See the pARMS function `parms_PCSetPermType()` for more information. 501f3161b27SJose E. Roman 502db781477SPatrick Sanan .seealso: `PCPARMS` 503f3161b27SJose E. Roman @*/ 5049371c9d4SSatish Balay PetscErrorCode PCPARMSSetNonsymPerm(PC pc, PetscBool nonsym) { 505f3161b27SJose E. Roman PetscFunctionBegin; 506f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 507cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetNonsymPerm_C", (PC, PetscBool), (pc, nonsym)); 508f3161b27SJose E. Roman PetscFunctionReturn(0); 509f3161b27SJose E. Roman } 510f3161b27SJose E. Roman 5119371c9d4SSatish Balay static PetscErrorCode PCPARMSSetFill_PARMS(PC pc, PetscInt lfil0, PetscInt lfil1, PetscInt lfil2) { 512f3161b27SJose E. Roman PC_PARMS *parms = (PC_PARMS *)pc->data; 513f3161b27SJose E. Roman 514f3161b27SJose E. Roman PetscFunctionBegin; 515f3161b27SJose E. Roman if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 516f3161b27SJose E. Roman parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 517f3161b27SJose E. Roman pc->setupcalled = 0; 518f3161b27SJose E. Roman } 519f3161b27SJose E. Roman if (lfil1 != parms->lfil[4]) { 520f3161b27SJose E. Roman parms->lfil[4] = lfil1; 521f3161b27SJose E. Roman pc->setupcalled = 0; 522f3161b27SJose E. Roman } 523f3161b27SJose E. Roman if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 524f3161b27SJose E. Roman parms->lfil[5] = parms->lfil[6] = lfil2; 525f3161b27SJose E. Roman pc->setupcalled = 0; 526f3161b27SJose E. Roman } 527f3161b27SJose E. Roman PetscFunctionReturn(0); 528f3161b27SJose E. Roman } 529f3161b27SJose E. Roman 5305de0dacdSJose E. Roman /*@ 531f3161b27SJose E. Roman PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 532f3161b27SJose E. Roman Consider the original matrix A = [B F; E C] and the approximate version 533f3161b27SJose E. Roman M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 534f3161b27SJose E. Roman 535*f1580f4eSBarry Smith Collective on pc 536f3161b27SJose E. Roman 537f3161b27SJose E. Roman Input Parameters: 538f3161b27SJose E. Roman + pc - the preconditioner context 539f3161b27SJose E. Roman . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 540f3161b27SJose E. Roman . fil1 - the level of fill-in kept in S 541f3161b27SJose E. Roman - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 542f3161b27SJose E. Roman 543f3161b27SJose E. Roman Options Database Keys: 544f3161b27SJose E. Roman + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 545f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 546f3161b27SJose E. Roman - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 547f3161b27SJose E. Roman 548f3161b27SJose E. Roman Level: intermediate 549f3161b27SJose E. Roman 550*f1580f4eSBarry Smith Note: 551*f1580f4eSBarry Smith See the pARMS function `parms_PCSetFill()` for more information. 552f3161b27SJose E. Roman 553db781477SPatrick Sanan .seealso: `PCPARMS` 554f3161b27SJose E. Roman @*/ 5559371c9d4SSatish Balay PetscErrorCode PCPARMSSetFill(PC pc, PetscInt lfil0, PetscInt lfil1, PetscInt lfil2) { 556f3161b27SJose E. Roman PetscFunctionBegin; 557f3161b27SJose E. Roman PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 558cac4c232SBarry Smith PetscTryMethod(pc, "PCPARMSSetFill_C", (PC, PetscInt, PetscInt, PetscInt), (pc, lfil0, lfil1, lfil2)); 559f3161b27SJose E. Roman PetscFunctionReturn(0); 560f3161b27SJose E. Roman } 561f3161b27SJose E. Roman 562f3161b27SJose E. Roman /*MC 563f3161b27SJose E. Roman PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 564f3161b27SJose E. Roman available in the package pARMS 565f3161b27SJose E. Roman 566f3161b27SJose E. Roman Options Database Keys: 567f3161b27SJose E. Roman + -pc_parms_global - one of ras, schur, bj 568f3161b27SJose E. Roman . -pc_parms_local - one of ilu0, iluk, ilut, arms 569f3161b27SJose E. Roman . -pc_parms_solve_tol - set the tolerance for local solve 570f3161b27SJose E. Roman . -pc_parms_levels - set the number of levels 571f3161b27SJose E. Roman . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 572f3161b27SJose E. Roman . -pc_parms_blocksize - set the block size 573f3161b27SJose E. Roman . -pc_parms_ind_tol - set the tolerance for independent sets 574f3161b27SJose E. Roman . -pc_parms_max_dim - set the inner krylov dimension 575f3161b27SJose E. Roman . -pc_parms_max_it - set the maximum number of inner iterations 5764b62eef3SJose E. Roman . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 577f3161b27SJose E. Roman . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 578f3161b27SJose E. Roman . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 579f3161b27SJose E. Roman . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 5804b62eef3SJose E. Roman . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 581f3161b27SJose E. Roman . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 582f3161b27SJose E. Roman . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 583f3161b27SJose E. Roman . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 584f3161b27SJose E. Roman . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 585f3161b27SJose E. Roman . -pc_parms_lfil_schur - set the amount of fill-in for schur 586f3161b27SJose E. Roman . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 587f3161b27SJose E. Roman . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 588f3161b27SJose E. Roman . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 589f3161b27SJose E. Roman - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 590f3161b27SJose E. Roman 591*f1580f4eSBarry Smith Note: 592f3161b27SJose E. Roman Unless configured appropriately, this preconditioner performs an inexact solve 593f3161b27SJose E. Roman as part of the preconditioner application. Therefore, it must be used in combination 594*f1580f4eSBarry Smith with flexible variants of iterative solvers, such as `KSPFGMRES` or `KSPGCR`. 595f3161b27SJose E. Roman 596f3161b27SJose E. Roman Level: intermediate 597f3161b27SJose E. Roman 598*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMG`, `PCGAMG`, `PCHYPRE`, `PCPARMSSetGlobal()`, 599*f1580f4eSBarry Smith `PCPARMSSetLocal()`, `PCPARMSSetSolveTolerances()`, `PCPARMSSetSolveRestart()`, `PCPARMSSetNonsymPerm()`, 600*f1580f4eSBarry Smith `PCPARMSSetFill()` 601f3161b27SJose E. Roman M*/ 602f3161b27SJose E. Roman 6039371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_PARMS(PC pc) { 604f3161b27SJose E. Roman PC_PARMS *parms; 605f3161b27SJose E. Roman 606f3161b27SJose E. Roman PetscFunctionBegin; 6079566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &parms)); 6082fa5cd67SKarl Rupp 609f3161b27SJose E. Roman parms->map = 0; 610f3161b27SJose E. Roman parms->A = 0; 611f3161b27SJose E. Roman parms->pc = 0; 612abd6c125SJose E. Roman parms->global = PC_PARMS_GLOBAL_RAS; 613f3161b27SJose E. Roman parms->local = PC_PARMS_LOCAL_ARMS; 614f3161b27SJose E. Roman parms->levels = 10; 615f3161b27SJose E. Roman parms->nonsymperm = PETSC_TRUE; 616f3161b27SJose E. Roman parms->blocksize = 250; 617abd6c125SJose E. Roman parms->maxdim = 0; 618abd6c125SJose E. Roman parms->maxits = 0; 6194b62eef3SJose E. Roman parms->meth[0] = PETSC_FALSE; 6204b62eef3SJose E. Roman parms->meth[1] = PETSC_FALSE; 6214b62eef3SJose E. Roman parms->meth[2] = PETSC_FALSE; 6224b62eef3SJose E. Roman parms->meth[3] = PETSC_FALSE; 6234b62eef3SJose E. Roman parms->meth[4] = PETSC_FALSE; 6244b62eef3SJose E. Roman parms->meth[5] = PETSC_FALSE; 6254b62eef3SJose E. Roman parms->meth[6] = PETSC_FALSE; 6264b62eef3SJose E. Roman parms->meth[7] = PETSC_FALSE; 627f3161b27SJose E. Roman parms->solvetol = 0.01; 628f3161b27SJose E. Roman parms->indtol = 0.4; 629f3161b27SJose E. Roman parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 630f3161b27SJose E. Roman parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 631f3161b27SJose E. Roman parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 632f3161b27SJose E. Roman parms->droptol[4] = 0.001; 633f3161b27SJose E. Roman parms->droptol[5] = parms->droptol[6] = 0.001; 6342fa5cd67SKarl Rupp 635f3161b27SJose E. Roman pc->data = parms; 636f3161b27SJose E. Roman pc->ops->destroy = PCDestroy_PARMS; 637f3161b27SJose E. Roman pc->ops->setfromoptions = PCSetFromOptions_PARMS; 638f3161b27SJose E. Roman pc->ops->setup = PCSetUp_PARMS; 639f3161b27SJose E. Roman pc->ops->apply = PCApply_PARMS; 640f3161b27SJose E. Roman pc->ops->view = PCView_PARMS; 6412fa5cd67SKarl Rupp 6429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetGlobal_C", PCPARMSSetGlobal_PARMS)); 6439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetLocal_C", PCPARMSSetLocal_PARMS)); 6449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveTolerances_C", PCPARMSSetSolveTolerances_PARMS)); 6459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetSolveRestart_C", PCPARMSSetSolveRestart_PARMS)); 6469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetNonsymPerm_C", PCPARMSSetNonsymPerm_PARMS)); 6479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCPARMSSetFill_C", PCPARMSSetFill_PARMS)); 648f3161b27SJose E. Roman PetscFunctionReturn(0); 649f3161b27SJose E. Roman } 650