1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 39b54502bSHong Zhang /* 49b54502bSHong Zhang Defines a ILU factorization preconditioner for any Mat implementation 59b54502bSHong Zhang */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 785b398ccSKris Buschelman #include "src/ksp/pc/impls/factor/ilu/ilu.h" 8b9147fbbSdalcinl #include "include/private/matimpl.h" 99b54502bSHong Zhang 109b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 11afaefe49SHong Zhang EXTERN_C_BEGIN 12afaefe49SHong Zhang #undef __FUNCT__ 132401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_ILU" 142401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag) 152401956bSBarry Smith { 162401956bSBarry Smith PC_ILU *lu; 172401956bSBarry Smith 182401956bSBarry Smith PetscFunctionBegin; 192401956bSBarry Smith lu = (PC_ILU*)pc->data; 202401956bSBarry Smith lu->reusefill = flag; 212401956bSBarry Smith PetscFunctionReturn(0); 222401956bSBarry Smith } 232401956bSBarry Smith EXTERN_C_END 242401956bSBarry Smith 252401956bSBarry Smith EXTERN_C_BEGIN 262401956bSBarry Smith #undef __FUNCT__ 27afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ILU" 28dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ILU(PC pc,PetscReal z) 29afaefe49SHong Zhang { 30afaefe49SHong Zhang PC_ILU *ilu; 31afaefe49SHong Zhang 32afaefe49SHong Zhang PetscFunctionBegin; 33afaefe49SHong Zhang ilu = (PC_ILU*)pc->data; 34afaefe49SHong Zhang ilu->info.zeropivot = z; 35afaefe49SHong Zhang PetscFunctionReturn(0); 36afaefe49SHong Zhang } 37afaefe49SHong Zhang EXTERN_C_END 38afaefe49SHong Zhang 39afaefe49SHong Zhang EXTERN_C_BEGIN 40afaefe49SHong Zhang #undef __FUNCT__ 41afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ILU" 42dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ILU(PC pc,PetscReal shift) 43afaefe49SHong Zhang { 44afaefe49SHong Zhang PC_ILU *dir; 45afaefe49SHong Zhang 46afaefe49SHong Zhang PetscFunctionBegin; 47afaefe49SHong Zhang dir = (PC_ILU*)pc->data; 48afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 49afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 50afaefe49SHong Zhang } else { 51afaefe49SHong Zhang dir->info.shiftnz = shift; 52afaefe49SHong Zhang } 53afaefe49SHong Zhang PetscFunctionReturn(0); 54afaefe49SHong Zhang } 55afaefe49SHong Zhang EXTERN_C_END 56afaefe49SHong Zhang 57afaefe49SHong Zhang EXTERN_C_BEGIN 58afaefe49SHong Zhang #undef __FUNCT__ 59afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ILU" 60dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ILU(PC pc,PetscTruth shift) 61afaefe49SHong Zhang { 62afaefe49SHong Zhang PC_ILU *dir; 63afaefe49SHong Zhang 64afaefe49SHong Zhang PetscFunctionBegin; 65afaefe49SHong Zhang dir = (PC_ILU*)pc->data; 66fbf22428SSatish Balay if (shift) { 67fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 68fbf22428SSatish Balay dir->info.shiftpd = 1.0; 69fbf22428SSatish Balay } else { 70fbf22428SSatish Balay dir->info.shiftpd = 0.0; 71fbf22428SSatish Balay } 72afaefe49SHong Zhang PetscFunctionReturn(0); 73afaefe49SHong Zhang } 74afaefe49SHong Zhang EXTERN_C_END 759b54502bSHong Zhang 769b54502bSHong Zhang EXTERN_C_BEGIN 779b54502bSHong Zhang #undef __FUNCT__ 782401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 792401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 809b54502bSHong Zhang { 819b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 829b54502bSHong Zhang 839b54502bSHong Zhang PetscFunctionBegin; 849b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 859b54502bSHong Zhang if (z == PETSC_DECIDE) { 869b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = 1.e-10; 879b54502bSHong Zhang } else { 889b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = z; 899b54502bSHong Zhang } 909b54502bSHong Zhang PetscFunctionReturn(0); 919b54502bSHong Zhang } 929b54502bSHong Zhang EXTERN_C_END 939b54502bSHong Zhang 949b54502bSHong Zhang #undef __FUNCT__ 959b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal" 969b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc) 979b54502bSHong Zhang { 989b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 999b54502bSHong Zhang PetscErrorCode ierr; 1009b54502bSHong Zhang 1019b54502bSHong Zhang PetscFunctionBegin; 1029b54502bSHong Zhang if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);} 1039b54502bSHong Zhang if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 1049b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 1059b54502bSHong Zhang PetscFunctionReturn(0); 1069b54502bSHong Zhang } 1079b54502bSHong Zhang 1089b54502bSHong Zhang EXTERN_C_BEGIN 1099b54502bSHong Zhang #undef __FUNCT__ 1102401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU" 1112401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 1129b54502bSHong Zhang { 1139b54502bSHong Zhang PC_ILU *ilu; 1149b54502bSHong Zhang PetscErrorCode ierr; 1159b54502bSHong Zhang 1169b54502bSHong Zhang PetscFunctionBegin; 1179b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1189b54502bSHong Zhang if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) { 1199b54502bSHong Zhang pc->setupcalled = 0; 1209b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1219b54502bSHong Zhang } 1229b54502bSHong Zhang ilu->usedt = PETSC_TRUE; 1239b54502bSHong Zhang ilu->info.dt = dt; 1249b54502bSHong Zhang ilu->info.dtcol = dtcol; 1259b54502bSHong Zhang ilu->info.dtcount = dtcount; 1269b54502bSHong Zhang ilu->info.fill = PETSC_DEFAULT; 1279b54502bSHong Zhang PetscFunctionReturn(0); 1289b54502bSHong Zhang } 1299b54502bSHong Zhang EXTERN_C_END 1309b54502bSHong Zhang 1319b54502bSHong Zhang EXTERN_C_BEGIN 1329b54502bSHong Zhang #undef __FUNCT__ 13355ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ILU" 13455ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ILU(PC pc,PetscReal fill) 1359b54502bSHong Zhang { 1369b54502bSHong Zhang PC_ILU *dir; 1379b54502bSHong Zhang 1389b54502bSHong Zhang PetscFunctionBegin; 1399b54502bSHong Zhang dir = (PC_ILU*)pc->data; 1409b54502bSHong Zhang dir->info.fill = fill; 1419b54502bSHong Zhang PetscFunctionReturn(0); 1429b54502bSHong Zhang } 1439b54502bSHong Zhang EXTERN_C_END 1449b54502bSHong Zhang 1459b54502bSHong Zhang EXTERN_C_BEGIN 1469b54502bSHong Zhang #undef __FUNCT__ 147*e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ILU" 148*e5a9bf91SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ILU(PC pc,MatOrderingType ordering) 1499b54502bSHong Zhang { 1509b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 1519b54502bSHong Zhang PetscErrorCode ierr; 1529b54502bSHong Zhang PetscTruth flg; 1539b54502bSHong Zhang 1549b54502bSHong Zhang PetscFunctionBegin; 1559b54502bSHong Zhang if (!pc->setupcalled) { 1569b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1579b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1589b54502bSHong Zhang } else { 1599b54502bSHong Zhang ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 1609b54502bSHong Zhang if (!flg) { 1619b54502bSHong Zhang pc->setupcalled = 0; 1629b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1639b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1649b54502bSHong Zhang /* free the data structures, then create them again */ 1659b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1669b54502bSHong Zhang } 1679b54502bSHong Zhang } 1689b54502bSHong Zhang PetscFunctionReturn(0); 1699b54502bSHong Zhang } 1709b54502bSHong Zhang EXTERN_C_END 1719b54502bSHong Zhang 1729b54502bSHong Zhang EXTERN_C_BEGIN 1739b54502bSHong Zhang #undef __FUNCT__ 1742401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 1752401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag) 1769b54502bSHong Zhang { 1779b54502bSHong Zhang PC_ILU *ilu; 1789b54502bSHong Zhang 1799b54502bSHong Zhang PetscFunctionBegin; 1809b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1819b54502bSHong Zhang ilu->reuseordering = flag; 1829b54502bSHong Zhang PetscFunctionReturn(0); 1839b54502bSHong Zhang } 1849b54502bSHong Zhang EXTERN_C_END 1859b54502bSHong Zhang 1869b54502bSHong Zhang EXTERN_C_BEGIN 1879b54502bSHong Zhang #undef __FUNCT__ 1882401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels_ILU" 1892401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ILU(PC pc,PetscInt levels) 1909b54502bSHong Zhang { 1919b54502bSHong Zhang PC_ILU *ilu; 1929b54502bSHong Zhang PetscErrorCode ierr; 1939b54502bSHong Zhang 1949b54502bSHong Zhang PetscFunctionBegin; 1959b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1969b54502bSHong Zhang 1979b54502bSHong Zhang if (!pc->setupcalled) { 1989b54502bSHong Zhang ilu->info.levels = levels; 1999b54502bSHong Zhang } else if (ilu->usedt || ilu->info.levels != levels) { 2009b54502bSHong Zhang ilu->info.levels = levels; 2019b54502bSHong Zhang pc->setupcalled = 0; 2029b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 2039b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 2049b54502bSHong Zhang } 2059b54502bSHong Zhang PetscFunctionReturn(0); 2069b54502bSHong Zhang } 2079b54502bSHong Zhang EXTERN_C_END 2089b54502bSHong Zhang 2099b54502bSHong Zhang EXTERN_C_BEGIN 2109b54502bSHong Zhang #undef __FUNCT__ 2112401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 2122401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc) 2139b54502bSHong Zhang { 2149b54502bSHong Zhang PC_ILU *dir; 2159b54502bSHong Zhang 2169b54502bSHong Zhang PetscFunctionBegin; 2179b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2189b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2199b54502bSHong Zhang PetscFunctionReturn(0); 2209b54502bSHong Zhang } 2219b54502bSHong Zhang EXTERN_C_END 2229b54502bSHong Zhang 2239b54502bSHong Zhang EXTERN_C_BEGIN 2249b54502bSHong Zhang #undef __FUNCT__ 22512a54f56SSatish Balay #define __FUNCT__ "PCFactorSetAllowDiagonalFill_ILU" 2262401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_ILU(PC pc) 2279b54502bSHong Zhang { 2289b54502bSHong Zhang PC_ILU *dir; 2299b54502bSHong Zhang 2309b54502bSHong Zhang PetscFunctionBegin; 2319b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2329b54502bSHong Zhang dir->info.diagonal_fill = 1; 2339b54502bSHong Zhang PetscFunctionReturn(0); 2349b54502bSHong Zhang } 2359b54502bSHong Zhang EXTERN_C_END 2369b54502bSHong Zhang 2379b54502bSHong Zhang #undef __FUNCT__ 2382401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance" 2399b54502bSHong Zhang /*@ 2402401956bSBarry Smith PCFactorSetUseDropTolerance - The preconditioner will use an ILU 2419b54502bSHong Zhang based on a drop tolerance. 2429b54502bSHong Zhang 2439b54502bSHong Zhang Collective on PC 2449b54502bSHong Zhang 2459b54502bSHong Zhang Input Parameters: 2469b54502bSHong Zhang + pc - the preconditioner context 2479b54502bSHong Zhang . dt - the drop tolerance, try from 1.e-10 to .1 2489b54502bSHong Zhang . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 2499b54502bSHong Zhang - maxrowcount - the max number of nonzeros allowed in a row, best value 2509b54502bSHong Zhang depends on the number of nonzeros in row of original matrix 2519b54502bSHong Zhang 2529b54502bSHong Zhang Options Database Key: 2532401956bSBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 2549b54502bSHong Zhang 2559b54502bSHong Zhang Level: intermediate 2569b54502bSHong Zhang 2579b54502bSHong Zhang Notes: 2589b54502bSHong Zhang This uses the iludt() code of Saad's SPARSKIT package 2599b54502bSHong Zhang 2609b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 2619b54502bSHong Zhang @*/ 2622401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 2639b54502bSHong Zhang { 2649b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 2659b54502bSHong Zhang 2669b54502bSHong Zhang PetscFunctionBegin; 2679b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2682401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 2699b54502bSHong Zhang if (f) { 2709b54502bSHong Zhang ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 2719b54502bSHong Zhang } 2729b54502bSHong Zhang PetscFunctionReturn(0); 2739b54502bSHong Zhang } 2749b54502bSHong Zhang 2759b54502bSHong Zhang #undef __FUNCT__ 2762401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels" 2779b54502bSHong Zhang /*@ 2782401956bSBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 2799b54502bSHong Zhang 2809b54502bSHong Zhang Collective on PC 2819b54502bSHong Zhang 2829b54502bSHong Zhang Input Parameters: 2839b54502bSHong Zhang + pc - the preconditioner context 2849b54502bSHong Zhang - levels - number of levels of fill 2859b54502bSHong Zhang 2869b54502bSHong Zhang Options Database Key: 2872401956bSBarry Smith . -pc_factor_levels <levels> - Sets fill level 2889b54502bSHong Zhang 2899b54502bSHong Zhang Level: intermediate 2909b54502bSHong Zhang 2919b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 2929b54502bSHong Zhang @*/ 2932401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 2949b54502bSHong Zhang { 2959b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 2969b54502bSHong Zhang 2979b54502bSHong Zhang PetscFunctionBegin; 2989b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2999b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 3002401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 3019b54502bSHong Zhang if (f) { 3029b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 3039b54502bSHong Zhang } 3049b54502bSHong Zhang PetscFunctionReturn(0); 3059b54502bSHong Zhang } 3069b54502bSHong Zhang 3079b54502bSHong Zhang #undef __FUNCT__ 3082401956bSBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 3099b54502bSHong Zhang /*@ 3102401956bSBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 3119b54502bSHong Zhang treated as level 0 fill even if there is no non-zero location. 3129b54502bSHong Zhang 3139b54502bSHong Zhang Collective on PC 3149b54502bSHong Zhang 3159b54502bSHong Zhang Input Parameters: 3169b54502bSHong Zhang + pc - the preconditioner context 3179b54502bSHong Zhang 3189b54502bSHong Zhang Options Database Key: 3192401956bSBarry Smith . -pc_factor_diagonal_fill 3209b54502bSHong Zhang 3219b54502bSHong Zhang Notes: 3229b54502bSHong Zhang Does not apply with 0 fill. 3239b54502bSHong Zhang 3249b54502bSHong Zhang Level: intermediate 3259b54502bSHong Zhang 3269b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 3279b54502bSHong Zhang @*/ 3282401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 3299b54502bSHong Zhang { 3309b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 3319b54502bSHong Zhang 3329b54502bSHong Zhang PetscFunctionBegin; 3339b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3342401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 3359b54502bSHong Zhang if (f) { 3369b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 3379b54502bSHong Zhang } 3389b54502bSHong Zhang PetscFunctionReturn(0); 3399b54502bSHong Zhang } 3409b54502bSHong Zhang 3419b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 3429b54502bSHong Zhang 3439b54502bSHong Zhang EXTERN_C_BEGIN 3449b54502bSHong Zhang #undef __FUNCT__ 3452401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU" 3462401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 3479b54502bSHong Zhang { 3489b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 3499b54502bSHong Zhang 3509b54502bSHong Zhang PetscFunctionBegin; 3519b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 3529b54502bSHong Zhang PetscFunctionReturn(0); 3539b54502bSHong Zhang } 3549b54502bSHong Zhang EXTERN_C_END 3559b54502bSHong Zhang 3569b54502bSHong Zhang #undef __FUNCT__ 3579b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 3589b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 3599b54502bSHong Zhang { 3609b54502bSHong Zhang PetscErrorCode ierr; 3619b54502bSHong Zhang PetscInt dtmax = 3,itmp; 3629b54502bSHong Zhang PetscTruth flg,set; 3639b54502bSHong Zhang PetscReal dt[3]; 3649b54502bSHong Zhang char tname[256]; 3659b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 3669b54502bSHong Zhang PetscFList ordlist; 3679b54502bSHong Zhang PetscReal tol; 3689b54502bSHong Zhang 3699b54502bSHong Zhang PetscFunctionBegin; 3709b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 3719b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 3722401956bSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 3739b54502bSHong Zhang if (flg) ilu->info.levels = itmp; 374edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 375edbb7feeSSatish Balay if (flg) ilu->inplace = PETSC_TRUE; 3762401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 3779b54502bSHong Zhang ilu->info.diagonal_fill = (double) flg; 378edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 379edbb7feeSSatish Balay if (flg) ilu->reusefill = PETSC_TRUE; 380edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 381edbb7feeSSatish Balay if (flg) ilu->reuseordering = PETSC_TRUE; 3829f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 3839b54502bSHong Zhang if (flg) { 384afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 3859b54502bSHong Zhang } 3869f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 3875e8efad8SHong Zhang 3889f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 3899b54502bSHong Zhang if (flg) { 3909f95998fSHong Zhang ierr = PetscOptionsInt("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr); 3919b54502bSHong Zhang if (flg && !itmp) { 392afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 3939b54502bSHong Zhang } else { 394afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 3959b54502bSHong Zhang } 3969b54502bSHong Zhang } 397ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 3989b54502bSHong Zhang 3999b54502bSHong Zhang dt[0] = ilu->info.dt; 4009b54502bSHong Zhang dt[1] = ilu->info.dtcol; 4019b54502bSHong Zhang dt[2] = ilu->info.dtcount; 4022401956bSBarry Smith ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 4039b54502bSHong Zhang if (flg) { 4042401956bSBarry Smith ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 4059b54502bSHong Zhang } 40655ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 4072401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 4089b54502bSHong Zhang if (flg) { 4099b54502bSHong Zhang tol = PETSC_DECIDE; 4102401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 4112401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 4129b54502bSHong Zhang } 4139b54502bSHong Zhang 4149b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 415*e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 4169b54502bSHong Zhang if (flg) { 417*e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 4189b54502bSHong Zhang } 4199b54502bSHong Zhang flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 4202401956bSBarry Smith ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 4219b54502bSHong Zhang if (set) { 4222401956bSBarry Smith ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 4239b54502bSHong Zhang } 4249b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4259b54502bSHong Zhang PetscFunctionReturn(0); 4269b54502bSHong Zhang } 4279b54502bSHong Zhang 4289b54502bSHong Zhang #undef __FUNCT__ 4299b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 4309b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 4319b54502bSHong Zhang { 4329b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 4339b54502bSHong Zhang PetscErrorCode ierr; 4349b54502bSHong Zhang PetscTruth isstring,iascii; 4359b54502bSHong Zhang 4369b54502bSHong Zhang PetscFunctionBegin; 4379b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 4389b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 4399b54502bSHong Zhang if (iascii) { 4409b54502bSHong Zhang if (ilu->usedt) { 441a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr); 4429b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 443a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr); 4449b54502bSHong Zhang } else if (ilu->info.levels == 1) { 4459b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 4469b54502bSHong Zhang } else { 4479b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 4489b54502bSHong Zhang } 449f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr); 450a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr); 4510a29876aSHong Zhang if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 4529b54502bSHong Zhang if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 4539b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 4549b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 4559b54502bSHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 4569b54502bSHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 4579b54502bSHong Zhang if (ilu->fact) { 458f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr); 4599b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 4609b54502bSHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 461f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 462f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4639b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 4649b54502bSHong Zhang ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 4659b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 4669b54502bSHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 467f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 468f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4699b54502bSHong Zhang } 4709b54502bSHong Zhang } else if (isstring) { 4719b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 4729b54502bSHong Zhang } else { 4739b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 4749b54502bSHong Zhang } 4759b54502bSHong Zhang PetscFunctionReturn(0); 4769b54502bSHong Zhang } 4779b54502bSHong Zhang 4789b54502bSHong Zhang #undef __FUNCT__ 4799b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 4809b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 4819b54502bSHong Zhang { 4829b54502bSHong Zhang PetscErrorCode ierr; 4839b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 484f3a39becSBarry Smith MatInfo info; 4859b54502bSHong Zhang 4869b54502bSHong Zhang PetscFunctionBegin; 4879b54502bSHong Zhang if (ilu->inplace) { 4889b54502bSHong Zhang if (!pc->setupcalled) { 4899b54502bSHong Zhang 4909b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 4919b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 4929b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 493efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 494efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 4959b54502bSHong Zhang } 4969b54502bSHong Zhang 4979b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 4989b54502bSHong Zhang cannot have levels of fill */ 4999b54502bSHong Zhang ilu->info.fill = 1.0; 5009b54502bSHong Zhang ilu->info.diagonal_fill = 0; 5019b54502bSHong Zhang ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 5029b54502bSHong Zhang ilu->fact = pc->pmat; 5039b54502bSHong Zhang } else if (ilu->usedt) { 5049b54502bSHong Zhang if (!pc->setupcalled) { 5059b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 506efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 507efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5087529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 50952e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5109b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 5119b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5129b54502bSHong Zhang if (!ilu->reuseordering) { 5139b54502bSHong Zhang if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 5149b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 5159b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 516efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 517efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5189b54502bSHong Zhang } 5197529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 52052e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5219b54502bSHong Zhang } else if (!ilu->reusefill) { 5229b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5237529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 52452e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5259b54502bSHong Zhang } else { 5269b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 5279b54502bSHong Zhang } 5289b54502bSHong Zhang } else { 5299b54502bSHong Zhang if (!pc->setupcalled) { 5309b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 5319b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 532efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 533efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5349b54502bSHong Zhang /* Remove zeros along diagonal? */ 5359b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5369b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5379b54502bSHong Zhang } 5389b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 539f3a39becSBarry Smith ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 540f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 54152e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5429b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 5439b54502bSHong Zhang if (!ilu->reuseordering) { 5449b54502bSHong Zhang /* compute a new ordering for the ILU */ 5459b54502bSHong Zhang ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 5469b54502bSHong Zhang ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 5479b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 548efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 549efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5509b54502bSHong Zhang /* Remove zeros along diagonal? */ 5519b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5529b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5539b54502bSHong Zhang } 5549b54502bSHong Zhang } 5559b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5569b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 557f3a39becSBarry Smith ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 558f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 55952e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5609b54502bSHong Zhang } 5619b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 5629b54502bSHong Zhang } 5639b54502bSHong Zhang PetscFunctionReturn(0); 5649b54502bSHong Zhang } 5659b54502bSHong Zhang 5669b54502bSHong Zhang #undef __FUNCT__ 5679b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 5689b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 5699b54502bSHong Zhang { 5709b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5719b54502bSHong Zhang PetscErrorCode ierr; 5729b54502bSHong Zhang 5739b54502bSHong Zhang PetscFunctionBegin; 5749b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 5759b54502bSHong Zhang ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 5769b54502bSHong Zhang ierr = PetscFree(ilu);CHKERRQ(ierr); 5779b54502bSHong Zhang PetscFunctionReturn(0); 5789b54502bSHong Zhang } 5799b54502bSHong Zhang 5809b54502bSHong Zhang #undef __FUNCT__ 5819b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 5829b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 5839b54502bSHong Zhang { 5849b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5859b54502bSHong Zhang PetscErrorCode ierr; 5869b54502bSHong Zhang 5879b54502bSHong Zhang PetscFunctionBegin; 5889b54502bSHong Zhang ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 5899b54502bSHong Zhang PetscFunctionReturn(0); 5909b54502bSHong Zhang } 5919b54502bSHong Zhang 5929b54502bSHong Zhang #undef __FUNCT__ 5939b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 5949b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 5959b54502bSHong Zhang { 5969b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5979b54502bSHong Zhang PetscErrorCode ierr; 5989b54502bSHong Zhang 5999b54502bSHong Zhang PetscFunctionBegin; 6009b54502bSHong Zhang ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 6019b54502bSHong Zhang PetscFunctionReturn(0); 6029b54502bSHong Zhang } 6039b54502bSHong Zhang 6049b54502bSHong Zhang #undef __FUNCT__ 6059b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ILU" 6069b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 6079b54502bSHong Zhang { 6089b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6099b54502bSHong Zhang 6109b54502bSHong Zhang PetscFunctionBegin; 6119b54502bSHong Zhang if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 6129b54502bSHong Zhang *mat = ilu->fact; 6139b54502bSHong Zhang PetscFunctionReturn(0); 6149b54502bSHong Zhang } 6159b54502bSHong Zhang 6169b54502bSHong Zhang /*MC 6179b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 6189b54502bSHong Zhang 6199b54502bSHong Zhang Options Database Keys: 6202401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 6212401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 6229b54502bSHong Zhang its factorization (overwrites original matrix) 6232401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 6242401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 6252401956bSBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 62655ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 6272401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 6289b54502bSHong Zhang this decreases the chance of getting a zero pivot 6292401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 6302401956bSBarry Smith . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 6319b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 6329b54502bSHong Zhang stability of the ILU factorization 633f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 634f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 635f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 6369b54502bSHong Zhang 6379b54502bSHong Zhang Level: beginner 6389b54502bSHong Zhang 6399b54502bSHong Zhang Concepts: incomplete factorization 6409b54502bSHong Zhang 6419b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 6429b54502bSHong Zhang 6439b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 6449b54502bSHong Zhang 6459b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 6462401956bSBarry Smith PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 647*e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 6482401956bSBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 649f251bdbdSHong Zhang PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 6509b54502bSHong Zhang 6519b54502bSHong Zhang M*/ 6529b54502bSHong Zhang 6539b54502bSHong Zhang EXTERN_C_BEGIN 6549b54502bSHong Zhang #undef __FUNCT__ 6559b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 656dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 6579b54502bSHong Zhang { 6589b54502bSHong Zhang PetscErrorCode ierr; 6599b54502bSHong Zhang PC_ILU *ilu; 6609b54502bSHong Zhang 6619b54502bSHong Zhang PetscFunctionBegin; 6629b54502bSHong Zhang ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr); 66352e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_ILU));CHKERRQ(ierr); 6649b54502bSHong Zhang 6659b54502bSHong Zhang ilu->fact = 0; 6669b54502bSHong Zhang ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 6679b54502bSHong Zhang ilu->info.levels = 0; 6689b54502bSHong Zhang ilu->info.fill = 1.0; 6699b54502bSHong Zhang ilu->col = 0; 6709b54502bSHong Zhang ilu->row = 0; 6719b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 6729b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 6739b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 6749b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 6759b54502bSHong Zhang ilu->info.dt = PETSC_DEFAULT; 6769b54502bSHong Zhang ilu->info.dtcount = PETSC_DEFAULT; 6779b54502bSHong Zhang ilu->info.dtcol = PETSC_DEFAULT; 6780a29876aSHong Zhang ilu->info.shiftnz = 0.0; 679fbf22428SSatish Balay ilu->info.shiftpd = 0.0; /* false */ 6809b54502bSHong Zhang ilu->info.shift_fraction = 0.0; 6819b54502bSHong Zhang ilu->info.zeropivot = 1.e-12; 6829b54502bSHong Zhang ilu->info.pivotinblocks = 1.0; 6839b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 6849b54502bSHong Zhang ilu->info.diagonal_fill = 0; 6859b54502bSHong Zhang pc->data = (void*)ilu; 6869b54502bSHong Zhang 6879b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 6889b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 6899b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 6909b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 6919b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 6929b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 6939b54502bSHong Zhang pc->ops->view = PCView_ILU; 6949b54502bSHong Zhang pc->ops->applyrichardson = 0; 6959b54502bSHong Zhang 696afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 697afaefe49SHong Zhang PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 698afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 699afaefe49SHong Zhang PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 700afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 701afaefe49SHong Zhang PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 702afaefe49SHong Zhang 7032401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 7042401956bSBarry Smith PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 70555ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 70655ba2a51SBarry Smith PCFactorSetFill_ILU);CHKERRQ(ierr); 707*e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ILU", 708*e5a9bf91SBarry Smith PCFactorSetMatOrderingType_ILU);CHKERRQ(ierr); 7092401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 7102401956bSBarry Smith PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 7112401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 7122401956bSBarry Smith PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 7132401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU", 7142401956bSBarry Smith PCFactorSetLevels_ILU);CHKERRQ(ierr); 7152401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 7162401956bSBarry Smith PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 7172401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU", 7182401956bSBarry Smith PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 7192401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU", 7202401956bSBarry Smith PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr); 7212401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 7222401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 7239b54502bSHong Zhang PetscFunctionReturn(0); 7249b54502bSHong Zhang } 7259b54502bSHong Zhang EXTERN_C_END 726