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