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 77087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 89b54502bSHong Zhang { 99b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 109b54502bSHong Zhang 119b54502bSHong Zhang PetscFunctionBegin; 129b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 132fa5cd67SKarl Rupp if (z == PETSC_DECIDE) ilu->nonzerosalongdiagonaltol = 1.e-10; 142fa5cd67SKarl Rupp else ilu->nonzerosalongdiagonaltol = z; 159b54502bSHong Zhang PetscFunctionReturn(0); 169b54502bSHong Zhang } 179b54502bSHong Zhang 18574deadeSBarry Smith PetscErrorCode PCReset_ILU(PC pc) 199b54502bSHong Zhang { 209b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 219b54502bSHong Zhang PetscErrorCode ierr; 229b54502bSHong Zhang 239b54502bSHong Zhang PetscFunctionBegin; 243d1c1ea0SBarry Smith if (!ilu->hdr.inplace) {ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 256bf464f9SBarry Smith if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(&ilu->row);CHKERRQ(ierr);} 26fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 279b54502bSHong Zhang PetscFunctionReturn(0); 289b54502bSHong Zhang } 299b54502bSHong Zhang 307087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 319b54502bSHong Zhang { 32075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 339b54502bSHong Zhang 349b54502bSHong Zhang PetscFunctionBegin; 354c9036c7SBarry Smith if (pc->setupcalled && (((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 36ce94432eSBarry Smith SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot change drop tolerance after using PC"); 379b54502bSHong Zhang } 38075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = dt; 39075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = dtcol; 40075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = dtcount; 414c9036c7SBarry Smith ((PC_Factor*)ilu)->info.usedt = 1.0; 429b54502bSHong Zhang PetscFunctionReturn(0); 439b54502bSHong Zhang } 449b54502bSHong Zhang 454416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ILU(PetscOptionItems *PetscOptionsObject,PC pc) 469b54502bSHong Zhang { 479b54502bSHong Zhang PetscErrorCode ierr; 4878fc6b22SHong Zhang PetscInt itmp; 498afaa268SBarry Smith PetscBool flg,set; 509b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 519b54502bSHong Zhang PetscReal tol; 529b54502bSHong Zhang 539b54502bSHong Zhang PetscFunctionBegin; 54e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"ILU Options");CHKERRQ(ierr); 55e55864a3SBarry Smith ierr = PCSetFromOptions_Factor(PetscOptionsObject,pc);CHKERRQ(ierr); 568ff23777SHong Zhang 57075768bcSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 58075768bcSBarry Smith if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 592fa5cd67SKarl Rupp 608afaa268SBarry 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); 618afaa268SBarry Smith if (set) ((PC_Factor*)ilu)->info.diagonal_fill = (PetscReal) flg; 622401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 639b54502bSHong Zhang if (flg) { 649b54502bSHong Zhang tol = PETSC_DECIDE; 650a545947SLisandro Dalcin ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,NULL);CHKERRQ(ierr); 662401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 679b54502bSHong Zhang } 689b54502bSHong Zhang 699b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 709b54502bSHong Zhang PetscFunctionReturn(0); 719b54502bSHong Zhang } 729b54502bSHong Zhang 739b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 749b54502bSHong Zhang { 759b54502bSHong Zhang PetscErrorCode ierr; 769b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 77f3a39becSBarry Smith MatInfo info; 78ace3abfcSBarry Smith PetscBool flg; 79ea799195SBarry Smith MatSolverType stype; 8000e125f8SBarry Smith MatFactorError err; 819b54502bSHong Zhang 829b54502bSHong Zhang PetscFunctionBegin; 83c6e4fdc6SHong Zhang pc->failedreason = PC_NOERROR; 8492927226SBarry Smith /* ugly hack to change default, since it is not support by some matrix types */ 8592927226SBarry Smith if (((PC_Factor*)ilu)->info.shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 8622d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 8792927226SBarry Smith if (!flg) { 8822d28d08SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATMPIAIJ,&flg);CHKERRQ(ierr); 8992927226SBarry Smith if (!flg) { 9092927226SBarry Smith ((PC_Factor*)ilu)->info.shifttype = (PetscReal)MAT_SHIFT_INBLOCKS; 91994fe344SLisandro Dalcin ierr = PetscInfo(pc,"Changing shift type from NONZERO to INBLOCKS because block matrices do not support NONZERO\n");CHKERRQ(ierr); 9292927226SBarry Smith } 9392927226SBarry Smith } 9492927226SBarry Smith } 9592927226SBarry Smith 9684d44b13SHong Zhang ierr = MatSetErrorIfFailure(pc->pmat,pc->erroriffailure);CHKERRQ(ierr); 973d1c1ea0SBarry Smith if (ilu->hdr.inplace) { 989b54502bSHong Zhang if (!pc->setupcalled) { 999b54502bSHong Zhang 1009b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 1019b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 1022c7c0729SBarry Smith /* Should not get the ordering if the factorization routine does not use it, but do not yet have access to the factor matrix */ 1034ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 1044ac6704cSBarry Smith ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 105075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1063bb1ff40SBarry Smith if (ilu->row) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr);} 1073bb1ff40SBarry Smith if (ilu->col) {ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr);} 1089b54502bSHong Zhang } 1099b54502bSHong Zhang 1109b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 1119b54502bSHong Zhang cannot have levels of fill */ 112075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 11375567043SBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0.0; 1142fa5cd67SKarl Rupp 115c3b366b1Sprj- ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 11600e125f8SBarry Smith ierr = MatFactorGetError(pc->pmat,&err);CHKERRQ(ierr); 11700e125f8SBarry Smith if (err) { /* Factor() fails */ 11800e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1196baea169SHong Zhang PetscFunctionReturn(0); 1206baea169SHong Zhang } 1216baea169SHong Zhang 122075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 123b33856dcSBarry Smith /* must update the pc record of the matrix state or the PC will attempt to run PCSetUp() yet again */ 124b33856dcSBarry Smith ierr = PetscObjectStateGet((PetscObject)pc->pmat,&pc->matstate);CHKERRQ(ierr); 1259b54502bSHong Zhang } else { 1269b54502bSHong Zhang if (!pc->setupcalled) { 1279b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 128f73b0415SBarry Smith PetscBool canuseordering; 1292c7c0729SBarry Smith if (!((PC_Factor*)ilu)->fact) { 1302c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1312c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1322c7c0729SBarry Smith } 133f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering);CHKERRQ(ierr); 134f73b0415SBarry Smith if (canuseordering) { 1354ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 136075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1372c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr); 1382c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr); 1399b54502bSHong Zhang /* Remove zeros along diagonal? */ 1409b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1419b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 1429b54502bSHong Zhang } 143a1f19f5aSHong Zhang } 144075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 145075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1463d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1479b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 1483d1c1ea0SBarry Smith if (!ilu->hdr.reuseordering) { 149f73b0415SBarry Smith PetscBool canuseordering; 1502c7c0729SBarry Smith ierr = MatDestroy(&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1512c7c0729SBarry Smith ierr = MatGetFactor(pc->pmat,((PC_Factor*)ilu)->solvertype,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 1522c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 153f73b0415SBarry Smith ierr = MatFactorGetCanUseOrdering(((PC_Factor*)ilu)->fact,&canuseordering);CHKERRQ(ierr); 154f73b0415SBarry Smith if (canuseordering) { 1559b54502bSHong Zhang /* compute a new ordering for the ILU */ 156fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->row);CHKERRQ(ierr); 157fcfd50ebSBarry Smith ierr = ISDestroy(&ilu->col);CHKERRQ(ierr); 1584ac6704cSBarry Smith ierr = PCFactorSetDefaultOrdering_Factor(pc);CHKERRQ(ierr); 159075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 1602c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->row);CHKERRQ(ierr); 1612c7c0729SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilu->col);CHKERRQ(ierr); 1629b54502bSHong Zhang /* Remove zeros along diagonal? */ 1639b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 1649b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 1659b54502bSHong Zhang } 1669b54502bSHong Zhang } 1672c7c0729SBarry Smith } 168075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 169075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 1703d1c1ea0SBarry Smith ilu->hdr.actualfill = info.fill_ratio_needed; 1719b54502bSHong Zhang } 17200e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 17300e125f8SBarry Smith if (err) { /* FactorSymbolic() fails */ 17400e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1756baea169SHong Zhang PetscFunctionReturn(0); 1766baea169SHong Zhang } 1776baea169SHong Zhang 178075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 17900e125f8SBarry Smith ierr = MatFactorGetError(((PC_Factor*)ilu)->fact,&err);CHKERRQ(ierr); 18000e125f8SBarry Smith if (err) { /* FactorNumeric() fails */ 18100e125f8SBarry Smith pc->failedreason = (PCFailedReason)err; 1826baea169SHong Zhang } 1839b54502bSHong Zhang } 18400c67f3bSHong Zhang 1853ca39a21SBarry Smith ierr = PCFactorGetMatSolverType(pc,&stype);CHKERRQ(ierr); 18600c67f3bSHong Zhang if (!stype) { 187ea799195SBarry Smith MatSolverType solverpackage; 1883ca39a21SBarry Smith ierr = MatFactorGetSolverType(((PC_Factor*)ilu)->fact,&solverpackage);CHKERRQ(ierr); 1893ca39a21SBarry Smith ierr = PCFactorSetMatSolverType(pc,solverpackage);CHKERRQ(ierr); 19000c67f3bSHong Zhang } 1919b54502bSHong Zhang PetscFunctionReturn(0); 1929b54502bSHong Zhang } 1939b54502bSHong Zhang 1949b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 1959b54502bSHong Zhang { 1969b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 1979b54502bSHong Zhang PetscErrorCode ierr; 1989b54502bSHong Zhang 1999b54502bSHong Zhang PetscFunctionBegin; 200574deadeSBarry Smith ierr = PCReset_ILU(pc);CHKERRQ(ierr); 201503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->solvertype);CHKERRQ(ierr); 202503cfb0cSBarry Smith ierr = PetscFree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 203c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2049b54502bSHong Zhang PetscFunctionReturn(0); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang 2079b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 2089b54502bSHong Zhang { 2099b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2109b54502bSHong Zhang PetscErrorCode ierr; 2119b54502bSHong Zhang 2129b54502bSHong Zhang PetscFunctionBegin; 213075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2149b54502bSHong Zhang PetscFunctionReturn(0); 2159b54502bSHong Zhang } 2169b54502bSHong Zhang 2177b6e2003SPierre Jolivet static PetscErrorCode PCMatApply_ILU(PC pc,Mat X,Mat Y) 2187b6e2003SPierre Jolivet { 2197b6e2003SPierre Jolivet PC_ILU *ilu = (PC_ILU*)pc->data; 2207b6e2003SPierre Jolivet PetscErrorCode ierr; 2217b6e2003SPierre Jolivet 2227b6e2003SPierre Jolivet PetscFunctionBegin; 2237b6e2003SPierre Jolivet ierr = MatMatSolve(((PC_Factor*)ilu)->fact,X,Y);CHKERRQ(ierr); 2247b6e2003SPierre Jolivet PetscFunctionReturn(0); 2257b6e2003SPierre Jolivet } 2267b6e2003SPierre Jolivet 2279b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 2289b54502bSHong Zhang { 2299b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 2309b54502bSHong Zhang PetscErrorCode ierr; 2319b54502bSHong Zhang 2329b54502bSHong Zhang PetscFunctionBegin; 233075768bcSBarry Smith ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 2349b54502bSHong Zhang PetscFunctionReturn(0); 2359b54502bSHong Zhang } 2369b54502bSHong Zhang 237f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricLeft_ILU(PC pc,Vec x,Vec y) 238f0b9ad6cSBarry Smith { 239f0b9ad6cSBarry Smith PetscErrorCode ierr; 240f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 241f0b9ad6cSBarry Smith 242f0b9ad6cSBarry Smith PetscFunctionBegin; 243f0b9ad6cSBarry Smith ierr = MatForwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 244f0b9ad6cSBarry Smith PetscFunctionReturn(0); 245f0b9ad6cSBarry Smith } 246f0b9ad6cSBarry Smith 247f0b9ad6cSBarry Smith static PetscErrorCode PCApplySymmetricRight_ILU(PC pc,Vec x,Vec y) 248f0b9ad6cSBarry Smith { 249f0b9ad6cSBarry Smith PetscErrorCode ierr; 250f0b9ad6cSBarry Smith PC_ILU *icc = (PC_ILU*)pc->data; 251f0b9ad6cSBarry Smith 252f0b9ad6cSBarry Smith PetscFunctionBegin; 253f0b9ad6cSBarry Smith ierr = MatBackwardSolve(((PC_Factor*)icc)->fact,x,y);CHKERRQ(ierr); 254f0b9ad6cSBarry Smith PetscFunctionReturn(0); 255f0b9ad6cSBarry Smith } 256f0b9ad6cSBarry Smith 2579b54502bSHong Zhang /*MC 2589b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 2599b54502bSHong Zhang 2609b54502bSHong Zhang Options Database Keys: 2612401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 2622401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 2639b54502bSHong Zhang its factorization (overwrites original matrix) 2642401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 2652401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 26655ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 2672401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 2689b54502bSHong Zhang this decreases the chance of getting a zero pivot 2692401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 270967c93d3SBarry Smith - -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 2719b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 2729b54502bSHong Zhang stability of the ILU factorization 2739b54502bSHong Zhang 2749b54502bSHong Zhang Level: beginner 2759b54502bSHong Zhang 27695452b02SPatrick Sanan Notes: 27795452b02SPatrick Sanan Only implemented for some matrix formats. (for parallel see PCHYPRE for hypre's ILU) 2789b54502bSHong Zhang 2799b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 2809b54502bSHong Zhang 281f0b9ad6cSBarry Smith The "symmetric" application of this preconditioner is not actually symmetric since L is not transpose(U) 282f0b9ad6cSBarry Smith even when the matrix is not symmetric since the U stores the diagonals of the factorization. 283f0b9ad6cSBarry Smith 2845ea7661aSPierre Jolivet If you are using MATSEQAIJCUSPARSE matrices (or MATMPIAIJCUSPARSE matrices with block Jacobi), factorization 285cd0a26f6SPaul Mullowney is never done on the GPU). 286ac793be5SBarry Smith 287c582cd25SBarry Smith References: 288*606c0280SSatish Balay + * - T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 28996a0c994SBarry Smith self adjoint elliptic difference equations. SIAM J. Numer. Anal., 5, 1968. 290*606c0280SSatish Balay . * - T.A. Oliphant. An implicit numerical method for solving two dimensional timedependent diffusion problems. Quart. Appl. Math., 19, 1961. 291*606c0280SSatish Balay - * - TONY F. CHAN AND HENK A. VAN DER VORST, APPROXIMATE AND INCOMPLETE FACTORIZATIONS, 29296a0c994SBarry Smith Chapter in Parallel Numerical 293c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 29496a0c994SBarry Smith Science and Engineering, Kluwer. 295c582cd25SBarry Smith 2969b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 297145b9266SHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftSetType(), PCFactorSetAmount(), 298b7c853c4SBarry Smith PCFactorSetDropTolerance(),PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 29992e9c092SBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 30092e9c092SBarry Smith PCFactorGetAllowDiagonalFill(), PCFactorGetUseInPlace() 3019b54502bSHong Zhang 3029b54502bSHong Zhang M*/ 3039b54502bSHong Zhang 3048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ILU(PC pc) 3059b54502bSHong Zhang { 3069b54502bSHong Zhang PetscErrorCode ierr; 3079b54502bSHong Zhang PC_ILU *ilu; 3089b54502bSHong Zhang 3099b54502bSHong Zhang PetscFunctionBegin; 310b00a9115SJed Brown ierr = PetscNewLog(pc,&ilu);CHKERRQ(ierr); 3113d1c1ea0SBarry Smith pc->data = (void*)ilu; 3124ac6704cSBarry Smith ierr = PCFactorInitialize(pc,MAT_FACTOR_ILU);CHKERRQ(ierr); 3139b54502bSHong Zhang 31475567043SBarry Smith ((PC_Factor*)ilu)->info.levels = 0.; 315075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 3160a545947SLisandro Dalcin ilu->col = NULL; 3170a545947SLisandro Dalcin ilu->row = NULL; 318075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 319075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 320075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 3219b54502bSHong Zhang 322574deadeSBarry Smith pc->ops->reset = PCReset_ILU; 3239b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 3249b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 3257b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_ILU; 3269b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 3279b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 3289b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 32992e08861SBarry Smith pc->ops->view = PCView_Factor; 330f0b9ad6cSBarry Smith pc->ops->applysymmetricleft = PCApplySymmetricLeft_ILU; 331f0b9ad6cSBarry Smith pc->ops->applysymmetricright = PCApplySymmetricRight_ILU; 3320a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 333bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",PCFactorSetDropTolerance_ILU);CHKERRQ(ierr); 334bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 3359b54502bSHong Zhang PetscFunctionReturn(0); 3369b54502bSHong Zhang } 337