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__ 147e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ILU" 148e5a9bf91SBarry 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 2603e32baedSBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 2613e32baedSBarry Smith matrix. We don't know how to compute reasonable values. 2623e32baedSBarry Smith 2639b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 2649b54502bSHong Zhang @*/ 2652401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 2669b54502bSHong Zhang { 2679b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 2689b54502bSHong Zhang 2699b54502bSHong Zhang PetscFunctionBegin; 2709b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2712401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 2729b54502bSHong Zhang if (f) { 2739b54502bSHong Zhang ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 2749b54502bSHong Zhang } 2759b54502bSHong Zhang PetscFunctionReturn(0); 2769b54502bSHong Zhang } 2779b54502bSHong Zhang 2789b54502bSHong Zhang #undef __FUNCT__ 2792401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels" 2809b54502bSHong Zhang /*@ 2812401956bSBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 2829b54502bSHong Zhang 2839b54502bSHong Zhang Collective on PC 2849b54502bSHong Zhang 2859b54502bSHong Zhang Input Parameters: 2869b54502bSHong Zhang + pc - the preconditioner context 2879b54502bSHong Zhang - levels - number of levels of fill 2889b54502bSHong Zhang 2899b54502bSHong Zhang Options Database Key: 2902401956bSBarry Smith . -pc_factor_levels <levels> - Sets fill level 2919b54502bSHong Zhang 2929b54502bSHong Zhang Level: intermediate 2939b54502bSHong Zhang 2949b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 2959b54502bSHong Zhang @*/ 2962401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 2979b54502bSHong Zhang { 2989b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 2999b54502bSHong Zhang 3009b54502bSHong Zhang PetscFunctionBegin; 3019b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3029b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 3032401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 3049b54502bSHong Zhang if (f) { 3059b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 3069b54502bSHong Zhang } 3079b54502bSHong Zhang PetscFunctionReturn(0); 3089b54502bSHong Zhang } 3099b54502bSHong Zhang 3109b54502bSHong Zhang #undef __FUNCT__ 3112401956bSBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 3129b54502bSHong Zhang /*@ 3132401956bSBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 3149b54502bSHong Zhang treated as level 0 fill even if there is no non-zero location. 3159b54502bSHong Zhang 3169b54502bSHong Zhang Collective on PC 3179b54502bSHong Zhang 3189b54502bSHong Zhang Input Parameters: 3199b54502bSHong Zhang + pc - the preconditioner context 3209b54502bSHong Zhang 3219b54502bSHong Zhang Options Database Key: 3222401956bSBarry Smith . -pc_factor_diagonal_fill 3239b54502bSHong Zhang 3249b54502bSHong Zhang Notes: 3259b54502bSHong Zhang Does not apply with 0 fill. 3269b54502bSHong Zhang 3279b54502bSHong Zhang Level: intermediate 3289b54502bSHong Zhang 3299b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 3309b54502bSHong Zhang @*/ 3312401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 3329b54502bSHong Zhang { 3339b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 3349b54502bSHong Zhang 3359b54502bSHong Zhang PetscFunctionBegin; 3369b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3372401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 3389b54502bSHong Zhang if (f) { 3399b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 3409b54502bSHong Zhang } 3419b54502bSHong Zhang PetscFunctionReturn(0); 3429b54502bSHong Zhang } 3439b54502bSHong Zhang 3449b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 3459b54502bSHong Zhang 3469b54502bSHong Zhang EXTERN_C_BEGIN 3479b54502bSHong Zhang #undef __FUNCT__ 3482401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU" 3492401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 3509b54502bSHong Zhang { 3519b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 3529b54502bSHong Zhang 3539b54502bSHong Zhang PetscFunctionBegin; 3549b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 3559b54502bSHong Zhang PetscFunctionReturn(0); 3569b54502bSHong Zhang } 3579b54502bSHong Zhang EXTERN_C_END 3589b54502bSHong Zhang 35962bba022SBarry Smith EXTERN_C_BEGIN 36062bba022SBarry Smith #undef __FUNCT__ 36162bba022SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks_ILU" 36262bba022SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_ILU(PC pc,PetscReal shift) 36362bba022SBarry Smith { 36462bba022SBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 36562bba022SBarry Smith 36662bba022SBarry Smith PetscFunctionBegin; 36762bba022SBarry Smith if (shift == PETSC_DEFAULT) { 36862bba022SBarry Smith dir->info.shiftinblocks = 1.e-12; 36962bba022SBarry Smith } else { 37062bba022SBarry Smith dir->info.shiftinblocks = shift; 37162bba022SBarry Smith } 37262bba022SBarry Smith PetscFunctionReturn(0); 37362bba022SBarry Smith } 37462bba022SBarry Smith EXTERN_C_END 37562bba022SBarry Smith 37662bba022SBarry Smith 3779b54502bSHong Zhang #undef __FUNCT__ 3789b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 3799b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 3809b54502bSHong Zhang { 3819b54502bSHong Zhang PetscErrorCode ierr; 3829b54502bSHong Zhang PetscInt dtmax = 3,itmp; 3839b54502bSHong Zhang PetscTruth flg,set; 3849b54502bSHong Zhang PetscReal dt[3]; 3859b54502bSHong Zhang char tname[256]; 3869b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 3879b54502bSHong Zhang PetscFList ordlist; 3889b54502bSHong Zhang PetscReal tol; 3899b54502bSHong Zhang 3909b54502bSHong Zhang PetscFunctionBegin; 3919b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 3929b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 3932401956bSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 3949b54502bSHong Zhang if (flg) ilu->info.levels = itmp; 395edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 396edbb7feeSSatish Balay if (flg) ilu->inplace = PETSC_TRUE; 3972401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 3989b54502bSHong Zhang ilu->info.diagonal_fill = (double) flg; 399edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 400edbb7feeSSatish Balay if (flg) ilu->reusefill = PETSC_TRUE; 401edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 402edbb7feeSSatish Balay if (flg) ilu->reuseordering = PETSC_TRUE; 4039f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 4049b54502bSHong Zhang if (flg) { 405afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 4069b54502bSHong Zhang } 4079f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 40862bba022SBarry Smith flg = (ilu->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE; 40962bba022SBarry Smith ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 41062bba022SBarry Smith ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr); 411ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 4129b54502bSHong Zhang 4139b54502bSHong Zhang dt[0] = ilu->info.dt; 4149b54502bSHong Zhang dt[1] = ilu->info.dtcol; 4159b54502bSHong Zhang dt[2] = ilu->info.dtcount; 4162401956bSBarry Smith ierr = PetscOptionsRealArray("-pc_factor_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCFactorSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 4179b54502bSHong Zhang if (flg) { 4182401956bSBarry Smith ierr = PCFactorSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 4199b54502bSHong Zhang } 42055ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 4212401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 4229b54502bSHong Zhang if (flg) { 4239b54502bSHong Zhang tol = PETSC_DECIDE; 4242401956bSBarry Smith ierr = PetscOptionsReal("-pc_factor_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCFactorReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 4252401956bSBarry Smith ierr = PCFactorReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 4269b54502bSHong Zhang } 4279b54502bSHong Zhang 4289b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 429e5a9bf91SBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 4309b54502bSHong Zhang if (flg) { 431e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 4329b54502bSHong Zhang } 4339b54502bSHong Zhang flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 4342401956bSBarry Smith ierr = PetscOptionsTruth("-pc_factor_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCFactorSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 4359b54502bSHong Zhang if (set) { 4362401956bSBarry Smith ierr = PCFactorSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 4379b54502bSHong Zhang } 43862bba022SBarry Smith ierr = PetscOptionsName("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",&flg);CHKERRQ(ierr); 43962bba022SBarry Smith if (flg) { 44062bba022SBarry Smith ierr = PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 44162bba022SBarry Smith } 44262bba022SBarry Smith ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",ilu->info.shiftinblocks,&ilu->info.shiftinblocks,0);CHKERRQ(ierr); 44362bba022SBarry Smith 4449b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4459b54502bSHong Zhang PetscFunctionReturn(0); 4469b54502bSHong Zhang } 4479b54502bSHong Zhang 4489b54502bSHong Zhang #undef __FUNCT__ 4499b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 4509b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 4519b54502bSHong Zhang { 4529b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 4539b54502bSHong Zhang PetscErrorCode ierr; 4549b54502bSHong Zhang PetscTruth isstring,iascii; 4559b54502bSHong Zhang 4569b54502bSHong Zhang PetscFunctionBegin; 4579b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 4589b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 4599b54502bSHong Zhang if (iascii) { 4609b54502bSHong Zhang if (ilu->usedt) { 461a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr); 4629b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 463a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr); 4649b54502bSHong Zhang } else if (ilu->info.levels == 1) { 4659b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 4669b54502bSHong Zhang } else { 4679b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 4689b54502bSHong Zhang } 469f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr); 470a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr); 4710a29876aSHong Zhang if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 47262bba022SBarry Smith if (ilu->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);} 47362bba022SBarry Smith if (ilu->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);} 4749b54502bSHong Zhang if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 4759b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 4769b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 4779b54502bSHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 4789b54502bSHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 4799b54502bSHong Zhang if (ilu->fact) { 480f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr); 4819b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 4829b54502bSHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 483f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 484f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4859b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 4869b54502bSHong Zhang ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 4879b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 4889b54502bSHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 489f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 490f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4919b54502bSHong Zhang } 4929b54502bSHong Zhang } else if (isstring) { 4939b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 4949b54502bSHong Zhang } else { 4959b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 4969b54502bSHong Zhang } 4979b54502bSHong Zhang PetscFunctionReturn(0); 4989b54502bSHong Zhang } 4999b54502bSHong Zhang 5009b54502bSHong Zhang #undef __FUNCT__ 5019b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 5029b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 5039b54502bSHong Zhang { 5049b54502bSHong Zhang PetscErrorCode ierr; 5059b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 506f3a39becSBarry Smith MatInfo info; 5079b54502bSHong Zhang 5089b54502bSHong Zhang PetscFunctionBegin; 5099b54502bSHong Zhang if (ilu->inplace) { 5109b54502bSHong Zhang if (!pc->setupcalled) { 5119b54502bSHong Zhang 5129b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 5139b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 5149b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 515efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 516efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5179b54502bSHong Zhang } 5189b54502bSHong Zhang 5199b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 5209b54502bSHong Zhang cannot have levels of fill */ 5219b54502bSHong Zhang ilu->info.fill = 1.0; 5229b54502bSHong Zhang ilu->info.diagonal_fill = 0; 5239b54502bSHong Zhang ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 5249b54502bSHong Zhang ilu->fact = pc->pmat; 5259b54502bSHong Zhang } else if (ilu->usedt) { 5269b54502bSHong Zhang if (!pc->setupcalled) { 5279b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 528efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 529efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5307529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 53152e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5329b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 5339b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5349b54502bSHong Zhang if (!ilu->reuseordering) { 5359b54502bSHong Zhang if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 5369b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 5379b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 538efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 539efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5409b54502bSHong Zhang } 5417529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 54252e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5439b54502bSHong Zhang } else if (!ilu->reusefill) { 5449b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5457529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 54652e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5479b54502bSHong Zhang } else { 5489b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 5499b54502bSHong Zhang } 5509b54502bSHong Zhang } else { 5519b54502bSHong Zhang if (!pc->setupcalled) { 5529b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 5539b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 554efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 555efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5569b54502bSHong Zhang /* Remove zeros along diagonal? */ 5579b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5589b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5599b54502bSHong Zhang } 5609b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 561f3a39becSBarry Smith ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 562f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 56352e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5649b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 5659b54502bSHong Zhang if (!ilu->reuseordering) { 5669b54502bSHong Zhang /* compute a new ordering for the ILU */ 5679b54502bSHong Zhang ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 5689b54502bSHong Zhang ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 5699b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 570efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 571efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5729b54502bSHong Zhang /* Remove zeros along diagonal? */ 5739b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5749b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5759b54502bSHong Zhang } 5769b54502bSHong Zhang } 5779b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 5789b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 579f3a39becSBarry Smith ierr = MatGetInfo(ilu->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 580f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 58152e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 5829b54502bSHong Zhang } 5839b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 5849b54502bSHong Zhang } 5859b54502bSHong Zhang PetscFunctionReturn(0); 5869b54502bSHong Zhang } 5879b54502bSHong Zhang 5889b54502bSHong Zhang #undef __FUNCT__ 5899b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 5909b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 5919b54502bSHong Zhang { 5929b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5939b54502bSHong Zhang PetscErrorCode ierr; 5949b54502bSHong Zhang 5959b54502bSHong Zhang PetscFunctionBegin; 5969b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 5979b54502bSHong Zhang ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 5989b54502bSHong Zhang ierr = PetscFree(ilu);CHKERRQ(ierr); 5999b54502bSHong Zhang PetscFunctionReturn(0); 6009b54502bSHong Zhang } 6019b54502bSHong Zhang 6029b54502bSHong Zhang #undef __FUNCT__ 6039b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 6049b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 6059b54502bSHong Zhang { 6069b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6079b54502bSHong Zhang PetscErrorCode ierr; 6089b54502bSHong Zhang 6099b54502bSHong Zhang PetscFunctionBegin; 6109b54502bSHong Zhang ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 6119b54502bSHong Zhang PetscFunctionReturn(0); 6129b54502bSHong Zhang } 6139b54502bSHong Zhang 6149b54502bSHong Zhang #undef __FUNCT__ 6159b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 6169b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 6179b54502bSHong Zhang { 6189b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6199b54502bSHong Zhang PetscErrorCode ierr; 6209b54502bSHong Zhang 6219b54502bSHong Zhang PetscFunctionBegin; 6229b54502bSHong Zhang ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 6239b54502bSHong Zhang PetscFunctionReturn(0); 6249b54502bSHong Zhang } 6259b54502bSHong Zhang 6269b54502bSHong Zhang #undef __FUNCT__ 627a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_ILU" 628a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_ILU(PC pc,Mat *mat) 6299b54502bSHong Zhang { 6309b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6319b54502bSHong Zhang 6329b54502bSHong Zhang PetscFunctionBegin; 6339b54502bSHong Zhang if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 6349b54502bSHong Zhang *mat = ilu->fact; 6359b54502bSHong Zhang PetscFunctionReturn(0); 6369b54502bSHong Zhang } 6379b54502bSHong Zhang 63862bba022SBarry Smith #undef __FUNCT__ 63962bba022SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks" 64062bba022SBarry Smith /*@ 64162bba022SBarry Smith PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot. 64262bba022SBarry Smith 64362bba022SBarry Smith Collective on PC 64462bba022SBarry Smith 64562bba022SBarry Smith Input Parameters: 64662bba022SBarry Smith + pc - the preconditioner context 64762bba022SBarry Smith - shift - amount to shift or PETSC_DECIDE 64862bba022SBarry Smith 64962bba022SBarry Smith Options Database Key: 65062bba022SBarry Smith . -pc_factor_shift_in_blocks <shift> 65162bba022SBarry Smith 65262bba022SBarry Smith Level: intermediate 65362bba022SBarry Smith 65462bba022SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks() 65562bba022SBarry Smith @*/ 65662bba022SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift) 65762bba022SBarry Smith { 65862bba022SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 65962bba022SBarry Smith 66062bba022SBarry Smith PetscFunctionBegin; 66162bba022SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 66262bba022SBarry Smith if (f) { 66362bba022SBarry Smith ierr = (*f)(pc,shift);CHKERRQ(ierr); 66462bba022SBarry Smith } 66562bba022SBarry Smith PetscFunctionReturn(0); 66662bba022SBarry Smith } 66762bba022SBarry Smith 6689b54502bSHong Zhang /*MC 6699b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 6709b54502bSHong Zhang 6719b54502bSHong Zhang Options Database Keys: 6722401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 6732401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 6749b54502bSHong Zhang its factorization (overwrites original matrix) 6752401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 6762401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 6772401956bSBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 67855ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 6792401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 6809b54502bSHong Zhang this decreases the chance of getting a zero pivot 6812401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 6822401956bSBarry Smith . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 6839b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 6849b54502bSHong Zhang stability of the ILU factorization 68562bba022SBarry Smith . -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization 686f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 687f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 688f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 6899b54502bSHong Zhang 6909b54502bSHong Zhang Level: beginner 6919b54502bSHong Zhang 6929b54502bSHong Zhang Concepts: incomplete factorization 6939b54502bSHong Zhang 694*c582cd25SBarry Smith Notes: Only implemented for some matrix formats. (for parallel use you 695*c582cd25SBarry Smith must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended 696*c582cd25SBarry Smith unless you really want a parallel ILU). 6979b54502bSHong Zhang 6989b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 6999b54502bSHong Zhang 700*c582cd25SBarry Smith References: 701*c582cd25SBarry Smith T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 702*c582cd25SBarry Smith self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 703*c582cd25SBarry Smith 704*c582cd25SBarry Smith T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 705*c582cd25SBarry Smith fusion problems. Quart. Appl. Math., 19:221--229, 1961. 706*c582cd25SBarry Smith 707*c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 708*c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 709*c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 710*c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 711*c582cd25SBarry Smith 712*c582cd25SBarry Smith 7139b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 7142401956bSBarry Smith PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 715e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 7162401956bSBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 717f251bdbdSHong Zhang PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 7189b54502bSHong Zhang 7199b54502bSHong Zhang M*/ 7209b54502bSHong Zhang 7219b54502bSHong Zhang EXTERN_C_BEGIN 7229b54502bSHong Zhang #undef __FUNCT__ 7239b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 724dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 7259b54502bSHong Zhang { 7269b54502bSHong Zhang PetscErrorCode ierr; 7279b54502bSHong Zhang PC_ILU *ilu; 7289b54502bSHong Zhang 7299b54502bSHong Zhang PetscFunctionBegin; 73038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 7319b54502bSHong Zhang 7329b54502bSHong Zhang ilu->fact = 0; 7339b54502bSHong Zhang ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 7349b54502bSHong Zhang ilu->info.levels = 0; 7359b54502bSHong Zhang ilu->info.fill = 1.0; 7369b54502bSHong Zhang ilu->col = 0; 7379b54502bSHong Zhang ilu->row = 0; 7389b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 7399b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 7409b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 7419b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 7429b54502bSHong Zhang ilu->info.dt = PETSC_DEFAULT; 7439b54502bSHong Zhang ilu->info.dtcount = PETSC_DEFAULT; 7449b54502bSHong Zhang ilu->info.dtcol = PETSC_DEFAULT; 74562bba022SBarry Smith ilu->info.shiftnz = 1.e-12; 746fbf22428SSatish Balay ilu->info.shiftpd = 0.0; /* false */ 7479b54502bSHong Zhang ilu->info.shift_fraction = 0.0; 7489b54502bSHong Zhang ilu->info.zeropivot = 1.e-12; 7499b54502bSHong Zhang ilu->info.pivotinblocks = 1.0; 75062bba022SBarry Smith ilu->info.shiftinblocks = 1.e-12; 7519b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 7529b54502bSHong Zhang ilu->info.diagonal_fill = 0; 7539b54502bSHong Zhang pc->data = (void*)ilu; 7549b54502bSHong Zhang 7559b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 7569b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 7579b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 7589b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 7599b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 760a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_ILU; 7619b54502bSHong Zhang pc->ops->view = PCView_ILU; 7629b54502bSHong Zhang pc->ops->applyrichardson = 0; 7639b54502bSHong Zhang 764afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 765afaefe49SHong Zhang PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 766afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 767afaefe49SHong Zhang PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 768afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 769afaefe49SHong Zhang PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 770afaefe49SHong Zhang 7712401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 7722401956bSBarry Smith PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 77355ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 77455ba2a51SBarry Smith PCFactorSetFill_ILU);CHKERRQ(ierr); 775e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ILU", 776e5a9bf91SBarry Smith PCFactorSetMatOrderingType_ILU);CHKERRQ(ierr); 7772401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 7782401956bSBarry Smith PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 7792401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 7802401956bSBarry Smith PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 7812401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU", 7822401956bSBarry Smith PCFactorSetLevels_ILU);CHKERRQ(ierr); 7832401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 7842401956bSBarry Smith PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 7852401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU", 7862401956bSBarry Smith PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 7872401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU", 7882401956bSBarry Smith PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr); 78962bba022SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_ILU", 79062bba022SBarry Smith PCFactorSetShiftInBlocks_ILU);CHKERRQ(ierr); 7912401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 7922401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 7939b54502bSHong Zhang PetscFunctionReturn(0); 7949b54502bSHong Zhang } 7959b54502bSHong Zhang EXTERN_C_END 796