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