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