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