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