1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 39b54502bSHong Zhang /* 49b54502bSHong Zhang Defines a ILU factorization preconditioner for any Mat implementation 59b54502bSHong Zhang */ 6*075768bcSBarry Smith #include "src/ksp/pc/impls/factor/ilu/ilu.h" /*I "petscpc.h" I*/ 79b54502bSHong Zhang 89b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 9afaefe49SHong Zhang EXTERN_C_BEGIN 10afaefe49SHong Zhang #undef __FUNCT__ 112401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseFill_ILU" 122401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseFill_ILU(PC pc,PetscTruth flag) 132401956bSBarry Smith { 14*075768bcSBarry Smith PC_ILU *lu = (PC_ILU*)pc->data; 152401956bSBarry Smith 162401956bSBarry Smith PetscFunctionBegin; 172401956bSBarry Smith lu->reusefill = flag; 182401956bSBarry Smith PetscFunctionReturn(0); 192401956bSBarry Smith } 202401956bSBarry Smith EXTERN_C_END 212401956bSBarry Smith 222401956bSBarry Smith EXTERN_C_BEGIN 232401956bSBarry Smith #undef __FUNCT__ 24afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ILU" 25dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ILU(PC pc,PetscReal z) 26afaefe49SHong Zhang { 27*075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 28afaefe49SHong Zhang 29afaefe49SHong Zhang PetscFunctionBegin; 30*075768bcSBarry Smith ((PC_Factor*)ilu)->info.zeropivot = z; 31afaefe49SHong Zhang PetscFunctionReturn(0); 32afaefe49SHong Zhang } 33afaefe49SHong Zhang EXTERN_C_END 34afaefe49SHong Zhang 35afaefe49SHong Zhang EXTERN_C_BEGIN 36afaefe49SHong Zhang #undef __FUNCT__ 37afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ILU" 38dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ILU(PC pc,PetscReal shift) 39afaefe49SHong Zhang { 40*075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 41afaefe49SHong Zhang 42afaefe49SHong Zhang PetscFunctionBegin; 43afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 44*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = 1.e-12; 45afaefe49SHong Zhang } else { 46*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftnz = shift; 47afaefe49SHong Zhang } 48afaefe49SHong Zhang PetscFunctionReturn(0); 49afaefe49SHong Zhang } 50afaefe49SHong Zhang EXTERN_C_END 51afaefe49SHong Zhang 52afaefe49SHong Zhang EXTERN_C_BEGIN 53afaefe49SHong Zhang #undef __FUNCT__ 54afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ILU" 55dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ILU(PC pc,PetscTruth shift) 56afaefe49SHong Zhang { 57*075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 58afaefe49SHong Zhang 59afaefe49SHong Zhang PetscFunctionBegin; 60fbf22428SSatish Balay if (shift) { 61*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 1.0; 62fbf22428SSatish Balay } else { 63*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftpd = 0.0; 64fbf22428SSatish Balay } 65afaefe49SHong Zhang PetscFunctionReturn(0); 66afaefe49SHong Zhang } 67afaefe49SHong Zhang EXTERN_C_END 689b54502bSHong Zhang 699b54502bSHong Zhang EXTERN_C_BEGIN 709b54502bSHong Zhang #undef __FUNCT__ 712401956bSBarry Smith #define __FUNCT__ "PCFactorReorderForNonzeroDiagonal_ILU" 722401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 739b54502bSHong Zhang { 749b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 759b54502bSHong Zhang 769b54502bSHong Zhang PetscFunctionBegin; 779b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 789b54502bSHong Zhang if (z == PETSC_DECIDE) { 799b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = 1.e-10; 809b54502bSHong Zhang } else { 819b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = z; 829b54502bSHong Zhang } 839b54502bSHong Zhang PetscFunctionReturn(0); 849b54502bSHong Zhang } 859b54502bSHong Zhang EXTERN_C_END 869b54502bSHong Zhang 879b54502bSHong Zhang #undef __FUNCT__ 889b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal" 899b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc) 909b54502bSHong Zhang { 919b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 929b54502bSHong Zhang PetscErrorCode ierr; 939b54502bSHong Zhang 949b54502bSHong Zhang PetscFunctionBegin; 95*075768bcSBarry Smith if (!ilu->inplace && ((PC_Factor*)ilu)->fact) {ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr);} 969b54502bSHong Zhang if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 979b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 989b54502bSHong Zhang PetscFunctionReturn(0); 999b54502bSHong Zhang } 1009b54502bSHong Zhang 1019b54502bSHong Zhang EXTERN_C_BEGIN 1029b54502bSHong Zhang #undef __FUNCT__ 1032401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance_ILU" 1042401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 1059b54502bSHong Zhang { 106*075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 1079b54502bSHong Zhang PetscErrorCode ierr; 1089b54502bSHong Zhang 1099b54502bSHong Zhang PetscFunctionBegin; 110*075768bcSBarry Smith if (pc->setupcalled && (!ilu->usedt || ((PC_Factor*)ilu)->info.dt != dt || ((PC_Factor*)ilu)->info.dtcol != dtcol || ((PC_Factor*)ilu)->info.dtcount != dtcount)) { 1119b54502bSHong Zhang pc->setupcalled = 0; 1129b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1139b54502bSHong Zhang } 1149b54502bSHong Zhang ilu->usedt = PETSC_TRUE; 115*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = dt; 116*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = dtcol; 117*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = dtcount; 118*075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = PETSC_DEFAULT; 1199b54502bSHong Zhang PetscFunctionReturn(0); 1209b54502bSHong Zhang } 1219b54502bSHong Zhang EXTERN_C_END 1229b54502bSHong Zhang 1239b54502bSHong Zhang EXTERN_C_BEGIN 1249b54502bSHong Zhang #undef __FUNCT__ 12555ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ILU" 12655ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ILU(PC pc,PetscReal fill) 1279b54502bSHong Zhang { 128*075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 1299b54502bSHong Zhang 1309b54502bSHong Zhang PetscFunctionBegin; 131*075768bcSBarry Smith ((PC_Factor*)dir)->info.fill = fill; 1329b54502bSHong Zhang PetscFunctionReturn(0); 1339b54502bSHong Zhang } 1349b54502bSHong Zhang EXTERN_C_END 1359b54502bSHong Zhang 1369b54502bSHong Zhang EXTERN_C_BEGIN 1379b54502bSHong Zhang #undef __FUNCT__ 138e5a9bf91SBarry Smith #define __FUNCT__ "PCFactorSetMatOrderingType_ILU" 139e36ff184SSatish Balay PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetMatOrderingType_ILU(PC pc,const MatOrderingType ordering) 1409b54502bSHong Zhang { 1419b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 1429b54502bSHong Zhang PetscErrorCode ierr; 1439b54502bSHong Zhang PetscTruth flg; 1449b54502bSHong Zhang 1459b54502bSHong Zhang PetscFunctionBegin; 1469b54502bSHong Zhang if (!pc->setupcalled) { 147*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 148*075768bcSBarry Smith ierr = PetscStrallocpy(ordering,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 1499b54502bSHong Zhang } else { 150*075768bcSBarry Smith ierr = PetscStrcmp(((PC_Factor*)dir)->ordering,ordering,&flg);CHKERRQ(ierr); 1519b54502bSHong Zhang if (!flg) { 1529b54502bSHong Zhang pc->setupcalled = 0; 153*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 154*075768bcSBarry Smith ierr = PetscStrallocpy(ordering,&((PC_Factor*)dir)->ordering);CHKERRQ(ierr); 1559b54502bSHong Zhang /* free the data structures, then create them again */ 1569b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1579b54502bSHong Zhang } 1589b54502bSHong Zhang } 1599b54502bSHong Zhang PetscFunctionReturn(0); 1609b54502bSHong Zhang } 1619b54502bSHong Zhang EXTERN_C_END 1629b54502bSHong Zhang 1639b54502bSHong Zhang EXTERN_C_BEGIN 1649b54502bSHong Zhang #undef __FUNCT__ 1652401956bSBarry Smith #define __FUNCT__ "PCFactorSetReuseOrdering_ILU" 1662401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetReuseOrdering_ILU(PC pc,PetscTruth flag) 1679b54502bSHong Zhang { 168*075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 1699b54502bSHong Zhang 1709b54502bSHong Zhang PetscFunctionBegin; 1719b54502bSHong Zhang ilu->reuseordering = flag; 1729b54502bSHong Zhang PetscFunctionReturn(0); 1739b54502bSHong Zhang } 1749b54502bSHong Zhang EXTERN_C_END 1759b54502bSHong Zhang 1769b54502bSHong Zhang EXTERN_C_BEGIN 1779b54502bSHong Zhang #undef __FUNCT__ 1782401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels_ILU" 1792401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels_ILU(PC pc,PetscInt levels) 1809b54502bSHong Zhang { 181*075768bcSBarry Smith PC_ILU *ilu = (PC_ILU*)pc->data; 1829b54502bSHong Zhang PetscErrorCode ierr; 1839b54502bSHong Zhang 1849b54502bSHong Zhang PetscFunctionBegin; 1859b54502bSHong Zhang if (!pc->setupcalled) { 186*075768bcSBarry Smith ((PC_Factor*)ilu)->info.levels = levels; 187*075768bcSBarry Smith } else if (ilu->usedt || ((PC_Factor*)ilu)->info.levels != levels) { 188*075768bcSBarry Smith ((PC_Factor*)ilu)->info.levels = levels; 1899b54502bSHong Zhang pc->setupcalled = 0; 1909b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 1919b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1929b54502bSHong Zhang } 1939b54502bSHong Zhang PetscFunctionReturn(0); 1949b54502bSHong Zhang } 1959b54502bSHong Zhang EXTERN_C_END 1969b54502bSHong Zhang 1979b54502bSHong Zhang EXTERN_C_BEGIN 1989b54502bSHong Zhang #undef __FUNCT__ 1992401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseInPlace_ILU" 2002401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseInPlace_ILU(PC pc) 2019b54502bSHong Zhang { 202*075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 2039b54502bSHong Zhang 2049b54502bSHong Zhang PetscFunctionBegin; 2059b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2069b54502bSHong Zhang PetscFunctionReturn(0); 2079b54502bSHong Zhang } 2089b54502bSHong Zhang EXTERN_C_END 2099b54502bSHong Zhang 2109b54502bSHong Zhang EXTERN_C_BEGIN 2119b54502bSHong Zhang #undef __FUNCT__ 21212a54f56SSatish Balay #define __FUNCT__ "PCFactorSetAllowDiagonalFill_ILU" 2132401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill_ILU(PC pc) 2149b54502bSHong Zhang { 215*075768bcSBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 2169b54502bSHong Zhang 2179b54502bSHong Zhang PetscFunctionBegin; 218*075768bcSBarry Smith ((PC_Factor*)dir)->info.diagonal_fill = 1; 2199b54502bSHong Zhang PetscFunctionReturn(0); 2209b54502bSHong Zhang } 2219b54502bSHong Zhang EXTERN_C_END 2229b54502bSHong Zhang 2239b54502bSHong Zhang #undef __FUNCT__ 2242401956bSBarry Smith #define __FUNCT__ "PCFactorSetUseDropTolerance" 2259b54502bSHong Zhang /*@ 2262401956bSBarry Smith PCFactorSetUseDropTolerance - The preconditioner will use an ILU 2279b54502bSHong Zhang based on a drop tolerance. 2289b54502bSHong Zhang 2299b54502bSHong Zhang Collective on PC 2309b54502bSHong Zhang 2319b54502bSHong Zhang Input Parameters: 2329b54502bSHong Zhang + pc - the preconditioner context 2339b54502bSHong Zhang . dt - the drop tolerance, try from 1.e-10 to .1 2349b54502bSHong Zhang . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 2359b54502bSHong Zhang - maxrowcount - the max number of nonzeros allowed in a row, best value 2369b54502bSHong Zhang depends on the number of nonzeros in row of original matrix 2379b54502bSHong Zhang 2389b54502bSHong Zhang Options Database Key: 2392401956bSBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 2409b54502bSHong Zhang 2419b54502bSHong Zhang Level: intermediate 2429b54502bSHong Zhang 2439b54502bSHong Zhang Notes: 2449b54502bSHong Zhang This uses the iludt() code of Saad's SPARSKIT package 2459b54502bSHong Zhang 2463e32baedSBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 2473e32baedSBarry Smith matrix. We don't know how to compute reasonable values. 2483e32baedSBarry Smith 2499b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 2509b54502bSHong Zhang @*/ 2512401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 2529b54502bSHong Zhang { 2539b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 2549b54502bSHong Zhang 2559b54502bSHong Zhang PetscFunctionBegin; 2569b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2572401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 2589b54502bSHong Zhang if (f) { 2599b54502bSHong Zhang ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 2609b54502bSHong Zhang } 2619b54502bSHong Zhang PetscFunctionReturn(0); 2629b54502bSHong Zhang } 2639b54502bSHong Zhang 2649b54502bSHong Zhang #undef __FUNCT__ 2652401956bSBarry Smith #define __FUNCT__ "PCFactorSetLevels" 2669b54502bSHong Zhang /*@ 2672401956bSBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 2689b54502bSHong Zhang 2699b54502bSHong Zhang Collective on PC 2709b54502bSHong Zhang 2719b54502bSHong Zhang Input Parameters: 2729b54502bSHong Zhang + pc - the preconditioner context 2739b54502bSHong Zhang - levels - number of levels of fill 2749b54502bSHong Zhang 2759b54502bSHong Zhang Options Database Key: 2762401956bSBarry Smith . -pc_factor_levels <levels> - Sets fill level 2779b54502bSHong Zhang 2789b54502bSHong Zhang Level: intermediate 2799b54502bSHong Zhang 2809b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 2819b54502bSHong Zhang @*/ 2822401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetLevels(PC pc,PetscInt levels) 2839b54502bSHong Zhang { 2849b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 2859b54502bSHong Zhang 2869b54502bSHong Zhang PetscFunctionBegin; 2879b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2889b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 2892401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 2909b54502bSHong Zhang if (f) { 2919b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 2929b54502bSHong Zhang } 2939b54502bSHong Zhang PetscFunctionReturn(0); 2949b54502bSHong Zhang } 2959b54502bSHong Zhang 2969b54502bSHong Zhang #undef __FUNCT__ 2972401956bSBarry Smith #define __FUNCT__ "PCFactorSetAllowDiagonalFill" 2989b54502bSHong Zhang /*@ 2992401956bSBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 3009b54502bSHong Zhang treated as level 0 fill even if there is no non-zero location. 3019b54502bSHong Zhang 3029b54502bSHong Zhang Collective on PC 3039b54502bSHong Zhang 3049b54502bSHong Zhang Input Parameters: 3059b54502bSHong Zhang + pc - the preconditioner context 3069b54502bSHong Zhang 3079b54502bSHong Zhang Options Database Key: 3082401956bSBarry Smith . -pc_factor_diagonal_fill 3099b54502bSHong Zhang 3109b54502bSHong Zhang Notes: 3119b54502bSHong Zhang Does not apply with 0 fill. 3129b54502bSHong Zhang 3139b54502bSHong Zhang Level: intermediate 3149b54502bSHong Zhang 3159b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 3169b54502bSHong Zhang @*/ 3172401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetAllowDiagonalFill(PC pc) 3189b54502bSHong Zhang { 3199b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 3209b54502bSHong Zhang 3219b54502bSHong Zhang PetscFunctionBegin; 3229b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3232401956bSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 3249b54502bSHong Zhang if (f) { 3259b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 3269b54502bSHong Zhang } 3279b54502bSHong Zhang PetscFunctionReturn(0); 3289b54502bSHong Zhang } 3299b54502bSHong Zhang 3309b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 3319b54502bSHong Zhang 3329b54502bSHong Zhang EXTERN_C_BEGIN 3339b54502bSHong Zhang #undef __FUNCT__ 3342401956bSBarry Smith #define __FUNCT__ "PCFactorSetPivotInBlocks_ILU" 3352401956bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 3369b54502bSHong Zhang { 3379b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 3389b54502bSHong Zhang 3399b54502bSHong Zhang PetscFunctionBegin; 340*075768bcSBarry Smith ((PC_Factor*)dir)->info.pivotinblocks = pivot ? 1.0 : 0.0; 3419b54502bSHong Zhang PetscFunctionReturn(0); 3429b54502bSHong Zhang } 3439b54502bSHong Zhang EXTERN_C_END 3449b54502bSHong Zhang 34562bba022SBarry Smith EXTERN_C_BEGIN 34662bba022SBarry Smith #undef __FUNCT__ 34762bba022SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks_ILU" 34862bba022SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks_ILU(PC pc,PetscReal shift) 34962bba022SBarry Smith { 35062bba022SBarry Smith PC_ILU *dir = (PC_ILU*)pc->data; 35162bba022SBarry Smith 35262bba022SBarry Smith PetscFunctionBegin; 35362bba022SBarry Smith if (shift == PETSC_DEFAULT) { 354*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftinblocks = 1.e-12; 35562bba022SBarry Smith } else { 356*075768bcSBarry Smith ((PC_Factor*)dir)->info.shiftinblocks = shift; 35762bba022SBarry Smith } 35862bba022SBarry Smith PetscFunctionReturn(0); 35962bba022SBarry Smith } 36062bba022SBarry Smith EXTERN_C_END 36162bba022SBarry Smith 36262bba022SBarry Smith 3639b54502bSHong Zhang #undef __FUNCT__ 3649b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 3659b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 3669b54502bSHong Zhang { 3679b54502bSHong Zhang PetscErrorCode ierr; 3689b54502bSHong Zhang PetscInt dtmax = 3,itmp; 3699b54502bSHong Zhang PetscTruth flg,set; 3709b54502bSHong Zhang PetscReal dt[3]; 3719b54502bSHong Zhang char tname[256]; 3729b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 3739b54502bSHong Zhang PetscFList ordlist; 3749b54502bSHong Zhang PetscReal tol; 3759b54502bSHong Zhang 3769b54502bSHong Zhang PetscFunctionBegin; 3779b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 3789b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 379*075768bcSBarry Smith ierr = PetscOptionsInt("-pc_factor_levels","levels of fill","PCFactorSetLevels",(PetscInt)((PC_Factor*)ilu)->info.levels,&itmp,&flg);CHKERRQ(ierr); 380*075768bcSBarry Smith if (flg) ((PC_Factor*)ilu)->info.levels = itmp; 381edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_in_place","do factorization in place","PCFactorSetUseInPlace",&flg);CHKERRQ(ierr); 382edbb7feeSSatish Balay if (flg) ilu->inplace = PETSC_TRUE; 3832401956bSBarry Smith ierr = PetscOptionsName("-pc_factor_diagonal_fill","Allow fill into empty diagonal entry","PCFactorSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 384*075768bcSBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = (double) flg; 385edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_fill","Reuse fill ratio from previous factorization","PCFactorSetReuseFill",&flg);CHKERRQ(ierr); 386edbb7feeSSatish Balay if (flg) ilu->reusefill = PETSC_TRUE; 387edbb7feeSSatish Balay ierr = PetscOptionsName("-pc_factor_reuse_ordering","Reuse previous reordering","PCFactorSetReuseOrdering",&flg);CHKERRQ(ierr); 388edbb7feeSSatish Balay if (flg) ilu->reuseordering = PETSC_TRUE; 3899f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 3909b54502bSHong Zhang if (flg) { 391afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 3929b54502bSHong Zhang } 393*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",((PC_Factor*)ilu)->info.shiftnz,&((PC_Factor*)ilu)->info.shiftnz,0);CHKERRQ(ierr); 394*075768bcSBarry Smith flg = (((PC_Factor*)ilu)->info.shiftpd > 0.0) ? PETSC_TRUE : PETSC_FALSE; 39562bba022SBarry Smith ierr = PetscOptionsTruth("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 39662bba022SBarry Smith ierr = PCFactorSetShiftPd(pc,flg);CHKERRQ(ierr); 397*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",((PC_Factor*)ilu)->info.zeropivot,&((PC_Factor*)ilu)->info.zeropivot,0);CHKERRQ(ierr); 3989b54502bSHong Zhang 399*075768bcSBarry Smith dt[0] = ((PC_Factor*)ilu)->info.dt; 400*075768bcSBarry Smith dt[1] = ((PC_Factor*)ilu)->info.dtcol; 401*075768bcSBarry Smith dt[2] = ((PC_Factor*)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 } 406*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",((PC_Factor*)ilu)->info.fill,&((PC_Factor*)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*075768bcSBarry Smith ierr = PetscOptionsList("-pc_factor_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCFactorSetMatOrderingType",ordlist,((PC_Factor*)ilu)->ordering,tname,256,&flg);CHKERRQ(ierr); 4169b54502bSHong Zhang if (flg) { 417e5a9bf91SBarry Smith ierr = PCFactorSetMatOrderingType(pc,tname);CHKERRQ(ierr); 4189b54502bSHong Zhang } 419*075768bcSBarry Smith flg = ((PC_Factor*)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 } 42462bba022SBarry Smith ierr = PetscOptionsName("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",&flg);CHKERRQ(ierr); 42562bba022SBarry Smith if (flg) { 42662bba022SBarry Smith ierr = PCFactorSetShiftInBlocks(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 42762bba022SBarry Smith } 428*075768bcSBarry Smith ierr = PetscOptionsReal("-pc_factor_shift_in_blocks","Shift added to diagonal of block","PCFactorSetShiftInBlocks",((PC_Factor*)ilu)->info.shiftinblocks,&((PC_Factor*)ilu)->info.shiftinblocks,0);CHKERRQ(ierr); 42962bba022SBarry Smith 4309b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4319b54502bSHong Zhang PetscFunctionReturn(0); 4329b54502bSHong Zhang } 4339b54502bSHong Zhang 4349b54502bSHong Zhang #undef __FUNCT__ 4359b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 4369b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 4379b54502bSHong Zhang { 4389b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 4399b54502bSHong Zhang PetscErrorCode ierr; 4409b54502bSHong Zhang PetscTruth isstring,iascii; 4419b54502bSHong Zhang 4429b54502bSHong Zhang PetscFunctionBegin; 4439b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 4449b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 4459b54502bSHong Zhang if (iascii) { 4469b54502bSHong Zhang if (ilu->usedt) { 447*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",((PC_Factor*)ilu)->info.dt);CHKERRQ(ierr); 448*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)((PC_Factor*)ilu)->info.dtcount);CHKERRQ(ierr); 449*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",((PC_Factor*)ilu)->info.dtcol);CHKERRQ(ierr); 450*075768bcSBarry Smith } else if (((PC_Factor*)ilu)->info.levels == 1) { 451*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr); 4529b54502bSHong Zhang } else { 453*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)((PC_Factor*)ilu)->info.levels);CHKERRQ(ierr); 4549b54502bSHong Zhang } 455*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio allocated %G\n",((PC_Factor*)ilu)->info.fill);CHKERRQ(ierr); 456*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",((PC_Factor*)ilu)->info.zeropivot);CHKERRQ(ierr); 457*075768bcSBarry Smith if (((PC_Factor*)ilu)->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 458*075768bcSBarry Smith if (((PC_Factor*)ilu)->info.shiftnz) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift to prevent zero pivot\n");CHKERRQ(ierr);} 459*075768bcSBarry Smith if (((PC_Factor*)ilu)->info.shiftinblocks) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using diagonal shift on blocks to prevent zero pivot\n");CHKERRQ(ierr);} 4609b54502bSHong Zhang if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 4619b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 462*075768bcSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 4639b54502bSHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 4649b54502bSHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 465*075768bcSBarry Smith if (((PC_Factor*)ilu)->fact) { 466f3a39becSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: factor fill ratio needed %G\n",ilu->actualfill);CHKERRQ(ierr); 4679b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 4689b54502bSHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 469f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 470f3a39becSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4719b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 472*075768bcSBarry Smith ierr = MatView(((PC_Factor*)ilu)->fact,viewer);CHKERRQ(ierr); 4739b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 4749b54502bSHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 475f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 476f3a39becSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4779b54502bSHong Zhang } 4789b54502bSHong Zhang } else if (isstring) { 479*075768bcSBarry Smith ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)((PC_Factor*)ilu)->info.levels,((PC_Factor*)ilu)->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 4809b54502bSHong Zhang } else { 4819b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 4829b54502bSHong Zhang } 4839b54502bSHong Zhang PetscFunctionReturn(0); 4849b54502bSHong Zhang } 4859b54502bSHong Zhang 4869b54502bSHong Zhang #undef __FUNCT__ 4879b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 4889b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 4899b54502bSHong Zhang { 4909b54502bSHong Zhang PetscErrorCode ierr; 4919b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 492f3a39becSBarry Smith MatInfo info; 4939b54502bSHong Zhang 4949b54502bSHong Zhang PetscFunctionBegin; 4959b54502bSHong Zhang if (ilu->inplace) { 496*075768bcSBarry Smith CHKMEMQ; 4979b54502bSHong Zhang if (!pc->setupcalled) { 4989b54502bSHong Zhang 4999b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 5009b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 501*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 502efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 503efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5049b54502bSHong Zhang } 5059b54502bSHong Zhang 5069b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 5079b54502bSHong Zhang cannot have levels of fill */ 508*075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 509*075768bcSBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0; 510*075768bcSBarry Smith ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 511*075768bcSBarry Smith CHKMEMQ; 512*075768bcSBarry Smith ((PC_Factor*)ilu)->fact = pc->pmat; 5139b54502bSHong Zhang } else if (ilu->usedt) { 5149b54502bSHong Zhang if (!pc->setupcalled) { 515*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 516*075768bcSBarry Smith CHKMEMQ; 517efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 518efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 519*075768bcSBarry Smith ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 520*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 5219b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 522*075768bcSBarry Smith CHKMEMQ; 523*075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 524*075768bcSBarry Smith CHKMEMQ; 5259b54502bSHong Zhang if (!ilu->reuseordering) { 5269b54502bSHong Zhang if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 5279b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 528*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 529efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 530efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5319b54502bSHong Zhang } 532*075768bcSBarry Smith ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 533*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 5349b54502bSHong Zhang } else if (!ilu->reusefill) { 535*075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 536*075768bcSBarry Smith ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 537*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 5389b54502bSHong Zhang } else { 539*075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 5409b54502bSHong Zhang } 5419b54502bSHong Zhang } else { 5429b54502bSHong Zhang if (!pc->setupcalled) { 5439b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 544*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 545efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 546efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5479b54502bSHong Zhang /* Remove zeros along diagonal? */ 5489b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5499b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5509b54502bSHong Zhang } 551*075768bcSBarry Smith CHKMEMQ; 552*075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 553*075768bcSBarry Smith CHKMEMQ; 554*075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 555*075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 556f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 557*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 5589b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 5599b54502bSHong Zhang if (!ilu->reuseordering) { 5609b54502bSHong Zhang /* compute a new ordering for the ILU */ 5619b54502bSHong Zhang ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 5629b54502bSHong Zhang ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 563*075768bcSBarry Smith ierr = MatGetOrdering(pc->pmat,((PC_Factor*)ilu)->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 564efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 565efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 5669b54502bSHong Zhang /* Remove zeros along diagonal? */ 5679b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 5689b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 5699b54502bSHong Zhang } 5709b54502bSHong Zhang } 571*075768bcSBarry Smith ierr = MatDestroy(((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 572*075768bcSBarry Smith ierr = MatGetFactor(pc->pmat,MAT_SOLVER_PETSC,MAT_FACTOR_ILU,&((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 573*075768bcSBarry Smith ierr = MatILUFactorSymbolic(((PC_Factor*)ilu)->fact,pc->pmat,ilu->row,ilu->col,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 574*075768bcSBarry Smith ierr = MatGetInfo(((PC_Factor*)ilu)->fact,MAT_LOCAL,&info);CHKERRQ(ierr); 575f3a39becSBarry Smith ilu->actualfill = info.fill_ratio_needed; 576*075768bcSBarry Smith ierr = PetscLogObjectParent(pc,((PC_Factor*)ilu)->fact);CHKERRQ(ierr); 5779b54502bSHong Zhang } 578*075768bcSBarry Smith CHKMEMQ; 579*075768bcSBarry Smith ierr = MatLUFactorNumeric(((PC_Factor*)ilu)->fact,pc->pmat,&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 580*075768bcSBarry Smith CHKMEMQ; 5819b54502bSHong Zhang } 5829b54502bSHong Zhang PetscFunctionReturn(0); 5839b54502bSHong Zhang } 5849b54502bSHong Zhang 5859b54502bSHong Zhang #undef __FUNCT__ 5869b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 5879b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 5889b54502bSHong Zhang { 5899b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5909b54502bSHong Zhang PetscErrorCode ierr; 5919b54502bSHong Zhang 5929b54502bSHong Zhang PetscFunctionBegin; 5939b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 594*075768bcSBarry Smith ierr = PetscStrfree(((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 5959b54502bSHong Zhang ierr = PetscFree(ilu);CHKERRQ(ierr); 5969b54502bSHong Zhang PetscFunctionReturn(0); 5979b54502bSHong Zhang } 5989b54502bSHong Zhang 5999b54502bSHong Zhang #undef __FUNCT__ 6009b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 6019b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 6029b54502bSHong Zhang { 6039b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6049b54502bSHong Zhang PetscErrorCode ierr; 6059b54502bSHong Zhang 6069b54502bSHong Zhang PetscFunctionBegin; 607*075768bcSBarry Smith ierr = MatSolve(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 6089b54502bSHong Zhang PetscFunctionReturn(0); 6099b54502bSHong Zhang } 6109b54502bSHong Zhang 6119b54502bSHong Zhang #undef __FUNCT__ 6129b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 6139b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 6149b54502bSHong Zhang { 6159b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6169b54502bSHong Zhang PetscErrorCode ierr; 6179b54502bSHong Zhang 6189b54502bSHong Zhang PetscFunctionBegin; 619*075768bcSBarry Smith ierr = MatSolveTranspose(((PC_Factor*)ilu)->fact,x,y);CHKERRQ(ierr); 6209b54502bSHong Zhang PetscFunctionReturn(0); 6219b54502bSHong Zhang } 6229b54502bSHong Zhang 6239b54502bSHong Zhang #undef __FUNCT__ 624a4fd02acSBarry Smith #define __FUNCT__ "PCFactorGetMatrix_ILU" 625a4fd02acSBarry Smith static PetscErrorCode PCFactorGetMatrix_ILU(PC pc,Mat *mat) 6269b54502bSHong Zhang { 6279b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6289b54502bSHong Zhang 6299b54502bSHong Zhang PetscFunctionBegin; 630*075768bcSBarry Smith if (!((PC_Factor*)ilu)->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 631*075768bcSBarry Smith *mat = ((PC_Factor*)ilu)->fact; 6329b54502bSHong Zhang PetscFunctionReturn(0); 6339b54502bSHong Zhang } 6349b54502bSHong Zhang 63562bba022SBarry Smith #undef __FUNCT__ 63662bba022SBarry Smith #define __FUNCT__ "PCFactorSetShiftInBlocks" 63762bba022SBarry Smith /*@ 63862bba022SBarry Smith PCFactorSetShiftInBlocks - Shifts the diagonal of any block if LU on that block detects a zero pivot. 63962bba022SBarry Smith 64062bba022SBarry Smith Collective on PC 64162bba022SBarry Smith 64262bba022SBarry Smith Input Parameters: 64362bba022SBarry Smith + pc - the preconditioner context 64462bba022SBarry Smith - shift - amount to shift or PETSC_DECIDE 64562bba022SBarry Smith 64662bba022SBarry Smith Options Database Key: 64762bba022SBarry Smith . -pc_factor_shift_in_blocks <shift> 64862bba022SBarry Smith 64962bba022SBarry Smith Level: intermediate 65062bba022SBarry Smith 65162bba022SBarry Smith .seealso: PCILUSetMatOrdering(), PCFactorSetPivoting(), PCFactorSetShiftInBlocks() 65262bba022SBarry Smith @*/ 65362bba022SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftInBlocks(PC pc,PetscReal shift) 65462bba022SBarry Smith { 65562bba022SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 65662bba022SBarry Smith 65762bba022SBarry Smith PetscFunctionBegin; 65862bba022SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCFactorSetShiftInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 65962bba022SBarry Smith if (f) { 66062bba022SBarry Smith ierr = (*f)(pc,shift);CHKERRQ(ierr); 66162bba022SBarry Smith } 66262bba022SBarry Smith PetscFunctionReturn(0); 66362bba022SBarry Smith } 66462bba022SBarry Smith 6659b54502bSHong Zhang /*MC 6669b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 6679b54502bSHong Zhang 6689b54502bSHong Zhang Options Database Keys: 6692401956bSBarry Smith + -pc_factor_levels <k> - number of levels of fill for ILU(k) 6702401956bSBarry Smith . -pc_factor_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 6719b54502bSHong Zhang its factorization (overwrites original matrix) 6722401956bSBarry Smith . -pc_factor_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 6732401956bSBarry Smith . -pc_factor_reuse_ordering - reuse ordering of factorized matrix from previous factorization 6742401956bSBarry Smith . -pc_factor_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 67555ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 6762401956bSBarry Smith . -pc_factor_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 6779b54502bSHong Zhang this decreases the chance of getting a zero pivot 6782401956bSBarry Smith . -pc_factor_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 6792401956bSBarry Smith . -pc_factor_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 6809b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 6819b54502bSHong Zhang stability of the ILU factorization 68262bba022SBarry Smith . -pc_factor_shift_in_blocks - adds a small diagonal to any block if it is singular during ILU factorization 683f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 684f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 685f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 6869b54502bSHong Zhang 6879b54502bSHong Zhang Level: beginner 6889b54502bSHong Zhang 6899b54502bSHong Zhang Concepts: incomplete factorization 6909b54502bSHong Zhang 691c582cd25SBarry Smith Notes: Only implemented for some matrix formats. (for parallel use you 692c582cd25SBarry Smith must use MATMPIROWBS, see MatCreateMPIRowbs(), this supports only ILU(0) and this is not recommended 693c582cd25SBarry Smith unless you really want a parallel ILU). 6949b54502bSHong Zhang 6959b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 6969b54502bSHong Zhang 697c582cd25SBarry Smith References: 698c582cd25SBarry Smith T. Dupont, R. Kendall, and H. Rachford. An approximate factorization procedure for solving 699c582cd25SBarry Smith self-adjoint elliptic difference equations. SIAM J. Numer. Anal., 5:559--573, 1968. 700c582cd25SBarry Smith 701c582cd25SBarry Smith T.A. Oliphant. An implicit numerical method for solving two-dimensional time-dependent dif- 702c582cd25SBarry Smith fusion problems. Quart. Appl. Math., 19:221--229, 1961. 703c582cd25SBarry Smith 704c582cd25SBarry Smith Review article: APPROXIMATE AND INCOMPLETE FACTORIZATIONS, TONY F. CHAN AND HENK A. VAN DER VORST 705c582cd25SBarry Smith http://igitur-archive.library.uu.nl/math/2001-0621-115821/proc.pdf chapter in Parallel Numerical 706c582cd25SBarry Smith Algorithms, edited by D. Keyes, A. Semah, V. Venkatakrishnan, ICASE/LaRC Interdisciplinary Series in 707c582cd25SBarry Smith Science and Engineering, Kluwer, pp. 167--202. 708c582cd25SBarry Smith 709c582cd25SBarry Smith 7109b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 7112401956bSBarry Smith PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCFactorSetUseDropTolerance(), 712e5a9bf91SBarry Smith PCFactorSetFill(), PCFactorSetMatOrderingType(), PCFactorSetReuseOrdering(), 7132401956bSBarry Smith PCFactorSetLevels(), PCFactorSetUseInPlace(), PCFactorSetAllowDiagonalFill(), PCFactorSetPivotInBlocks(), 714f251bdbdSHong Zhang PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 7159b54502bSHong Zhang 7169b54502bSHong Zhang M*/ 7179b54502bSHong Zhang 7189b54502bSHong Zhang EXTERN_C_BEGIN 7199b54502bSHong Zhang #undef __FUNCT__ 7209b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 721dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 7229b54502bSHong Zhang { 7239b54502bSHong Zhang PetscErrorCode ierr; 7249b54502bSHong Zhang PC_ILU *ilu; 7259b54502bSHong Zhang 7269b54502bSHong Zhang PetscFunctionBegin; 72738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ILU,&ilu);CHKERRQ(ierr); 7289b54502bSHong Zhang 729*075768bcSBarry Smith ((PC_Factor*)ilu)->fact = 0; 730*075768bcSBarry Smith ierr = MatFactorInfoInitialize(&((PC_Factor*)ilu)->info);CHKERRQ(ierr); 731*075768bcSBarry Smith ((PC_Factor*)ilu)->info.levels = 0; 732*075768bcSBarry Smith ((PC_Factor*)ilu)->info.fill = 1.0; 7339b54502bSHong Zhang ilu->col = 0; 7349b54502bSHong Zhang ilu->row = 0; 7359b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 736*075768bcSBarry Smith ierr = PetscStrallocpy(MATORDERING_NATURAL,&((PC_Factor*)ilu)->ordering);CHKERRQ(ierr); 7379b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 7389b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 739*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dt = PETSC_DEFAULT; 740*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcount = PETSC_DEFAULT; 741*075768bcSBarry Smith ((PC_Factor*)ilu)->info.dtcol = PETSC_DEFAULT; 742*075768bcSBarry Smith ((PC_Factor*)ilu)->info.shiftnz = 1.e-12; 743*075768bcSBarry Smith ((PC_Factor*)ilu)->info.shiftpd = 0.0; /* false */ 744*075768bcSBarry Smith ((PC_Factor*)ilu)->info.zeropivot = 1.e-12; 745*075768bcSBarry Smith ((PC_Factor*)ilu)->info.pivotinblocks = 1.0; 746*075768bcSBarry Smith ((PC_Factor*)ilu)->info.shiftinblocks = 1.e-12; 7479b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 748*075768bcSBarry Smith ((PC_Factor*)ilu)->info.diagonal_fill = 0; 7499b54502bSHong Zhang pc->data = (void*)ilu; 7509b54502bSHong Zhang 7519b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 7529b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 7539b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 7549b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 7559b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 756a4fd02acSBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_ILU; 7579b54502bSHong Zhang pc->ops->view = PCView_ILU; 7589b54502bSHong Zhang pc->ops->applyrichardson = 0; 7599b54502bSHong Zhang 760afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 761afaefe49SHong Zhang PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 762afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 763afaefe49SHong Zhang PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 764afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 765afaefe49SHong Zhang PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 766afaefe49SHong Zhang 7672401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseDropTolerance_C","PCFactorSetUseDropTolerance_ILU", 7682401956bSBarry Smith PCFactorSetUseDropTolerance_ILU);CHKERRQ(ierr); 76955ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 77055ba2a51SBarry Smith PCFactorSetFill_ILU);CHKERRQ(ierr); 771e5a9bf91SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetMatOrderingType_C","PCFactorSetMatOrderingType_ILU", 772e5a9bf91SBarry Smith PCFactorSetMatOrderingType_ILU);CHKERRQ(ierr); 7732401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseOrdering_C","PCFactorSetReuseOrdering_ILU", 7742401956bSBarry Smith PCFactorSetReuseOrdering_ILU);CHKERRQ(ierr); 7752401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetReuseFill_C","PCFactorSetReuseFill_ILU", 7762401956bSBarry Smith PCFactorSetReuseFill_ILU);CHKERRQ(ierr); 7772401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetLevels_C","PCFactorSetLevels_ILU", 7782401956bSBarry Smith PCFactorSetLevels_ILU);CHKERRQ(ierr); 7792401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetUseInPlace_C","PCFactorSetUseInPlace_ILU", 7802401956bSBarry Smith PCFactorSetUseInPlace_ILU);CHKERRQ(ierr); 7812401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C","PCFactorSetAllowDiagonalFill_ILU", 7822401956bSBarry Smith PCFactorSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 7832401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetPivotInBlocks_C","PCFactorSetPivotInBlocks_ILU", 7842401956bSBarry Smith PCFactorSetPivotInBlocks_ILU);CHKERRQ(ierr); 78562bba022SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftInBlocks_C","PCFactorSetShiftInBlocks_ILU", 78662bba022SBarry Smith PCFactorSetShiftInBlocks_ILU);CHKERRQ(ierr); 7872401956bSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C","PCFactorReorderForNonzeroDiagonal_ILU", 7882401956bSBarry Smith PCFactorReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 7899b54502bSHong Zhang PetscFunctionReturn(0); 7909b54502bSHong Zhang } 7919b54502bSHong Zhang EXTERN_C_END 792