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++) 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 = PETSC_NULL; 78 } 79 80 /* create pARMS map object */ 81 parms_MapCreateFromPtr(&parms->map,(int)n,maptmp,mapptr,((PetscObject)pmat)->comm,1,NONINTERLACED); 82 83 /* if created, destroy the previous pARMS matrix */ 84 if (parms->A) { 85 parms_MatFree(&parms->A); 86 parms->A = PETSC_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 = PetscMalloc((lsize+1)*sizeof(int),&ia);CHKERRQ(ierr); 94 ia[0] = 1; 95 ierr = MatGetInfo(pmat,MAT_LOCAL,&matinfo);CHKERRQ(ierr); 96 length = matinfo.nz_used; 97 ierr = PetscMalloc(length*sizeof(int),&ja);CHKERRQ(ierr); 98 ierr = PetscMalloc(length*sizeof(PetscScalar),&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 = PetscMalloc(length*sizeof(int),&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 = PetscMalloc(length*sizeof(PetscScalar),&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 = PetscMalloc(lsize*sizeof(int),&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 = PETSC_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 = PetscMalloc(lsize*sizeof(PetscScalar), &parms->lvec0);CHKERRQ(ierr); 181 if (parms->lvec1) { ierr = PetscFree(parms->lvec1);CHKERRQ(ierr); } 182 ierr = PetscMalloc(lsize*sizeof(PetscScalar), &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 = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetGlobal_C","",PETSC_NULL);CHKERRQ(ierr); 269 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetLocal_C","",PETSC_NULL);CHKERRQ(ierr); 270 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveTolerances_C","",PETSC_NULL);CHKERRQ(ierr); 271 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveRestart_C","",PETSC_NULL);CHKERRQ(ierr); 272 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetNonsymPerm_C","",PETSC_NULL);CHKERRQ(ierr); 273 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetFill_C","",PETSC_NULL);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 #undef __FUNCT__ 278 #define __FUNCT__ "PCSetFromOptions_PARMS" 279 static PetscErrorCode PCSetFromOptions_PARMS(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("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,&flag);CHKERRQ(ierr); 294 ierr = PetscOptionsInt("-pc_parms_levels","Number of levels","None",parms->levels,&parms->levels,&flag);CHKERRQ(ierr); 295 ierr = PetscOptionsBool("-pc_parms_nonsymmetric_perm","Use nonsymmetric permutation","PCPARMSSetNonsymPerm",parms->nonsymperm,&parms->nonsymperm,&flag);CHKERRQ(ierr); 296 ierr = PetscOptionsInt("-pc_parms_blocksize","Block size","None",parms->blocksize,&parms->blocksize,&flag);CHKERRQ(ierr); 297 ierr = PetscOptionsReal("-pc_parms_ind_tol","Tolerance for independent sets","None",parms->indtol,&parms->indtol,&flag);CHKERRQ(ierr); 298 ierr = PetscOptionsInt("-pc_parms_max_dim","Inner Krylov dimension","PCPARMSSetSolveRestart",parms->maxdim,&parms->maxdim,&flag);CHKERRQ(ierr); 299 ierr = PetscOptionsInt("-pc_parms_max_it","Maximum number of inner iterations","PCPARMSSetSolveTolerances",parms->maxits,&parms->maxits,&flag);CHKERRQ(ierr); 300 ierr = PetscOptionsBool("-pc_parms_inter_nonsymmetric_perm","nonsymmetric permutation for interlevel blocks","None",parms->meth[0],&parms->meth[0],&flag);CHKERRQ(ierr); 301 ierr = PetscOptionsBool("-pc_parms_inter_column_perm","column permutation for interlevel blocks","None",parms->meth[1],&parms->meth[1],&flag);CHKERRQ(ierr); 302 ierr = PetscOptionsBool("-pc_parms_inter_row_scaling","row scaling for interlevel blocks","None",parms->meth[2],&parms->meth[2],&flag);CHKERRQ(ierr); 303 ierr = PetscOptionsBool("-pc_parms_inter_column_scaling","column scaling for interlevel blocks","None",parms->meth[3],&parms->meth[3],&flag);CHKERRQ(ierr); 304 ierr = PetscOptionsBool("-pc_parms_last_nonsymmetric_perm","nonsymmetric permutation for last level blocks","None",parms->meth[4],&parms->meth[4],&flag);CHKERRQ(ierr); 305 ierr = PetscOptionsBool("-pc_parms_last_column_perm","column permutation for last level blocks","None",parms->meth[5],&parms->meth[5],&flag);CHKERRQ(ierr); 306 ierr = PetscOptionsBool("-pc_parms_last_row_scaling","row scaling for last level blocks","None",parms->meth[6],&parms->meth[6],&flag);CHKERRQ(ierr); 307 ierr = PetscOptionsBool("-pc_parms_last_column_scaling","column scaling for last level blocks","None",parms->meth[7],&parms->meth[7],&flag);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],&flag);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],&flag);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 EXTERN_C_BEGIN 343 #undef __FUNCT__ 344 #define __FUNCT__ "PCPARMSSetGlobal_PARMS" 345 PetscErrorCode PCPARMSSetGlobal_PARMS(PC pc,PCPARMSGlobalType type) 346 { 347 PC_PARMS *parms = (PC_PARMS*)pc->data; 348 349 PetscFunctionBegin; 350 if (type != parms->global) { 351 parms->global = type; 352 pc->setupcalled = 0; 353 } 354 PetscFunctionReturn(0); 355 } 356 EXTERN_C_END 357 358 #undef __FUNCT__ 359 #define __FUNCT__ "PCPARMSSetGlobal" 360 /*@ 361 PCPARMSSetGlobal - Sets the global preconditioner to be used in PARMS. 362 363 Collective on PC 364 365 Input Parameters: 366 + pc - the preconditioner context 367 - type - the global preconditioner type, one of 368 .vb 369 PC_PARMS_GLOBAL_RAS - Restricted additive Schwarz 370 PC_PARMS_GLOBAL_SCHUR - Schur complement 371 PC_PARMS_GLOBAL_BJ - Block Jacobi 372 .ve 373 374 Options Database Keys: 375 -pc_parms_global [ras,schur,bj] - Sets global preconditioner 376 377 Level: intermediate 378 379 Notes: 380 See the pARMS function parms_PCSetType for more information. 381 382 .seealso: PCPARMS, PCPARMSSetLocal() 383 @*/ 384 PetscErrorCode PCPARMSSetGlobal(PC pc,PCPARMSGlobalType type) 385 { 386 PetscErrorCode ierr; 387 388 PetscFunctionBegin; 389 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390 PetscValidLogicalCollectiveEnum(pc,type,2); 391 ierr = PetscTryMethod(pc,"PCPARMSSetGlobal_C",(PC,PCPARMSGlobalType),(pc,type));CHKERRQ(ierr); 392 PetscFunctionReturn(0); 393 } 394 395 EXTERN_C_BEGIN 396 #undef __FUNCT__ 397 #define __FUNCT__ "PCPARMSSetLocal_PARMS" 398 PetscErrorCode PCPARMSSetLocal_PARMS(PC pc,PCPARMSLocalType type) 399 { 400 PC_PARMS *parms = (PC_PARMS*)pc->data; 401 402 PetscFunctionBegin; 403 if (type != parms->local) { 404 parms->local = type; 405 pc->setupcalled = 0; 406 } 407 PetscFunctionReturn(0); 408 } 409 EXTERN_C_END 410 411 #undef __FUNCT__ 412 #define __FUNCT__ "PCPARMSSetLocal" 413 /*@ 414 PCPARMSSetLocal - Sets the local preconditioner to be used in PARMS. 415 416 Collective on PC 417 418 Input Parameters: 419 + pc - the preconditioner context 420 - type - the local preconditioner type, one of 421 .vb 422 PC_PARMS_LOCAL_ILU0 - ILU0 preconditioner 423 PC_PARMS_LOCAL_ILUK - ILU(k) preconditioner 424 PC_PARMS_LOCAL_ILUT - ILUT preconditioner 425 PC_PARMS_LOCAL_ARMS - ARMS preconditioner 426 .ve 427 428 Options Database Keys: 429 -pc_parms_local [ilu0,iluk,ilut,arms] - Sets local preconditioner 430 431 Level: intermediate 432 433 Notes: 434 For the ARMS preconditioner, one can use either the symmetric ARMS or the non-symmetric 435 variant (ARMS-ddPQ) by setting the permutation type with PCPARMSSetNonsymPerm(). 436 437 See the pARMS function parms_PCILUSetType for more information. 438 439 .seealso: PCPARMS, PCPARMSSetGlobal(), PCPARMSSetNonsymPerm() 440 441 @*/ 442 PetscErrorCode PCPARMSSetLocal(PC pc,PCPARMSLocalType type) 443 { 444 PetscErrorCode ierr; 445 446 PetscFunctionBegin; 447 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 448 PetscValidLogicalCollectiveEnum(pc,type,2); 449 ierr = PetscTryMethod(pc,"PCPARMSSetLocal_C",(PC,PCPARMSLocalType),(pc,type));CHKERRQ(ierr); 450 PetscFunctionReturn(0); 451 } 452 453 EXTERN_C_BEGIN 454 #undef __FUNCT__ 455 #define __FUNCT__ "PCPARMSSetSolveTolerances_PARMS" 456 PetscErrorCode PCPARMSSetSolveTolerances_PARMS(PC pc,PetscReal tol,PetscInt maxits) 457 { 458 PC_PARMS *parms = (PC_PARMS*)pc->data; 459 460 PetscFunctionBegin; 461 if (tol != parms->solvetol) { 462 parms->solvetol = tol; 463 pc->setupcalled = 0; 464 } 465 if (maxits != parms->maxits) { 466 parms->maxits = maxits; 467 pc->setupcalled = 0; 468 } 469 PetscFunctionReturn(0); 470 } 471 EXTERN_C_END 472 473 #undef __FUNCT__ 474 #define __FUNCT__ "PCPARMSSetSolveTolerances" 475 /*@ 476 PCPARMSSetSolveTolerances - Sets the convergence tolerance and the maximum iterations for the 477 inner GMRES solver, when the Schur global preconditioner is used. 478 479 Collective on PC 480 481 Input Parameters: 482 + pc - the preconditioner context 483 . tol - the convergence tolerance 484 - maxits - the maximum number of iterations to use 485 486 Options Database Keys: 487 + -pc_parms_solve_tol - set the tolerance for local solve 488 - -pc_parms_max_it - set the maximum number of inner iterations 489 490 Level: intermediate 491 492 Notes: 493 See the pARMS functions parms_PCSetInnerEps and parms_PCSetInnerMaxits for more information. 494 495 .seealso: PCPARMS, PCPARMSSetSolveRestart() 496 @*/ 497 PetscErrorCode PCPARMSSetSolveTolerances(PC pc,PetscReal tol,PetscInt maxits) 498 { 499 PetscErrorCode ierr; 500 501 PetscFunctionBegin; 502 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 503 ierr = PetscTryMethod(pc,"PCPARMSSetSolveTolerances_C",(PC,PetscReal,PetscInt),(pc,tol,maxits));CHKERRQ(ierr); 504 PetscFunctionReturn(0); 505 } 506 507 EXTERN_C_BEGIN 508 #undef __FUNCT__ 509 #define __FUNCT__ "PCPARMSSetSolveRestart_PARMS" 510 PetscErrorCode PCPARMSSetSolveRestart_PARMS(PC pc,PetscInt restart) 511 { 512 PC_PARMS *parms = (PC_PARMS*)pc->data; 513 514 PetscFunctionBegin; 515 if (restart != parms->maxdim) { 516 parms->maxdim = restart; 517 pc->setupcalled = 0; 518 } 519 PetscFunctionReturn(0); 520 } 521 EXTERN_C_END 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "PCPARMSSetSolveRestart" 525 /*@ 526 PCPARMSSetSolveRestart - Sets the number of iterations at which the 527 inner GMRES solver restarts. 528 529 Collective on PC 530 531 Input Parameters: 532 + pc - the preconditioner context 533 - restart - maximum dimension of the Krylov subspace 534 535 Options Database Keys: 536 . -pc_parms_max_dim - sets the inner Krylov dimension 537 538 Level: intermediate 539 540 Notes: 541 See the pARMS function parms_PCSetInnerKSize for more information. 542 543 .seealso: PCPARMS, PCPARMSSetSolveTolerances() 544 @*/ 545 PetscErrorCode PCPARMSSetSolveRestart(PC pc,PetscInt restart) 546 { 547 PetscErrorCode ierr; 548 549 PetscFunctionBegin; 550 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 551 ierr = PetscTryMethod(pc,"PCPARMSSetSolveRestart_C",(PC,PetscInt),(pc,restart));CHKERRQ(ierr); 552 PetscFunctionReturn(0); 553 } 554 555 EXTERN_C_BEGIN 556 #undef __FUNCT__ 557 #define __FUNCT__ "PCPARMSSetNonsymPerm_PARMS" 558 PetscErrorCode PCPARMSSetNonsymPerm_PARMS(PC pc,PetscBool nonsym) 559 { 560 PC_PARMS *parms = (PC_PARMS*)pc->data; 561 562 PetscFunctionBegin; 563 if ((nonsym && !parms->nonsymperm) || (!nonsym && parms->nonsymperm)) { 564 parms->nonsymperm = nonsym; 565 pc->setupcalled = 0; 566 } 567 PetscFunctionReturn(0); 568 } 569 EXTERN_C_END 570 571 #undef __FUNCT__ 572 #define __FUNCT__ "PCPARMSSetNonsymPerm" 573 /*@ 574 PCPARMSSetNonsymPerm - Sets the type of permutation for the ARMS preconditioner: the standard 575 symmetric ARMS or the non-symmetric ARMS (ARMS-ddPQ). 576 577 Collective on PC 578 579 Input Parameters: 580 + pc - the preconditioner context 581 - nonsym - PETSC_TRUE indicates the non-symmetric ARMS is used; 582 PETSC_FALSE indicates the symmetric ARMS is used 583 584 Options Database Keys: 585 . -pc_parms_nonsymmetric_perm - sets the use of nonsymmetric permutation 586 587 Level: intermediate 588 589 Notes: 590 See the pARMS function parms_PCSetPermType for more information. 591 592 .seealso: PCPARMS 593 @*/ 594 PetscErrorCode PCPARMSSetNonsymPerm(PC pc,PetscBool nonsym) 595 { 596 PetscErrorCode ierr; 597 598 PetscFunctionBegin; 599 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 600 ierr = PetscTryMethod(pc,"PCPARMSSetNonsymPerm_C",(PC,PetscBool),(pc,nonsym));CHKERRQ(ierr); 601 PetscFunctionReturn(0); 602 } 603 604 EXTERN_C_BEGIN 605 #undef __FUNCT__ 606 #define __FUNCT__ "PCPARMSSetFill_PARMS" 607 PetscErrorCode PCPARMSSetFill_PARMS(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 608 { 609 PC_PARMS *parms = (PC_PARMS*)pc->data; 610 611 PetscFunctionBegin; 612 if (lfil0 != parms->lfil[0] || lfil0 != parms->lfil[1] || lfil0 != parms->lfil[2] || lfil0 != parms->lfil[3]) { 613 parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = parms->lfil[0] = lfil0; 614 pc->setupcalled = 0; 615 } 616 if (lfil1 != parms->lfil[4]) { 617 parms->lfil[4] = lfil1; 618 pc->setupcalled = 0; 619 } 620 if (lfil2 != parms->lfil[5] || lfil2 != parms->lfil[6]) { 621 parms->lfil[5] = parms->lfil[6] = lfil2; 622 pc->setupcalled = 0; 623 } 624 PetscFunctionReturn(0); 625 } 626 EXTERN_C_END 627 628 #undef __FUNCT__ 629 #define __FUNCT__ "PCPARMSSetFill" 630 /*@ 631 PCPARMSSetFill - Sets the fill-in parameters for ILUT, ILUK and ARMS preconditioners. 632 Consider the original matrix A = [B F; E C] and the approximate version 633 M = [LB 0; E/UB I]*[UB LB\F; 0 S]. 634 635 Collective on PC 636 637 Input Parameters: 638 + pc - the preconditioner context 639 . fil0 - the level of fill-in kept in LB, UB, E/UB and LB\F 640 . fil1 - the level of fill-in kept in S 641 - fil2 - the level of fill-in kept in the L and U parts of the LU factorization of S 642 643 Options Database Keys: 644 + -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 645 . -pc_parms_lfil_schur - set the amount of fill-in for schur 646 - -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 647 648 Level: intermediate 649 650 Notes: 651 See the pARMS function parms_PCSetFill for more information. 652 653 .seealso: PCPARMS 654 @*/ 655 PetscErrorCode PCPARMSSetFill(PC pc,PetscInt lfil0,PetscInt lfil1,PetscInt lfil2) 656 { 657 PetscErrorCode ierr; 658 659 PetscFunctionBegin; 660 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 661 ierr = PetscTryMethod(pc,"PCPARMSSetFill_C",(PC,PetscInt,PetscInt,PetscInt),(pc,lfil0,lfil1,lfil2));CHKERRQ(ierr); 662 PetscFunctionReturn(0); 663 } 664 665 /*MC 666 PCPARMS - Allows the use of the parallel Algebraic Recursive Multilevel Solvers 667 available in the package pARMS 668 669 Options Database Keys: 670 + -pc_parms_global - one of ras, schur, bj 671 . -pc_parms_local - one of ilu0, iluk, ilut, arms 672 . -pc_parms_solve_tol - set the tolerance for local solve 673 . -pc_parms_levels - set the number of levels 674 . -pc_parms_nonsymmetric_perm - set the use of nonsymmetric permutation 675 . -pc_parms_blocksize - set the block size 676 . -pc_parms_ind_tol - set the tolerance for independent sets 677 . -pc_parms_max_dim - set the inner krylov dimension 678 . -pc_parms_max_it - set the maximum number of inner iterations 679 . -pc_parms_inter_nonsymmetric_perm - set the use of nonsymmetric permutation for interlevel blocks 680 . -pc_parms_inter_column_perm - set the use of column permutation for interlevel blocks 681 . -pc_parms_inter_row_scaling - set the use of row scaling for interlevel blocks 682 . -pc_parms_inter_column_scaling - set the use of column scaling for interlevel blocks 683 . -pc_parms_last_nonsymmetric_perm - set the use of nonsymmetric permutation for last level blocks 684 . -pc_parms_last_column_perm - set the use of column permutation for last level blocks 685 . -pc_parms_last_row_scaling - set the use of row scaling for last level blocks 686 . -pc_parms_last_column_scaling - set the use of column scaling for last level blocks 687 . -pc_parms_lfil_ilu_arms - set the amount of fill-in for ilut, iluk and arms 688 . -pc_parms_lfil_schur - set the amount of fill-in for schur 689 . -pc_parms_lfil_ilut_L_U - set the amount of fill-in for ILUT L and U 690 . -pc_parms_droptol_factors - set the drop tolerance for L, U, L^{-1}F and EU^{-1} 691 . -pc_parms_droptol_schur_compl - set the drop tolerance for schur complement at each level 692 - -pc_parms_droptol_last_schur - set the drop tolerance for ILUT in last level schur complement 693 694 IMPORTANT: 695 Unless configured appropriately, this preconditioner performs an inexact solve 696 as part of the preconditioner application. Therefore, it must be used in combination 697 with flexible variants of iterative solvers, such as KSPFGMRES or KSPCGR. 698 699 Level: intermediate 700 701 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC 702 M*/ 703 704 EXTERN_C_BEGIN 705 #undef __FUNCT__ 706 #define __FUNCT__ "PCCreate_PARMS" 707 PetscErrorCode PCCreate_PARMS(PC pc) 708 { 709 PC_PARMS *parms; 710 PetscErrorCode ierr; 711 712 PetscFunctionBegin; 713 ierr = PetscNewLog(pc,PC_PARMS,&parms);CHKERRQ(ierr); 714 715 parms->map = 0; 716 parms->A = 0; 717 parms->pc = 0; 718 parms->global = PC_PARMS_GLOBAL_RAS; 719 parms->local = PC_PARMS_LOCAL_ARMS; 720 parms->levels = 10; 721 parms->nonsymperm = PETSC_TRUE; 722 parms->blocksize = 250; 723 parms->maxdim = 0; 724 parms->maxits = 0; 725 parms->meth[0] = PETSC_FALSE; 726 parms->meth[1] = PETSC_FALSE; 727 parms->meth[2] = PETSC_FALSE; 728 parms->meth[3] = PETSC_FALSE; 729 parms->meth[4] = PETSC_FALSE; 730 parms->meth[5] = PETSC_FALSE; 731 parms->meth[6] = PETSC_FALSE; 732 parms->meth[7] = PETSC_FALSE; 733 parms->solvetol = 0.01; 734 parms->indtol = 0.4; 735 parms->lfil[0] = parms->lfil[1] = parms->lfil[2] = parms->lfil[3] = 20; 736 parms->lfil[4] = parms->lfil[5] = parms->lfil[6] = 20; 737 parms->droptol[0] = parms->droptol[1] = parms->droptol[2] = parms->droptol[3] = 0.00001; 738 parms->droptol[4] = 0.001; 739 parms->droptol[5] = parms->droptol[6] = 0.001; 740 741 pc->data = parms; 742 pc->ops->destroy = PCDestroy_PARMS; 743 pc->ops->setfromoptions = PCSetFromOptions_PARMS; 744 pc->ops->setup = PCSetUp_PARMS; 745 pc->ops->apply = PCApply_PARMS; 746 pc->ops->view = PCView_PARMS; 747 748 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetGlobal_C","PCPARMSSetGlobal_PARMS",PCPARMSSetGlobal_PARMS);CHKERRQ(ierr); 749 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetLocal_C","PCPARMSSetLocal_PARMS",PCPARMSSetLocal_PARMS);CHKERRQ(ierr); 750 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveTolerances_C","PCPARMSSetSolveTolerances_PARMS",PCPARMSSetSolveTolerances_PARMS);CHKERRQ(ierr); 751 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetSolveRestart_C","PCPARMSSetSolveRestart_PARMS",PCPARMSSetSolveRestart_PARMS);CHKERRQ(ierr); 752 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetNonsymPerm_C","PCPARMSSetNonsymPerm_PARMS",PCPARMSSetNonsymPerm_PARMS);CHKERRQ(ierr); 753 ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCPARMSSetFill_C","PCPARMSSetFill_PARMS",PCPARMSSetFill_PARMS);CHKERRQ(ierr); 754 PetscFunctionReturn(0); 755 } 756 EXTERN_C_END 757