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" 89b54502bSHong Zhang #include "src/mat/matimpl.h" 99b54502bSHong Zhang 109b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 11afaefe49SHong Zhang EXTERN_C_BEGIN 12afaefe49SHong Zhang #undef __FUNCT__ 13afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetZeroPivot_ILU" 14dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetZeroPivot_ILU(PC pc,PetscReal z) 15afaefe49SHong Zhang { 16afaefe49SHong Zhang PC_ILU *ilu; 17afaefe49SHong Zhang 18afaefe49SHong Zhang PetscFunctionBegin; 19afaefe49SHong Zhang ilu = (PC_ILU*)pc->data; 20afaefe49SHong Zhang ilu->info.zeropivot = z; 21afaefe49SHong Zhang PetscFunctionReturn(0); 22afaefe49SHong Zhang } 23afaefe49SHong Zhang EXTERN_C_END 24afaefe49SHong Zhang 25afaefe49SHong Zhang EXTERN_C_BEGIN 26afaefe49SHong Zhang #undef __FUNCT__ 27afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftNonzero_ILU" 28dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftNonzero_ILU(PC pc,PetscReal shift) 29afaefe49SHong Zhang { 30afaefe49SHong Zhang PC_ILU *dir; 31afaefe49SHong Zhang 32afaefe49SHong Zhang PetscFunctionBegin; 33afaefe49SHong Zhang dir = (PC_ILU*)pc->data; 34afaefe49SHong Zhang if (shift == (PetscReal) PETSC_DECIDE) { 35afaefe49SHong Zhang dir->info.shiftnz = 1.e-12; 36afaefe49SHong Zhang } else { 37afaefe49SHong Zhang dir->info.shiftnz = shift; 38afaefe49SHong Zhang } 39afaefe49SHong Zhang PetscFunctionReturn(0); 40afaefe49SHong Zhang } 41afaefe49SHong Zhang EXTERN_C_END 42afaefe49SHong Zhang 43afaefe49SHong Zhang EXTERN_C_BEGIN 44afaefe49SHong Zhang #undef __FUNCT__ 45afaefe49SHong Zhang #define __FUNCT__ "PCFactorSetShiftPd_ILU" 46dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetShiftPd_ILU(PC pc,PetscTruth shift) 47afaefe49SHong Zhang { 48afaefe49SHong Zhang PC_ILU *dir; 49afaefe49SHong Zhang 50afaefe49SHong Zhang PetscFunctionBegin; 51afaefe49SHong Zhang dir = (PC_ILU*)pc->data; 52fbf22428SSatish Balay if (shift) { 53fbf22428SSatish Balay dir->info.shift_fraction = 0.0; 54fbf22428SSatish Balay dir->info.shiftpd = 1.0; 55fbf22428SSatish Balay } else { 56fbf22428SSatish Balay dir->info.shiftpd = 0.0; 57fbf22428SSatish Balay } 58afaefe49SHong Zhang PetscFunctionReturn(0); 59afaefe49SHong Zhang } 60afaefe49SHong Zhang EXTERN_C_END 619b54502bSHong Zhang 629b54502bSHong Zhang EXTERN_C_BEGIN 639b54502bSHong Zhang #undef __FUNCT__ 649b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal_ILU" 65dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 669b54502bSHong Zhang { 679b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 689b54502bSHong Zhang 699b54502bSHong Zhang PetscFunctionBegin; 709b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 719b54502bSHong Zhang if (z == PETSC_DECIDE) { 729b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = 1.e-10; 739b54502bSHong Zhang } else { 749b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = z; 759b54502bSHong Zhang } 769b54502bSHong Zhang PetscFunctionReturn(0); 779b54502bSHong Zhang } 789b54502bSHong Zhang EXTERN_C_END 799b54502bSHong Zhang 809b54502bSHong Zhang #undef __FUNCT__ 819b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal" 829b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc) 839b54502bSHong Zhang { 849b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 859b54502bSHong Zhang PetscErrorCode ierr; 869b54502bSHong Zhang 879b54502bSHong Zhang PetscFunctionBegin; 889b54502bSHong Zhang if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);} 899b54502bSHong Zhang if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 909b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 919b54502bSHong Zhang PetscFunctionReturn(0); 929b54502bSHong Zhang } 939b54502bSHong Zhang 949b54502bSHong Zhang EXTERN_C_BEGIN 959b54502bSHong Zhang #undef __FUNCT__ 969b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance_ILU" 97dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 989b54502bSHong Zhang { 999b54502bSHong Zhang PC_ILU *ilu; 1009b54502bSHong Zhang PetscErrorCode ierr; 1019b54502bSHong Zhang 1029b54502bSHong Zhang PetscFunctionBegin; 1039b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1049b54502bSHong Zhang if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) { 1059b54502bSHong Zhang pc->setupcalled = 0; 1069b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1079b54502bSHong Zhang } 1089b54502bSHong Zhang ilu->usedt = PETSC_TRUE; 1099b54502bSHong Zhang ilu->info.dt = dt; 1109b54502bSHong Zhang ilu->info.dtcol = dtcol; 1119b54502bSHong Zhang ilu->info.dtcount = dtcount; 1129b54502bSHong Zhang ilu->info.fill = PETSC_DEFAULT; 1139b54502bSHong Zhang PetscFunctionReturn(0); 1149b54502bSHong Zhang } 1159b54502bSHong Zhang EXTERN_C_END 1169b54502bSHong Zhang 1179b54502bSHong Zhang EXTERN_C_BEGIN 1189b54502bSHong Zhang #undef __FUNCT__ 119*55ba2a51SBarry Smith #define __FUNCT__ "PCFactorSetFill_ILU" 120*55ba2a51SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCFactorSetFill_ILU(PC pc,PetscReal fill) 1219b54502bSHong Zhang { 1229b54502bSHong Zhang PC_ILU *dir; 1239b54502bSHong Zhang 1249b54502bSHong Zhang PetscFunctionBegin; 1259b54502bSHong Zhang dir = (PC_ILU*)pc->data; 1269b54502bSHong Zhang dir->info.fill = fill; 1279b54502bSHong Zhang PetscFunctionReturn(0); 1289b54502bSHong Zhang } 1299b54502bSHong Zhang EXTERN_C_END 1309b54502bSHong Zhang 1319b54502bSHong Zhang EXTERN_C_BEGIN 1329b54502bSHong Zhang #undef __FUNCT__ 1339b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering_ILU" 134dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering) 1359b54502bSHong Zhang { 1369b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 1379b54502bSHong Zhang PetscErrorCode ierr; 1389b54502bSHong Zhang PetscTruth flg; 1399b54502bSHong Zhang 1409b54502bSHong Zhang PetscFunctionBegin; 1419b54502bSHong Zhang if (!pc->setupcalled) { 1429b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1439b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1449b54502bSHong Zhang } else { 1459b54502bSHong Zhang ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 1469b54502bSHong Zhang if (!flg) { 1479b54502bSHong Zhang pc->setupcalled = 0; 1489b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1499b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1509b54502bSHong Zhang /* free the data structures, then create them again */ 1519b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1529b54502bSHong Zhang } 1539b54502bSHong Zhang } 1549b54502bSHong Zhang PetscFunctionReturn(0); 1559b54502bSHong Zhang } 1569b54502bSHong Zhang EXTERN_C_END 1579b54502bSHong Zhang 1589b54502bSHong Zhang EXTERN_C_BEGIN 1599b54502bSHong Zhang #undef __FUNCT__ 1609b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering_ILU" 161dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag) 1629b54502bSHong Zhang { 1639b54502bSHong Zhang PC_ILU *ilu; 1649b54502bSHong Zhang 1659b54502bSHong Zhang PetscFunctionBegin; 1669b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1679b54502bSHong Zhang ilu->reuseordering = flag; 1689b54502bSHong Zhang PetscFunctionReturn(0); 1699b54502bSHong Zhang } 1709b54502bSHong Zhang EXTERN_C_END 1719b54502bSHong Zhang 1729b54502bSHong Zhang EXTERN_C_BEGIN 1739b54502bSHong Zhang #undef __FUNCT__ 1749b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill_ILUDT" 175dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUDTSetReuseFill_ILUDT(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->reusefill = flag; 1829b54502bSHong Zhang if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported"); 1839b54502bSHong Zhang PetscFunctionReturn(0); 1849b54502bSHong Zhang } 1859b54502bSHong Zhang EXTERN_C_END 1869b54502bSHong Zhang 1879b54502bSHong Zhang EXTERN_C_BEGIN 1889b54502bSHong Zhang #undef __FUNCT__ 1899b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels_ILU" 190dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetLevels_ILU(PC pc,PetscInt levels) 1919b54502bSHong Zhang { 1929b54502bSHong Zhang PC_ILU *ilu; 1939b54502bSHong Zhang PetscErrorCode ierr; 1949b54502bSHong Zhang 1959b54502bSHong Zhang PetscFunctionBegin; 1969b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1979b54502bSHong Zhang 1989b54502bSHong Zhang if (!pc->setupcalled) { 1999b54502bSHong Zhang ilu->info.levels = levels; 2009b54502bSHong Zhang } else if (ilu->usedt || ilu->info.levels != levels) { 2019b54502bSHong Zhang ilu->info.levels = levels; 2029b54502bSHong Zhang pc->setupcalled = 0; 2039b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 2049b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 2059b54502bSHong Zhang } 2069b54502bSHong Zhang PetscFunctionReturn(0); 2079b54502bSHong Zhang } 2089b54502bSHong Zhang EXTERN_C_END 2099b54502bSHong Zhang 2109b54502bSHong Zhang EXTERN_C_BEGIN 2119b54502bSHong Zhang #undef __FUNCT__ 2129b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace_ILU" 213dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetUseInPlace_ILU(PC pc) 2149b54502bSHong Zhang { 2159b54502bSHong Zhang PC_ILU *dir; 2169b54502bSHong Zhang 2179b54502bSHong Zhang PetscFunctionBegin; 2189b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2199b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2209b54502bSHong Zhang PetscFunctionReturn(0); 2219b54502bSHong Zhang } 2229b54502bSHong Zhang EXTERN_C_END 2239b54502bSHong Zhang 2249b54502bSHong Zhang EXTERN_C_BEGIN 2259b54502bSHong Zhang #undef __FUNCT__ 2269b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill" 227dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetAllowDiagonalFill_ILU(PC pc) 2289b54502bSHong Zhang { 2299b54502bSHong Zhang PC_ILU *dir; 2309b54502bSHong Zhang 2319b54502bSHong Zhang PetscFunctionBegin; 2329b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2339b54502bSHong Zhang dir->info.diagonal_fill = 1; 2349b54502bSHong Zhang PetscFunctionReturn(0); 2359b54502bSHong Zhang } 2369b54502bSHong Zhang EXTERN_C_END 2379b54502bSHong Zhang 2389b54502bSHong Zhang #undef __FUNCT__ 2399b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance" 2409b54502bSHong Zhang /*@ 2419b54502bSHong Zhang PCILUSetUseDropTolerance - The preconditioner will use an ILU 2429b54502bSHong Zhang based on a drop tolerance. 2439b54502bSHong Zhang 2449b54502bSHong Zhang Collective on PC 2459b54502bSHong Zhang 2469b54502bSHong Zhang Input Parameters: 2479b54502bSHong Zhang + pc - the preconditioner context 2489b54502bSHong Zhang . dt - the drop tolerance, try from 1.e-10 to .1 2499b54502bSHong Zhang . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 2509b54502bSHong Zhang - maxrowcount - the max number of nonzeros allowed in a row, best value 2519b54502bSHong Zhang depends on the number of nonzeros in row of original matrix 2529b54502bSHong Zhang 2539b54502bSHong Zhang Options Database Key: 2549b54502bSHong Zhang . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 2559b54502bSHong Zhang 2569b54502bSHong Zhang Level: intermediate 2579b54502bSHong Zhang 2589b54502bSHong Zhang Notes: 2599b54502bSHong Zhang This uses the iludt() code of Saad's SPARSKIT package 2609b54502bSHong Zhang 2619b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 2629b54502bSHong Zhang @*/ 263dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 2649b54502bSHong Zhang { 2659b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 2669b54502bSHong Zhang 2679b54502bSHong Zhang PetscFunctionBegin; 2689b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2699b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 2709b54502bSHong Zhang if (f) { 2719b54502bSHong Zhang ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 2729b54502bSHong Zhang } 2739b54502bSHong Zhang PetscFunctionReturn(0); 2749b54502bSHong Zhang } 2759b54502bSHong Zhang 2769b54502bSHong Zhang #undef __FUNCT__ 2779b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering" 2789b54502bSHong Zhang /*@C 2799b54502bSHong Zhang PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to 2809b54502bSHong Zhang be used in the ILU factorization. 2819b54502bSHong Zhang 2829b54502bSHong Zhang Collective on PC 2839b54502bSHong Zhang 2849b54502bSHong Zhang Input Parameters: 2859b54502bSHong Zhang + pc - the preconditioner context 2869b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 2879b54502bSHong Zhang 2889b54502bSHong Zhang Options Database Key: 2899b54502bSHong Zhang . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 2909b54502bSHong Zhang 2919b54502bSHong Zhang Level: intermediate 2929b54502bSHong Zhang 2939b54502bSHong Zhang Notes: natural ordering is used by default 2949b54502bSHong Zhang 2959b54502bSHong Zhang .seealso: PCLUSetMatOrdering() 2969b54502bSHong Zhang 2979b54502bSHong Zhang .keywords: PC, ILU, set, matrix, reordering 2989b54502bSHong Zhang 2999b54502bSHong Zhang @*/ 300dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetMatOrdering(PC pc,MatOrderingType ordering) 3019b54502bSHong Zhang { 3029b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 3039b54502bSHong Zhang 3049b54502bSHong Zhang PetscFunctionBegin; 3059b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3069b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 3079b54502bSHong Zhang if (f) { 3089b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 3099b54502bSHong Zhang } 3109b54502bSHong Zhang PetscFunctionReturn(0); 3119b54502bSHong Zhang } 3129b54502bSHong Zhang 3139b54502bSHong Zhang #undef __FUNCT__ 3149b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering" 3159b54502bSHong Zhang /*@ 3169b54502bSHong Zhang PCILUSetReuseOrdering - When similar matrices are factored, this 3179b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 3189b54502bSHong Zhang following factors; applies to both fill and drop tolerance ILUs. 3199b54502bSHong Zhang 3209b54502bSHong Zhang Collective on PC 3219b54502bSHong Zhang 3229b54502bSHong Zhang Input Parameters: 3239b54502bSHong Zhang + pc - the preconditioner context 3249b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 3259b54502bSHong Zhang 3269b54502bSHong Zhang Options Database Key: 3279b54502bSHong Zhang . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering() 3289b54502bSHong Zhang 3299b54502bSHong Zhang Level: intermediate 3309b54502bSHong Zhang 3319b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 3329b54502bSHong Zhang 3339b54502bSHong Zhang .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 3349b54502bSHong Zhang @*/ 335dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetReuseOrdering(PC pc,PetscTruth flag) 3369b54502bSHong Zhang { 3379b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 3389b54502bSHong Zhang 3399b54502bSHong Zhang PetscFunctionBegin; 3409b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3419b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 3429b54502bSHong Zhang if (f) { 3439b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 3449b54502bSHong Zhang } 3459b54502bSHong Zhang PetscFunctionReturn(0); 3469b54502bSHong Zhang } 3479b54502bSHong Zhang 3489b54502bSHong Zhang #undef __FUNCT__ 3499b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill" 3509b54502bSHong Zhang /*@ 3519b54502bSHong Zhang PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored, 3529b54502bSHong Zhang this causes later ones to use the fill computed in the initial factorization. 3539b54502bSHong Zhang 3549b54502bSHong Zhang Collective on PC 3559b54502bSHong Zhang 3569b54502bSHong Zhang Input Parameters: 3579b54502bSHong Zhang + pc - the preconditioner context 3589b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 3599b54502bSHong Zhang 3609b54502bSHong Zhang Options Database Key: 3619b54502bSHong Zhang . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill() 3629b54502bSHong Zhang 3639b54502bSHong Zhang Level: intermediate 3649b54502bSHong Zhang 3659b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 3669b54502bSHong Zhang 3679b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 3689b54502bSHong Zhang @*/ 369dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUDTSetReuseFill(PC pc,PetscTruth flag) 3709b54502bSHong Zhang { 3719b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 3729b54502bSHong Zhang 3739b54502bSHong Zhang PetscFunctionBegin; 3749b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3759b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 3769b54502bSHong Zhang if (f) { 3779b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 3789b54502bSHong Zhang } 3799b54502bSHong Zhang PetscFunctionReturn(0); 3809b54502bSHong Zhang } 3819b54502bSHong Zhang 3829b54502bSHong Zhang #undef __FUNCT__ 3839b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels" 3849b54502bSHong Zhang /*@ 3859b54502bSHong Zhang PCILUSetLevels - Sets the number of levels of fill to use. 3869b54502bSHong Zhang 3879b54502bSHong Zhang Collective on PC 3889b54502bSHong Zhang 3899b54502bSHong Zhang Input Parameters: 3909b54502bSHong Zhang + pc - the preconditioner context 3919b54502bSHong Zhang - levels - number of levels of fill 3929b54502bSHong Zhang 3939b54502bSHong Zhang Options Database Key: 3949b54502bSHong Zhang . -pc_ilu_levels <levels> - Sets fill level 3959b54502bSHong Zhang 3969b54502bSHong Zhang Level: intermediate 3979b54502bSHong Zhang 3989b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 3999b54502bSHong Zhang @*/ 400dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetLevels(PC pc,PetscInt levels) 4019b54502bSHong Zhang { 4029b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 4039b54502bSHong Zhang 4049b54502bSHong Zhang PetscFunctionBegin; 4059b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4069b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 4079b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 4089b54502bSHong Zhang if (f) { 4099b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 4109b54502bSHong Zhang } 4119b54502bSHong Zhang PetscFunctionReturn(0); 4129b54502bSHong Zhang } 4139b54502bSHong Zhang 4149b54502bSHong Zhang #undef __FUNCT__ 4159b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill" 4169b54502bSHong Zhang /*@ 4179b54502bSHong Zhang PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be 4189b54502bSHong Zhang treated as level 0 fill even if there is no non-zero location. 4199b54502bSHong Zhang 4209b54502bSHong Zhang Collective on PC 4219b54502bSHong Zhang 4229b54502bSHong Zhang Input Parameters: 4239b54502bSHong Zhang + pc - the preconditioner context 4249b54502bSHong Zhang 4259b54502bSHong Zhang Options Database Key: 4269b54502bSHong Zhang . -pc_ilu_diagonal_fill 4279b54502bSHong Zhang 4289b54502bSHong Zhang Notes: 4299b54502bSHong Zhang Does not apply with 0 fill. 4309b54502bSHong Zhang 4319b54502bSHong Zhang Level: intermediate 4329b54502bSHong Zhang 4339b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 4349b54502bSHong Zhang @*/ 435dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetAllowDiagonalFill(PC pc) 4369b54502bSHong Zhang { 4379b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 4389b54502bSHong Zhang 4399b54502bSHong Zhang PetscFunctionBegin; 4409b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4419b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 4429b54502bSHong Zhang if (f) { 4439b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 4449b54502bSHong Zhang } 4459b54502bSHong Zhang PetscFunctionReturn(0); 4469b54502bSHong Zhang } 4479b54502bSHong Zhang 4489b54502bSHong Zhang #undef __FUNCT__ 4499b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace" 4509b54502bSHong Zhang /*@ 4519b54502bSHong Zhang PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization. 4529b54502bSHong Zhang Collective on PC 4539b54502bSHong Zhang 4549b54502bSHong Zhang Input Parameters: 4559b54502bSHong Zhang . pc - the preconditioner context 4569b54502bSHong Zhang 4579b54502bSHong Zhang Options Database Key: 4589b54502bSHong Zhang . -pc_ilu_in_place - Activates in-place factorization 4599b54502bSHong Zhang 4609b54502bSHong Zhang Notes: 4619b54502bSHong Zhang PCILUSetUseInPlace() is intended for use with matrix-free variants of 4629b54502bSHong Zhang Krylov methods, or when a different matrices are employed for the linear 4639b54502bSHong Zhang system and preconditioner, or with ASM preconditioning. Do NOT use 4649b54502bSHong Zhang this option if the linear system 4659b54502bSHong Zhang matrix also serves as the preconditioning matrix, since the factored 4669b54502bSHong Zhang matrix would then overwrite the original matrix. 4679b54502bSHong Zhang 4689b54502bSHong Zhang Only works well with ILU(0). 4699b54502bSHong Zhang 4709b54502bSHong Zhang Level: intermediate 4719b54502bSHong Zhang 4729b54502bSHong Zhang .keywords: PC, set, factorization, inplace, in-place, ILU 4739b54502bSHong Zhang 4749b54502bSHong Zhang .seealso: PCLUSetUseInPlace() 4759b54502bSHong Zhang @*/ 476dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetUseInPlace(PC pc) 4779b54502bSHong Zhang { 4789b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 4799b54502bSHong Zhang 4809b54502bSHong Zhang PetscFunctionBegin; 4819b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4829b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 4839b54502bSHong Zhang if (f) { 4849b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 4859b54502bSHong Zhang } 4869b54502bSHong Zhang PetscFunctionReturn(0); 4879b54502bSHong Zhang } 4889b54502bSHong Zhang 4899b54502bSHong Zhang #undef __FUNCT__ 4909b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal" 4919b54502bSHong Zhang /*@ 4929b54502bSHong Zhang PCILUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 4939b54502bSHong Zhang 4949b54502bSHong Zhang Collective on PC 4959b54502bSHong Zhang 4969b54502bSHong Zhang Input Parameters: 4979b54502bSHong Zhang + pc - the preconditioner context 4989b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 4999b54502bSHong Zhang 5009b54502bSHong Zhang Options Database Key: 5019b54502bSHong Zhang . -pc_lu_nonzeros_along_diagonal 5029b54502bSHong Zhang 5039b54502bSHong Zhang Level: intermediate 5049b54502bSHong Zhang 5059b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 5069b54502bSHong Zhang 507*55ba2a51SBarry Smith .seealso: PCFactorSetFill(), MatReorderForNonzeroDiagonal() 5089b54502bSHong Zhang @*/ 509dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 5109b54502bSHong Zhang { 5119b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 5129b54502bSHong Zhang 5139b54502bSHong Zhang PetscFunctionBegin; 5149b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5159b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 5169b54502bSHong Zhang if (f) { 5179b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 5189b54502bSHong Zhang } 5199b54502bSHong Zhang PetscFunctionReturn(0); 5209b54502bSHong Zhang } 5219b54502bSHong Zhang 5229b54502bSHong Zhang #undef __FUNCT__ 5239b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks" 5249b54502bSHong Zhang /*@ 5259b54502bSHong Zhang PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block 5269b54502bSHong Zhang with BAIJ or SBAIJ matrices 5279b54502bSHong Zhang 5289b54502bSHong Zhang Collective on PC 5299b54502bSHong Zhang 5309b54502bSHong Zhang Input Parameters: 5319b54502bSHong Zhang + pc - the preconditioner context 5329b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 5339b54502bSHong Zhang 5349b54502bSHong Zhang Options Database Key: 5359b54502bSHong Zhang . -pc_ilu_pivot_in_blocks <true,false> 5369b54502bSHong Zhang 5379b54502bSHong Zhang Level: intermediate 5389b54502bSHong Zhang 5399b54502bSHong Zhang .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting() 5409b54502bSHong Zhang @*/ 541dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetPivotInBlocks(PC pc,PetscTruth pivot) 5429b54502bSHong Zhang { 5439b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 5449b54502bSHong Zhang 5459b54502bSHong Zhang PetscFunctionBegin; 5469b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 5479b54502bSHong Zhang if (f) { 5489b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 5499b54502bSHong Zhang } 5509b54502bSHong Zhang PetscFunctionReturn(0); 5519b54502bSHong Zhang } 5529b54502bSHong Zhang 5539b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 5549b54502bSHong Zhang 5559b54502bSHong Zhang EXTERN_C_BEGIN 5569b54502bSHong Zhang #undef __FUNCT__ 5579b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks_ILU" 558dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 5599b54502bSHong Zhang { 5609b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 5619b54502bSHong Zhang 5629b54502bSHong Zhang PetscFunctionBegin; 5639b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 5649b54502bSHong Zhang PetscFunctionReturn(0); 5659b54502bSHong Zhang } 5669b54502bSHong Zhang EXTERN_C_END 5679b54502bSHong Zhang 5689b54502bSHong Zhang #undef __FUNCT__ 5699b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 5709b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 5719b54502bSHong Zhang { 5729b54502bSHong Zhang PetscErrorCode ierr; 5739b54502bSHong Zhang PetscInt dtmax = 3,itmp; 5749b54502bSHong Zhang PetscTruth flg,set; 5759b54502bSHong Zhang PetscReal dt[3]; 5769b54502bSHong Zhang char tname[256]; 5779b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 5789b54502bSHong Zhang PetscFList ordlist; 5799b54502bSHong Zhang PetscReal tol; 5809b54502bSHong Zhang 5819b54502bSHong Zhang PetscFunctionBegin; 5829b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 5839b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 5849b54502bSHong Zhang ierr = PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 5859b54502bSHong Zhang if (flg) ilu->info.levels = itmp; 5869b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);CHKERRQ(ierr); 5879b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 5889b54502bSHong Zhang ilu->info.diagonal_fill = (double) flg; 5899b54502bSHong Zhang ierr = PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);CHKERRQ(ierr); 5909b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr); 5915e8efad8SHong Zhang 5929f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 5939b54502bSHong Zhang if (flg) { 594afaefe49SHong Zhang ierr = PCFactorSetShiftNonzero(pc,(PetscReal)PETSC_DECIDE);CHKERRQ(ierr); 5959b54502bSHong Zhang } 5969f95998fSHong Zhang ierr = PetscOptionsReal("-pc_factor_shift_nonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.shiftnz,&ilu->info.shiftnz,0);CHKERRQ(ierr); 5975e8efad8SHong Zhang 5989f95998fSHong Zhang ierr = PetscOptionsName("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",&flg);CHKERRQ(ierr); 5999b54502bSHong Zhang if (flg) { 6009f95998fSHong Zhang ierr = PetscOptionsInt("-pc_factor_shift_positive_definite","Manteuffel shift applied to diagonal","PCFactorSetShiftPd",(PetscInt)ilu->info.shiftpd,&itmp,&flg); CHKERRQ(ierr); 6019b54502bSHong Zhang if (flg && !itmp) { 602afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_FALSE);CHKERRQ(ierr); 6039b54502bSHong Zhang } else { 604afaefe49SHong Zhang ierr = PCFactorSetShiftPd(pc,PETSC_TRUE);CHKERRQ(ierr); 6059b54502bSHong Zhang } 6069b54502bSHong Zhang } 607ee45ca4aSHong Zhang ierr = PetscOptionsReal("-pc_factor_zeropivot","Pivot is considered zero if less than","PCFactorSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 6089b54502bSHong Zhang 6099b54502bSHong Zhang dt[0] = ilu->info.dt; 6109b54502bSHong Zhang dt[1] = ilu->info.dtcol; 6119b54502bSHong Zhang dt[2] = ilu->info.dtcount; 6129b54502bSHong Zhang ierr = PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 6139b54502bSHong Zhang if (flg) { 6149b54502bSHong Zhang ierr = PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 6159b54502bSHong Zhang } 616*55ba2a51SBarry Smith ierr = PetscOptionsReal("-pc_factor_fill","Expected fill in factorization","PCFactorSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 6179b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 6189b54502bSHong Zhang if (flg) { 6199b54502bSHong Zhang tol = PETSC_DECIDE; 6209b54502bSHong Zhang ierr = PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 6219b54502bSHong Zhang ierr = PCILUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 6229b54502bSHong Zhang } 6239b54502bSHong Zhang 6249b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 6259b54502bSHong Zhang ierr = PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 6269b54502bSHong Zhang if (flg) { 6279b54502bSHong Zhang ierr = PCILUSetMatOrdering(pc,tname);CHKERRQ(ierr); 6289b54502bSHong Zhang } 6299b54502bSHong Zhang flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 6304dc4c822SBarry Smith ierr = PetscOptionsTruth("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 6319b54502bSHong Zhang if (set) { 6329b54502bSHong Zhang ierr = PCILUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 6339b54502bSHong Zhang } 6349b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 6359b54502bSHong Zhang PetscFunctionReturn(0); 6369b54502bSHong Zhang } 6379b54502bSHong Zhang 6389b54502bSHong Zhang #undef __FUNCT__ 6399b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 6409b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 6419b54502bSHong Zhang { 6429b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6439b54502bSHong Zhang PetscErrorCode ierr; 6449b54502bSHong Zhang PetscTruth isstring,iascii; 6459b54502bSHong Zhang 6469b54502bSHong Zhang PetscFunctionBegin; 6479b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 6489b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 6499b54502bSHong Zhang if (iascii) { 6509b54502bSHong Zhang if (ilu->usedt) { 651a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %G\n",ilu->info.dt);CHKERRQ(ierr); 6529b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 653a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %G\n",ilu->info.dtcol);CHKERRQ(ierr); 6549b54502bSHong Zhang } else if (ilu->info.levels == 1) { 6559b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 6569b54502bSHong Zhang } else { 6579b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 6589b54502bSHong Zhang } 659a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %G\n",ilu->info.fill);CHKERRQ(ierr); 660a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %G\n",ilu->info.zeropivot);CHKERRQ(ierr); 6610a29876aSHong Zhang if (ilu->info.shiftpd) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 6629b54502bSHong Zhang if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 6639b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 6649b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 6659b54502bSHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 6669b54502bSHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 6679b54502bSHong Zhang if (ilu->fact) { 6689b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 6699b54502bSHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6709b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 6719b54502bSHong Zhang ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 6729b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 6739b54502bSHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 6749b54502bSHong Zhang } 6759b54502bSHong Zhang } else if (isstring) { 6769b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 6779b54502bSHong Zhang } else { 6789b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 6799b54502bSHong Zhang } 6809b54502bSHong Zhang PetscFunctionReturn(0); 6819b54502bSHong Zhang } 6829b54502bSHong Zhang 6839b54502bSHong Zhang #undef __FUNCT__ 6849b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 6859b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 6869b54502bSHong Zhang { 6879b54502bSHong Zhang PetscErrorCode ierr; 6889b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 6899b54502bSHong Zhang 6909b54502bSHong Zhang PetscFunctionBegin; 6919b54502bSHong Zhang if (ilu->inplace) { 6929b54502bSHong Zhang if (!pc->setupcalled) { 6939b54502bSHong Zhang 6949b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 6959b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 6969b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 697efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 698efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 6999b54502bSHong Zhang } 7009b54502bSHong Zhang 7019b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 7029b54502bSHong Zhang cannot have levels of fill */ 7039b54502bSHong Zhang ilu->info.fill = 1.0; 7049b54502bSHong Zhang ilu->info.diagonal_fill = 0; 7059b54502bSHong Zhang ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 7069b54502bSHong Zhang ilu->fact = pc->pmat; 7079b54502bSHong Zhang } else if (ilu->usedt) { 7089b54502bSHong Zhang if (!pc->setupcalled) { 7099b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 710efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 711efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 7127529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 71352e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 7149b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 7159b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 7169b54502bSHong Zhang if (!ilu->reuseordering) { 7179b54502bSHong Zhang if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 7189b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 7199b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 720efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 721efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 7229b54502bSHong Zhang } 7237529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 72452e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 7259b54502bSHong Zhang } else if (!ilu->reusefill) { 7269b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 7277529efd4SKris Buschelman ierr = MatILUDTFactor(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 72852e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 7299b54502bSHong Zhang } else { 7309b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 7319b54502bSHong Zhang } 7329b54502bSHong Zhang } else { 7339b54502bSHong Zhang if (!pc->setupcalled) { 7349b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 7359b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 736efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 737efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 7389b54502bSHong Zhang /* Remove zeros along diagonal? */ 7399b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 7409b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 7419b54502bSHong Zhang } 7429b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 74352e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 7449b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 7459b54502bSHong Zhang if (!ilu->reuseordering) { 7469b54502bSHong Zhang /* compute a new ordering for the ILU */ 7479b54502bSHong Zhang ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 7489b54502bSHong Zhang ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 7499b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 750efee365bSSatish Balay if (ilu->row) {ierr = PetscLogObjectParent(pc,ilu->row);CHKERRQ(ierr);} 751efee365bSSatish Balay if (ilu->col) {ierr = PetscLogObjectParent(pc,ilu->col);CHKERRQ(ierr);} 7529b54502bSHong Zhang /* Remove zeros along diagonal? */ 7539b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 7549b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 7559b54502bSHong Zhang } 7569b54502bSHong Zhang } 7579b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 7589b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 75952e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,ilu->fact);CHKERRQ(ierr); 7609b54502bSHong Zhang } 7619b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 7629b54502bSHong Zhang } 7639b54502bSHong Zhang PetscFunctionReturn(0); 7649b54502bSHong Zhang } 7659b54502bSHong Zhang 7669b54502bSHong Zhang #undef __FUNCT__ 7679b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 7689b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 7699b54502bSHong Zhang { 7709b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 7719b54502bSHong Zhang PetscErrorCode ierr; 7729b54502bSHong Zhang 7739b54502bSHong Zhang PetscFunctionBegin; 7749b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 7759b54502bSHong Zhang ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 7769b54502bSHong Zhang ierr = PetscFree(ilu);CHKERRQ(ierr); 7779b54502bSHong Zhang PetscFunctionReturn(0); 7789b54502bSHong Zhang } 7799b54502bSHong Zhang 7809b54502bSHong Zhang #undef __FUNCT__ 7819b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 7829b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 7839b54502bSHong Zhang { 7849b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 7859b54502bSHong Zhang PetscErrorCode ierr; 7869b54502bSHong Zhang 7879b54502bSHong Zhang PetscFunctionBegin; 7889b54502bSHong Zhang ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 7899b54502bSHong Zhang PetscFunctionReturn(0); 7909b54502bSHong Zhang } 7919b54502bSHong Zhang 7929b54502bSHong Zhang #undef __FUNCT__ 7939b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 7949b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 7959b54502bSHong Zhang { 7969b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 7979b54502bSHong Zhang PetscErrorCode ierr; 7989b54502bSHong Zhang 7999b54502bSHong Zhang PetscFunctionBegin; 8009b54502bSHong Zhang ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 8019b54502bSHong Zhang PetscFunctionReturn(0); 8029b54502bSHong Zhang } 8039b54502bSHong Zhang 8049b54502bSHong Zhang #undef __FUNCT__ 8059b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ILU" 8069b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 8079b54502bSHong Zhang { 8089b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 8099b54502bSHong Zhang 8109b54502bSHong Zhang PetscFunctionBegin; 8119b54502bSHong Zhang if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 8129b54502bSHong Zhang *mat = ilu->fact; 8139b54502bSHong Zhang PetscFunctionReturn(0); 8149b54502bSHong Zhang } 8159b54502bSHong Zhang 8169b54502bSHong Zhang /*MC 8179b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 8189b54502bSHong Zhang 8199b54502bSHong Zhang Options Database Keys: 8209b54502bSHong Zhang + -pc_ilu_levels <k> - number of levels of fill for ILU(k) 8219b54502bSHong Zhang . -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 8229b54502bSHong Zhang its factorization (overwrites original matrix) 8239b54502bSHong Zhang . -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 8249b54502bSHong Zhang . -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization 8259b54502bSHong Zhang . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 826*55ba2a51SBarry Smith . -pc_factor_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 8279b54502bSHong Zhang . -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 8289b54502bSHong Zhang this decreases the chance of getting a zero pivot 8299b54502bSHong Zhang . -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 830f251bdbdSHong Zhang . -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 8319b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 8329b54502bSHong Zhang stability of the ILU factorization 833f251bdbdSHong Zhang . -pc_factor_shift_nonzero <shift> - Sets shift amount or PETSC_DECIDE for the default 834f251bdbdSHong Zhang - -pc_factor_shift_positive_definite [PETSC_TRUE/PETSC_FALSE] - Activate/Deactivate PCFactorSetShiftPd(); the value 835f251bdbdSHong Zhang is optional with PETSC_TRUE being the default 8369b54502bSHong Zhang 8379b54502bSHong Zhang Level: beginner 8389b54502bSHong Zhang 8399b54502bSHong Zhang Concepts: incomplete factorization 8409b54502bSHong Zhang 8419b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 8429b54502bSHong Zhang 8439b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 8449b54502bSHong Zhang 8459b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 846ee45ca4aSHong Zhang PCFactorSetZeroPivot(), PCFactorSetShiftNonzero(), PCFactorSetShiftPd(), PCILUSetUseDropTolerance(), 847*55ba2a51SBarry Smith PCFactorSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(), 8489b54502bSHong Zhang PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(), 849f251bdbdSHong Zhang PCFactorSetShiftNonzero(),PCFactorSetShiftPd() 8509b54502bSHong Zhang 8519b54502bSHong Zhang M*/ 8529b54502bSHong Zhang 8539b54502bSHong Zhang EXTERN_C_BEGIN 8549b54502bSHong Zhang #undef __FUNCT__ 8559b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 856dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ILU(PC pc) 8579b54502bSHong Zhang { 8589b54502bSHong Zhang PetscErrorCode ierr; 8599b54502bSHong Zhang PC_ILU *ilu; 8609b54502bSHong Zhang 8619b54502bSHong Zhang PetscFunctionBegin; 8629b54502bSHong Zhang ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr); 86352e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_ILU));CHKERRQ(ierr); 8649b54502bSHong Zhang 8659b54502bSHong Zhang ilu->fact = 0; 8669b54502bSHong Zhang ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 8679b54502bSHong Zhang ilu->info.levels = 0; 8689b54502bSHong Zhang ilu->info.fill = 1.0; 8699b54502bSHong Zhang ilu->col = 0; 8709b54502bSHong Zhang ilu->row = 0; 8719b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 8729b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 8739b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 8749b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 8759b54502bSHong Zhang ilu->info.dt = PETSC_DEFAULT; 8769b54502bSHong Zhang ilu->info.dtcount = PETSC_DEFAULT; 8779b54502bSHong Zhang ilu->info.dtcol = PETSC_DEFAULT; 8780a29876aSHong Zhang ilu->info.shiftnz = 0.0; 879fbf22428SSatish Balay ilu->info.shiftpd = 0.0; /* false */ 8809b54502bSHong Zhang ilu->info.shift_fraction = 0.0; 8819b54502bSHong Zhang ilu->info.zeropivot = 1.e-12; 8829b54502bSHong Zhang ilu->info.pivotinblocks = 1.0; 8839b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 8849b54502bSHong Zhang ilu->info.diagonal_fill = 0; 8859b54502bSHong Zhang pc->data = (void*)ilu; 8869b54502bSHong Zhang 8879b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 8889b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 8899b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 8909b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 8919b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 8929b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 8939b54502bSHong Zhang pc->ops->view = PCView_ILU; 8949b54502bSHong Zhang pc->ops->applyrichardson = 0; 8959b54502bSHong Zhang 896afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetZeroPivot_C","PCFactorSetZeroPivot_ILU", 897afaefe49SHong Zhang PCFactorSetZeroPivot_ILU);CHKERRQ(ierr); 898afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftNonzero_C","PCFactorSetShiftNonzero_ILU", 899afaefe49SHong Zhang PCFactorSetShiftNonzero_ILU);CHKERRQ(ierr); 900afaefe49SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetShiftPd_C","PCFactorSetShiftPd_ILU", 901afaefe49SHong Zhang PCFactorSetShiftPd_ILU);CHKERRQ(ierr); 902afaefe49SHong Zhang 9039b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU", 9049b54502bSHong Zhang PCILUSetUseDropTolerance_ILU);CHKERRQ(ierr); 905*55ba2a51SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCFactorSetFill_C","PCFactorSetFill_ILU", 906*55ba2a51SBarry Smith PCFactorSetFill_ILU);CHKERRQ(ierr); 9079b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU", 9089b54502bSHong Zhang PCILUSetMatOrdering_ILU);CHKERRQ(ierr); 9099b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU", 9109b54502bSHong Zhang PCILUSetReuseOrdering_ILU);CHKERRQ(ierr); 9119b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT", 9129b54502bSHong Zhang PCILUDTSetReuseFill_ILUDT);CHKERRQ(ierr); 9139b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU", 9149b54502bSHong Zhang PCILUSetLevels_ILU);CHKERRQ(ierr); 9159b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU", 9169b54502bSHong Zhang PCILUSetUseInPlace_ILU);CHKERRQ(ierr); 9179b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU", 9189b54502bSHong Zhang PCILUSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 9199b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU", 9209b54502bSHong Zhang PCILUSetPivotInBlocks_ILU);CHKERRQ(ierr); 9219b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C","PCILUReorderForNonzeroDiagonal_ILU", 9229b54502bSHong Zhang PCILUReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 9239b54502bSHong Zhang PetscFunctionReturn(0); 9249b54502bSHong Zhang } 9259b54502bSHong Zhang EXTERN_C_END 926