1dba47a55SKris Buschelman 29b54502bSHong Zhang /* 39b54502bSHong Zhang Defines a ILU factorization preconditioner for any Mat implementation 49b54502bSHong Zhang */ 5c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/ilu/ilu.h> /*I "petscpc.h" I*/ 69b54502bSHong Zhang 79b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 8afaefe49SHong Zhang #undef __FUNCT__ 92401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_ILU" 107087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill_ILU(PC pc,PetscBool flag) 112401956bSBarry Smith { 12075768bcSBarry Smith PC_ILU *lu = (PC_ILU*)pc->data; 132401956bSBarry Smith 142401956bSBarry Smith PetscFunctionBegin; 152401956bSBarry Smith lu->reusefill = flag; 162401956bSBarry Smith PetscFunctionReturn(0); 172401956bSBarry Smith } 182401956bSBarry Smith 192401956bSBarry Smith #undef __FUNCT__ 202401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 217087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 229b54502bSHong Zhang { 239b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 249b54502bSHong Zhang 259b54502bSHong Zhang PetscFunctionBegin; 269b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 272fa5cd67SKarl Rupp if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10; 282fa5cd67SKarl Rupp else ilu->nonzerosalongdiagonaltol = z; 299b54502bSHong Zhang PetscFunctionReturn(0); 309b54502bSHong Zhang } 319b54502bSHong Zhang 329b54502bSHong Zhang #undef __FUNCT__ 33574deadeSBarry Smith #define __FUNCT__ "PCReset_ILU" 34574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc) 359b54502bSHong Zhang { 369b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 379b54502bSHong Zhang PetscErrorCode ierr; 389b54502bSHong Zhang 399b54502bSHong Zhang PetscFunctionBegin; 40298cc208SBarry Smith if (!ilu->inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 416bf464f9SBarry Smith if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);} 42fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 439b54502bSHong Zhang PetscFunctionReturn(0); 449b54502bSHong Zhang } 459b54502bSHong Zhang 469b54502bSHong Zhang #undef __FUNCT__ 47b7c853c4SBarry Smith #define __FUNCT__ "PCFactorSetDropTolerance_ILU" 487087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 499b54502bSHong Zhang { 50075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 519b54502bSHong Zhang 529b54502bSHong Zhang PetscFunctionBegin; 534c9036c7SBarry Smith if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 54ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 559b54502bSHong Zhang } 56075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = dt; 57075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = dtcol; 58075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = dtcount; 594c9036c7SBarry Smith ((PC_Factor*)ilu)->info.usedt = 1.0; 609b54502bSHong Zhang PetscFunctionReturn(0); 619b54502bSHong Zhang } 629b54502bSHong Zhang 639b54502bSHong Zhang #undef __FUNCT__ 642401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 657087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseOrdering_ILU(PC pc,PetscBool flag) 669b54502bSHong Zhang { 67075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 689b54502bSHong Zhang 699b54502bSHong Zhang PetscFunctionBegin; 709b54502bSHong Zhang ilu->reuseordering = flag; 719b54502bSHong Zhang PetscFunctionReturn(0); 729b54502bSHong Zhang } 739b54502bSHong Zhang 749b54502bSHong Zhang #undef __FUNCT__ 752401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 767087cfbeSBarry Smith PetscErrorCode PCFactorSetUseInPlace_ILU(PC pc) 779b54502bSHong Zhang { 78075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 799b54502bSHong Zhang 809b54502bSHong Zhang PetscFunctionBegin; 819b54502bSHong Zhang dir->inplace = PETSC_TRUE; 829b54502bSHong Zhang PetscFunctionReturn(0); 839b54502bSHong Zhang } 849b54502bSHong Zhang 859b54502bSHong Zhang #undef __FUNCT__ 869b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 879b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 889b54502bSHong Zhang { 899b54502bSHong Zhang PetscErrorCode ierr; 9078fc6b22SHong Zhang PetscInt itmp; 918afaa268SBarry Smith PetscBool flg,set; 929b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 939b54502bSHong Zhang PetscReal tol; 949b54502bSHong Zhang 959b54502bSHong Zhang PetscFunctionBegin; 969b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 978ff23777SHong Zhang ierr = PCSetFromOptions_Factor(pc);CHKERRQ(ierr); 988ff23777SHong Zhang 99075768bcSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 100075768bcSBarry Smith if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 1012fa5cd67SKarl Rupp 1028afaa268SBarry Smith ierr = PetscOptionsBool("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",((PC_Factor*)ilu)->info.diagonal_fill ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 1038afaa268SBarry Smith if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg; 1042401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 1059b54502bSHong Zhang if (flg) { 1069b54502bSHong Zhang tol = PETSC_DECIDE; 1072401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 1082401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 1099b54502bSHong Zhang } 1109b54502bSHong Zhang 1119b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 1129b54502bSHong Zhang PetscFunctionReturn(0); 1139b54502bSHong Zhang } 1149b54502bSHong Zhang 1159b54502bSHong Zhang #undef __FUNCT__ 1169b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 1179b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 1189b54502bSHong Zhang { 1199b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 1209b54502bSHong Zhang PetscErrorCode ierr; 121ace3abfcSBarry Smith PetscBool iascii; 1229b54502bSHong Zhang 1239b54502bSHong Zhang PetscFunctionBegin; 124251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 1259b54502bSHong Zhang if (iascii) { 126914a5d51SHong Zhang if (ilu->inplace) { 127914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: in-place factorization\n");CHKERRQ(ierr); 1289b54502bSHong Zhang } else { 129914a5d51SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: out-of-place factorization\n");CHKERRQ(ierr); 130145b9266SHong Zhang } 131145b9266SHong Zhang 132914a5d51SHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," ILU: Reusing fill from past factorization\n");CHKERRQ(ierr);} 133914a5d51SHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," ILU: Reusing reordering from past factorization\n");CHKERRQ(ierr);} 1349b54502bSHong Zhang } 135914a5d51SHong Zhang ierr = PCView_Factor(pc,viewer);CHKERRQ(ierr); 1369b54502bSHong Zhang PetscFunctionReturn(0); 1379b54502bSHong Zhang } 1389b54502bSHong Zhang 1399b54502bSHong Zhang #undef __FUNCT__ 1409b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 1419b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 1429b54502bSHong Zhang { 1439b54502bSHong Zhang PetscErrorCode ierr; 1449b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 145f3a39becSBarry Smith MatInfo info; 146ace3abfcSBarry Smith PetscBool flg; 1479b54502bSHong Zhang 1489b54502bSHong Zhang PetscFunctionBegin; 14992927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 15092927226SBarry Smith if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 15122d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 15292927226SBarry Smith if (!flg) { 15322d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr); 15492927226SBarry Smith if (!flg) { 15592927226SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 15622d28d08SBarry Smith PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO");CHKERRQ(ierr); 15792927226SBarry Smith } 15892927226SBarry Smith } 15992927226SBarry Smith } 16092927226SBarry Smith 1619b54502bSHong Zhang if (ilu->inplace) { 1629b54502bSHong Zhang if (!pc->setupcalled) { 1639b54502bSHong Zhang 1649b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 1659b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 166075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1673bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1683bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1699b54502bSHong Zhang } 1709b54502bSHong Zhang 1719b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1729b54502bSHong Zhang cannot have levels of fill */ 173075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 17475567043SBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 1752fa5cd67SKarl Rupp 176365a8a9eSBarry Smith ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr);CHKERRQ(ierr); 177075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 178b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 179b33856dcSBarry Smith ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr); 1809b54502bSHong Zhang } else { 1819b54502bSHong Zhang if (!pc->setupcalled) { 1829b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 183075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1843bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1853bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1869b54502bSHong Zhang /* Remove zeros along diagonal? */ 1879b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1889b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 1899b54502bSHong Zhang } 19086daa42cSBarry Smith if (!((PC_Factor*)ilu)->fact) { 19111bfe483SHong Zhang ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 192a1f19f5aSHong Zhang } 193075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 194075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1952fa5cd67SKarl Rupp 196f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 1972fa5cd67SKarl Rupp 1983bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1999b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 2009b54502bSHong Zhang if (!ilu->reuseordering) { 2019b54502bSHong Zhang /* compute a new ordering for the ILU */ 202fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->row);CHKERRQ(ierr); 203fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 204075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 2053bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 2063bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 2079b54502bSHong Zhang /* Remove zeros along diagonal? */ 2089b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 2099b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 2109b54502bSHong Zhang } 2119b54502bSHong Zhang } 2126bf464f9SBarry Smith ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 21314798fb4SJed Brown ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 214075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 215075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 2162fa5cd67SKarl Rupp 217f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 2182fa5cd67SKarl Rupp 2193bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 2209b54502bSHong Zhang } 221075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 2229b54502bSHong Zhang } 2239b54502bSHong Zhang PetscFunctionReturn(0); 2249b54502bSHong Zhang } 2259b54502bSHong Zhang 2269b54502bSHong Zhang #undef __FUNCT__ 2279b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 2289b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 2299b54502bSHong Zhang { 2309b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2319b54502bSHong Zhang PetscErrorCode ierr; 2329b54502bSHong Zhang 2339b54502bSHong Zhang PetscFunctionBegin; 234574deadeSBarry Smith ierr = PCReset_ILU(pc);CHKERRQ(ierr); 235503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 236503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 237c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2389b54502bSHong Zhang PetscFunctionReturn(0); 2399b54502bSHong Zhang } 2409b54502bSHong Zhang 2419b54502bSHong Zhang #undef __FUNCT__ 2429b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 2439b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 2449b54502bSHong Zhang { 2459b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2469b54502bSHong Zhang PetscErrorCode ierr; 2479b54502bSHong Zhang 2489b54502bSHong Zhang PetscFunctionBegin; 249075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2509b54502bSHong Zhang PetscFunctionReturn(0); 2519b54502bSHong Zhang } 2529b54502bSHong Zhang 2539b54502bSHong Zhang #undef __FUNCT__ 2549b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 2559b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 2569b54502bSHong Zhang { 2579b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2589b54502bSHong Zhang PetscErrorCode ierr; 2599b54502bSHong Zhang 2609b54502bSHong Zhang PetscFunctionBegin; 261075768bcSBarry Smith ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2629b54502bSHong Zhang PetscFunctionReturn(0); 2639b54502bSHong Zhang } 2649b54502bSHong Zhang 265f0b9ad6cSBarry Smith #undef __FUNCT__ 266f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricLeft_ILU" 267f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 268f0b9ad6cSBarry Smith { 269f0b9ad6cSBarry Smith PetscErrorCode ierr; 270f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 271f0b9ad6cSBarry Smith 272f0b9ad6cSBarry Smith PetscFunctionBegin; 273f0b9ad6cSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 274f0b9ad6cSBarry Smith PetscFunctionReturn(0); 275f0b9ad6cSBarry Smith } 276f0b9ad6cSBarry Smith 277f0b9ad6cSBarry Smith #undef __FUNCT__ 278f0b9ad6cSBarry Smith #define __FUNCT__ "PCApplySymmetricRight_ILU" 279f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 280f0b9ad6cSBarry Smith { 281f0b9ad6cSBarry Smith PetscErrorCode ierr; 282f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 283f0b9ad6cSBarry Smith 284f0b9ad6cSBarry Smith PetscFunctionBegin; 285f0b9ad6cSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 286f0b9ad6cSBarry Smith PetscFunctionReturn(0); 287f0b9ad6cSBarry Smith } 288f0b9ad6cSBarry Smith 2899b54502bSHong Zhang /*MC 2909b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2919b54502bSHong Zhang 2929b54502bSHong Zhang Options Database Keys: 2932401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2942401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2959b54502bSHong Zhang its factorization (overwrites original matrix) 2962401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2972401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 29855ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2992401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 3009b54502bSHong Zhang this decreases the chance of getting a zero pivot 3012401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 3022401956bSBarry Smith . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 3039b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 3049b54502bSHong Zhang stability of the ILU factorization 3059b54502bSHong Zhang 3069b54502bSHong Zhang Level: beginner 3079b54502bSHong Zhang 3089b54502bSHong Zhang Concepts: incomplete factorization 3099b54502bSHong Zhang 310d9ba8547SBarry Smith Notes: Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 3119b54502bSHong Zhang 3129b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 3139b54502bSHong Zhang 314f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 315f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 316f0b9ad6cSBarry Smith 317cd0a26f6SPaul Mullowney If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARESE matrices with block Jacobi), factorization 318cd0a26f6SPaul Mullowney is never done on the GPU). 319ac793be5SBarry Smith 320c582cd25SBarry Smith References: 321c582cd25SBarry Smith T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 322c582cd25SBarry Smith self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 323c582cd25SBarry Smith 324c582cd25SBarry Smith T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 325c582cd25SBarry Smith fusion problems. Quart. Appl. Math., 19:221--229, 1961. 326c582cd25SBarry Smith 327c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 328c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 329c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 330c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 331c582cd25SBarry Smith 332c582cd25SBarry Smith 3339b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 334145b9266SHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 335b7c853c4SBarry Smith PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 3368ff23777SHong Zhang PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks() 3379b54502bSHong Zhang 3389b54502bSHong Zhang M*/ 3399b54502bSHong Zhang 3409b54502bSHong Zhang #undef __FUNCT__ 3419b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 3428cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 3439b54502bSHong Zhang { 3449b54502bSHong Zhang PetscErrorCode ierr; 3459b54502bSHong Zhang PC_ILU *ilu; 3469b54502bSHong Zhang 3479b54502bSHong Zhang PetscFunctionBegin; 348b00a9115SJed Brown ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr); 3499b54502bSHong Zhang 350075768bcSBarry Smith ((PC_Factor*)ilu)->fact = 0; 351075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 352879e8a4dSBarry Smith ((PC_Factor*)ilu)->factortype = MAT_FACTOR_ILU; 35375567043SBarry Smith ((PC_Factor*)ilu)->info.levels = 0.; 354075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 3559b54502bSHong Zhang ilu->col = 0; 3569b54502bSHong Zhang ilu->row = 0; 3579b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 3582692d6eeSBarry Smith ierr = PetscStrallocpy(MATSOLVERPETSC,&((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 35919fd82e9SBarry Smith ierr = PetscStrallocpy(MATORDERINGNATURAL,(char**)&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 3609b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 361075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 362075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 363075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 364*c378d667SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 3650ed735ceSHong Zhang ((PC_Factor*)ilu)->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 3660ed735ceSHong Zhang ((PC_Factor*)ilu)->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 367075768bcSBarry Smith ((PC_Factor*)ilu)->info.pivotinblocks = 1.0; 3689b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 3690ed735ceSHong Zhang ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 3709b54502bSHong Zhang pc->data = (void*)ilu; 3719b54502bSHong Zhang 372574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3739b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3749b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3759b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3769b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3779b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 37885317021SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 3799b54502bSHong Zhang pc->ops->view = PCView_ILU; 380f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 381f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3829b54502bSHong Zhang pc->ops->applyrichardson = 0; 3839b54502bSHong Zhang 384bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor);CHKERRQ(ierr); 385bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor);CHKERRQ(ierr); 386bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor);CHKERRQ(ierr); 387bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverPackage_C",PCFactorGetMatSolverPackage_Factor);CHKERRQ(ierr); 388bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverPackage_C",PCFactorSetMatSolverPackage_Factor);CHKERRQ(ierr); 389bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverPackage_C",PCFactorSetUpMatSolverPackage_Factor);CHKERRQ(ierr); 390bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 391bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor);CHKERRQ(ierr); 392bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor);CHKERRQ(ierr); 393bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 394bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 395bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor);CHKERRQ(ierr); 3962591b318SToby Isaac ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor);CHKERRQ(ierr); 397bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 398bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor);CHKERRQ(ierr); 399bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor);CHKERRQ(ierr); 400bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 4019b54502bSHong Zhang PetscFunctionReturn(0); 4029b54502bSHong Zhang } 403