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 /* ------------------------------------------------------------------------------------------*/ 99b54502bSHong Zhang 109b54502bSHong Zhang EXTERN_C_BEGIN 119b54502bSHong Zhang #undef __FUNCT__ 129b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal_ILU" 139b54502bSHong Zhang PetscErrorCode PCILUReorderForNonzeroDiagonal_ILU(PC pc,PetscReal z) 149b54502bSHong Zhang { 159b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 169b54502bSHong Zhang 179b54502bSHong Zhang PetscFunctionBegin; 189b54502bSHong Zhang ilu->nonzerosalongdiagonal = PETSC_TRUE; 199b54502bSHong Zhang if (z == PETSC_DECIDE) { 209b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = 1.e-10; 219b54502bSHong Zhang } else { 229b54502bSHong Zhang ilu->nonzerosalongdiagonaltol = z; 239b54502bSHong Zhang } 249b54502bSHong Zhang PetscFunctionReturn(0); 259b54502bSHong Zhang } 269b54502bSHong Zhang EXTERN_C_END 279b54502bSHong Zhang 289b54502bSHong Zhang EXTERN_C_BEGIN 299b54502bSHong Zhang #undef __FUNCT__ 309b54502bSHong Zhang #define __FUNCT__ "PCILUSetSetZeroPivot_ILU" 319b54502bSHong Zhang PetscErrorCode PCILUSetZeroPivot_ILU(PC pc,PetscReal z) 329b54502bSHong Zhang { 339b54502bSHong Zhang PC_ILU *lu; 349b54502bSHong Zhang 359b54502bSHong Zhang PetscFunctionBegin; 369b54502bSHong Zhang lu = (PC_ILU*)pc->data; 379b54502bSHong Zhang lu->info.zeropivot = z; 389b54502bSHong Zhang PetscFunctionReturn(0); 399b54502bSHong Zhang } 409b54502bSHong Zhang EXTERN_C_END 419b54502bSHong Zhang 429b54502bSHong Zhang #undef __FUNCT__ 439b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU_Internal" 449b54502bSHong Zhang PetscErrorCode PCDestroy_ILU_Internal(PC pc) 459b54502bSHong Zhang { 469b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 479b54502bSHong Zhang PetscErrorCode ierr; 489b54502bSHong Zhang 499b54502bSHong Zhang PetscFunctionBegin; 509b54502bSHong Zhang if (!ilu->inplace && ilu->fact) {ierr = MatDestroy(ilu->fact);CHKERRQ(ierr);} 519b54502bSHong Zhang if (ilu->row && ilu->col && ilu->row != ilu->col) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 529b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 539b54502bSHong Zhang PetscFunctionReturn(0); 549b54502bSHong Zhang } 559b54502bSHong Zhang 569b54502bSHong Zhang EXTERN_C_BEGIN 579b54502bSHong Zhang #undef __FUNCT__ 589b54502bSHong Zhang #define __FUNCT__ "PCILUSetDamping_ILU" 599b54502bSHong Zhang PetscErrorCode PCILUSetDamping_ILU(PC pc,PetscReal damping) 609b54502bSHong Zhang { 619b54502bSHong Zhang PC_ILU *dir; 629b54502bSHong Zhang 639b54502bSHong Zhang PetscFunctionBegin; 649b54502bSHong Zhang dir = (PC_ILU*)pc->data; 659b54502bSHong Zhang if (damping == (PetscReal) PETSC_DECIDE) { 669b54502bSHong Zhang dir->info.damping = 1.e-12; 679b54502bSHong Zhang } else { 689b54502bSHong Zhang dir->info.damping = damping; 699b54502bSHong Zhang } 709b54502bSHong Zhang PetscFunctionReturn(0); 719b54502bSHong Zhang } 729b54502bSHong Zhang EXTERN_C_END 739b54502bSHong Zhang 749b54502bSHong Zhang EXTERN_C_BEGIN 759b54502bSHong Zhang #undef __FUNCT__ 769b54502bSHong Zhang #define __FUNCT__ "PCILUSetShift_ILU" 779b54502bSHong Zhang PetscErrorCode PCILUSetShift_ILU(PC pc,PetscTruth shift) 789b54502bSHong Zhang { 799b54502bSHong Zhang PC_ILU *dir; 809b54502bSHong Zhang 819b54502bSHong Zhang PetscFunctionBegin; 829b54502bSHong Zhang dir = (PC_ILU*)pc->data; 839b54502bSHong Zhang dir->info.shift = shift; 849b54502bSHong Zhang if (shift) dir->info.shift_fraction = 0.0; 859b54502bSHong Zhang PetscFunctionReturn(0); 869b54502bSHong Zhang } 879b54502bSHong Zhang EXTERN_C_END 889b54502bSHong Zhang 899b54502bSHong Zhang EXTERN_C_BEGIN 909b54502bSHong Zhang #undef __FUNCT__ 919b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance_ILU" 929b54502bSHong Zhang PetscErrorCode PCILUSetUseDropTolerance_ILU(PC pc,PetscReal dt,PetscReal dtcol,PetscInt dtcount) 939b54502bSHong Zhang { 949b54502bSHong Zhang PC_ILU *ilu; 959b54502bSHong Zhang PetscErrorCode ierr; 969b54502bSHong Zhang 979b54502bSHong Zhang PetscFunctionBegin; 989b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 999b54502bSHong Zhang if (pc->setupcalled && (!ilu->usedt || ilu->info.dt != dt || ilu->info.dtcol != dtcol || ilu->info.dtcount != dtcount)) { 1009b54502bSHong Zhang pc->setupcalled = 0; 1019b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1029b54502bSHong Zhang } 1039b54502bSHong Zhang ilu->usedt = PETSC_TRUE; 1049b54502bSHong Zhang ilu->info.dt = dt; 1059b54502bSHong Zhang ilu->info.dtcol = dtcol; 1069b54502bSHong Zhang ilu->info.dtcount = dtcount; 1079b54502bSHong Zhang ilu->info.fill = PETSC_DEFAULT; 1089b54502bSHong Zhang PetscFunctionReturn(0); 1099b54502bSHong Zhang } 1109b54502bSHong Zhang EXTERN_C_END 1119b54502bSHong Zhang 1129b54502bSHong Zhang EXTERN_C_BEGIN 1139b54502bSHong Zhang #undef __FUNCT__ 1149b54502bSHong Zhang #define __FUNCT__ "PCILUSetFill_ILU" 1159b54502bSHong Zhang PetscErrorCode PCILUSetFill_ILU(PC pc,PetscReal fill) 1169b54502bSHong Zhang { 1179b54502bSHong Zhang PC_ILU *dir; 1189b54502bSHong Zhang 1199b54502bSHong Zhang PetscFunctionBegin; 1209b54502bSHong Zhang dir = (PC_ILU*)pc->data; 1219b54502bSHong Zhang dir->info.fill = fill; 1229b54502bSHong Zhang PetscFunctionReturn(0); 1239b54502bSHong Zhang } 1249b54502bSHong Zhang EXTERN_C_END 1259b54502bSHong Zhang 1269b54502bSHong Zhang EXTERN_C_BEGIN 1279b54502bSHong Zhang #undef __FUNCT__ 1289b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering_ILU" 1299b54502bSHong Zhang PetscErrorCode PCILUSetMatOrdering_ILU(PC pc,MatOrderingType ordering) 1309b54502bSHong Zhang { 1319b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 1329b54502bSHong Zhang PetscErrorCode ierr; 1339b54502bSHong Zhang PetscTruth flg; 1349b54502bSHong Zhang 1359b54502bSHong Zhang PetscFunctionBegin; 1369b54502bSHong Zhang if (!pc->setupcalled) { 1379b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1389b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1399b54502bSHong Zhang } else { 1409b54502bSHong Zhang ierr = PetscStrcmp(dir->ordering,ordering,&flg);CHKERRQ(ierr); 1419b54502bSHong Zhang if (!flg) { 1429b54502bSHong Zhang pc->setupcalled = 0; 1439b54502bSHong Zhang ierr = PetscStrfree(dir->ordering);CHKERRQ(ierr); 1449b54502bSHong Zhang ierr = PetscStrallocpy(ordering,&dir->ordering);CHKERRQ(ierr); 1459b54502bSHong Zhang /* free the data structures, then create them again */ 1469b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 1479b54502bSHong Zhang } 1489b54502bSHong Zhang } 1499b54502bSHong Zhang PetscFunctionReturn(0); 1509b54502bSHong Zhang } 1519b54502bSHong Zhang EXTERN_C_END 1529b54502bSHong Zhang 1539b54502bSHong Zhang EXTERN_C_BEGIN 1549b54502bSHong Zhang #undef __FUNCT__ 1559b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering_ILU" 1569b54502bSHong Zhang PetscErrorCode PCILUSetReuseOrdering_ILU(PC pc,PetscTruth flag) 1579b54502bSHong Zhang { 1589b54502bSHong Zhang PC_ILU *ilu; 1599b54502bSHong Zhang 1609b54502bSHong Zhang PetscFunctionBegin; 1619b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1629b54502bSHong Zhang ilu->reuseordering = flag; 1639b54502bSHong Zhang PetscFunctionReturn(0); 1649b54502bSHong Zhang } 1659b54502bSHong Zhang EXTERN_C_END 1669b54502bSHong Zhang 1679b54502bSHong Zhang EXTERN_C_BEGIN 1689b54502bSHong Zhang #undef __FUNCT__ 1699b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill_ILUDT" 1709b54502bSHong Zhang PetscErrorCode PCILUDTSetReuseFill_ILUDT(PC pc,PetscTruth flag) 1719b54502bSHong Zhang { 1729b54502bSHong Zhang PC_ILU *ilu; 1739b54502bSHong Zhang 1749b54502bSHong Zhang PetscFunctionBegin; 1759b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1769b54502bSHong Zhang ilu->reusefill = flag; 1779b54502bSHong Zhang if (flag) SETERRQ(PETSC_ERR_SUP,"Not yet supported"); 1789b54502bSHong Zhang PetscFunctionReturn(0); 1799b54502bSHong Zhang } 1809b54502bSHong Zhang EXTERN_C_END 1819b54502bSHong Zhang 1829b54502bSHong Zhang EXTERN_C_BEGIN 1839b54502bSHong Zhang #undef __FUNCT__ 1849b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels_ILU" 1859b54502bSHong Zhang PetscErrorCode PCILUSetLevels_ILU(PC pc,PetscInt levels) 1869b54502bSHong Zhang { 1879b54502bSHong Zhang PC_ILU *ilu; 1889b54502bSHong Zhang PetscErrorCode ierr; 1899b54502bSHong Zhang 1909b54502bSHong Zhang PetscFunctionBegin; 1919b54502bSHong Zhang ilu = (PC_ILU*)pc->data; 1929b54502bSHong Zhang 1939b54502bSHong Zhang if (!pc->setupcalled) { 1949b54502bSHong Zhang ilu->info.levels = levels; 1959b54502bSHong Zhang } else if (ilu->usedt || ilu->info.levels != levels) { 1969b54502bSHong Zhang ilu->info.levels = levels; 1979b54502bSHong Zhang pc->setupcalled = 0; 1989b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 1999b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 2009b54502bSHong Zhang } 2019b54502bSHong Zhang PetscFunctionReturn(0); 2029b54502bSHong Zhang } 2039b54502bSHong Zhang EXTERN_C_END 2049b54502bSHong Zhang 2059b54502bSHong Zhang EXTERN_C_BEGIN 2069b54502bSHong Zhang #undef __FUNCT__ 2079b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace_ILU" 2089b54502bSHong Zhang PetscErrorCode PCILUSetUseInPlace_ILU(PC pc) 2099b54502bSHong Zhang { 2109b54502bSHong Zhang PC_ILU *dir; 2119b54502bSHong Zhang 2129b54502bSHong Zhang PetscFunctionBegin; 2139b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2149b54502bSHong Zhang dir->inplace = PETSC_TRUE; 2159b54502bSHong Zhang PetscFunctionReturn(0); 2169b54502bSHong Zhang } 2179b54502bSHong Zhang EXTERN_C_END 2189b54502bSHong Zhang 2199b54502bSHong Zhang EXTERN_C_BEGIN 2209b54502bSHong Zhang #undef __FUNCT__ 2219b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill" 2229b54502bSHong Zhang PetscErrorCode PCILUSetAllowDiagonalFill_ILU(PC pc) 2239b54502bSHong Zhang { 2249b54502bSHong Zhang PC_ILU *dir; 2259b54502bSHong Zhang 2269b54502bSHong Zhang PetscFunctionBegin; 2279b54502bSHong Zhang dir = (PC_ILU*)pc->data; 2289b54502bSHong Zhang dir->info.diagonal_fill = 1; 2299b54502bSHong Zhang PetscFunctionReturn(0); 2309b54502bSHong Zhang } 2319b54502bSHong Zhang EXTERN_C_END 2329b54502bSHong Zhang 2339b54502bSHong Zhang #undef __FUNCT__ 2349b54502bSHong Zhang #define __FUNCT__ "PCILUSetZeroPivot" 2359b54502bSHong Zhang /*@ 2369b54502bSHong Zhang PCILUSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 2379b54502bSHong Zhang 2389b54502bSHong Zhang Collective on PC 2399b54502bSHong Zhang 2409b54502bSHong Zhang Input Parameters: 2419b54502bSHong Zhang + pc - the preconditioner context 2429b54502bSHong Zhang - zero - all pivots smaller than this will be considered zero 2439b54502bSHong Zhang 2449b54502bSHong Zhang Options Database Key: 2459b54502bSHong Zhang . -pc_ilu_zeropivot <zero> - Sets the zero pivot size 2469b54502bSHong Zhang 2479b54502bSHong Zhang Level: intermediate 2489b54502bSHong Zhang 2499b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 2509b54502bSHong Zhang 2519b54502bSHong Zhang .seealso: PCILUSetFill(), PCLUSetDamp(), PCLUSetZeroPivot() 2529b54502bSHong Zhang @*/ 2539b54502bSHong Zhang PetscErrorCode PCILUSetZeroPivot(PC pc,PetscReal zero) 2549b54502bSHong Zhang { 2559b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 2569b54502bSHong Zhang 2579b54502bSHong Zhang PetscFunctionBegin; 2589b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2599b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetZeroPivot_C",(void (**)(void))&f);CHKERRQ(ierr); 2609b54502bSHong Zhang if (f) { 2619b54502bSHong Zhang ierr = (*f)(pc,zero);CHKERRQ(ierr); 2629b54502bSHong Zhang } 2639b54502bSHong Zhang PetscFunctionReturn(0); 2649b54502bSHong Zhang } 2659b54502bSHong Zhang 2669b54502bSHong Zhang #undef __FUNCT__ 2679b54502bSHong Zhang #define __FUNCT__ "PCILUSetDamping" 2689b54502bSHong Zhang /*@ 2699b54502bSHong Zhang PCILUSetDamping - adds this quantity to the diagonal of the matrix during the 2709b54502bSHong Zhang ILU numerical factorization 2719b54502bSHong Zhang 2729b54502bSHong Zhang Collective on PC 2739b54502bSHong Zhang 2749b54502bSHong Zhang Input Parameters: 2759b54502bSHong Zhang + pc - the preconditioner context 2769b54502bSHong Zhang - damping - amount of damping 2779b54502bSHong Zhang 2789b54502bSHong Zhang Options Database Key: 2799b54502bSHong Zhang . -pc_ilu_damping <damping> - Sets damping amount or PETSC_DECIDE for the default 2809b54502bSHong Zhang 2819b54502bSHong Zhang Note: If 0.0 is given, then no damping is used. If a diagonal element is classified as a zero 2829b54502bSHong Zhang pivot, then the damping is doubled until this is alleviated. 2839b54502bSHong Zhang 2849b54502bSHong Zhang Level: intermediate 2859b54502bSHong Zhang 2869b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 2879b54502bSHong Zhang 2889b54502bSHong Zhang .seealso: PCILUSetFill(), PCLUSetDamping(), PCILUSetShift() 2899b54502bSHong Zhang @*/ 2909b54502bSHong Zhang PetscErrorCode PCILUSetDamping(PC pc,PetscReal damping) 2919b54502bSHong Zhang { 2929b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 2939b54502bSHong Zhang 2949b54502bSHong Zhang PetscFunctionBegin; 2959b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2969b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetDamping_C",(void (**)(void))&f);CHKERRQ(ierr); 2979b54502bSHong Zhang if (f) { 2989b54502bSHong Zhang ierr = (*f)(pc,damping);CHKERRQ(ierr); 2999b54502bSHong Zhang } 3009b54502bSHong Zhang PetscFunctionReturn(0); 3019b54502bSHong Zhang } 3029b54502bSHong Zhang 3039b54502bSHong Zhang #undef __FUNCT__ 3049b54502bSHong Zhang #define __FUNCT__ "PCILUSetShift" 3059b54502bSHong Zhang /*@ 3069b54502bSHong Zhang PCILUSetShift - specify whether to use Manteuffel shifting of ILU. 3079b54502bSHong Zhang If an ILU factorisation breaks down because of nonpositive pivots, 3089b54502bSHong Zhang adding sufficient identity to the diagonal will remedy this. 3099b54502bSHong Zhang Setting this causes a bisection method to find the minimum shift that 3109b54502bSHong Zhang will lead to a well-defined ILU. 3119b54502bSHong Zhang 3129b54502bSHong Zhang Input parameters: 3139b54502bSHong Zhang + pc - the preconditioner context 3149b54502bSHong Zhang - shifting - PETSC_TRUE to set shift else PETSC_FALSE 3159b54502bSHong Zhang 3169b54502bSHong Zhang Options Database Key: 3179b54502bSHong Zhang . -pc_ilu_shift [1/0] - Activate/Deactivate PCILUSetShift(); the value 3189b54502bSHong Zhang is optional with 1 being the default 3199b54502bSHong Zhang 3209b54502bSHong Zhang Level: intermediate 3219b54502bSHong Zhang 3229b54502bSHong Zhang .keywords: PC, indefinite, factorization, incomplete, ILU 3239b54502bSHong Zhang 3249b54502bSHong Zhang .seealso: PCILUSetDamping() 3259b54502bSHong Zhang @*/ 3269b54502bSHong Zhang PetscErrorCode PCILUSetShift(PC pc,PetscTruth shifting) 3279b54502bSHong Zhang { 3289b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 3299b54502bSHong Zhang 3309b54502bSHong Zhang PetscFunctionBegin; 3319b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3329b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetShift_C",(void (**)(void))&f);CHKERRQ(ierr); 3339b54502bSHong Zhang if (f) { 3349b54502bSHong Zhang ierr = (*f)(pc,shifting);CHKERRQ(ierr); 3359b54502bSHong Zhang } 3369b54502bSHong Zhang PetscFunctionReturn(0); 3379b54502bSHong Zhang } 3389b54502bSHong Zhang 3399b54502bSHong Zhang 3409b54502bSHong Zhang #undef __FUNCT__ 3419b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseDropTolerance" 3429b54502bSHong Zhang /*@ 3439b54502bSHong Zhang PCILUSetUseDropTolerance - The preconditioner will use an ILU 3449b54502bSHong Zhang based on a drop tolerance. 3459b54502bSHong Zhang 3469b54502bSHong Zhang Collective on PC 3479b54502bSHong Zhang 3489b54502bSHong Zhang Input Parameters: 3499b54502bSHong Zhang + pc - the preconditioner context 3509b54502bSHong Zhang . dt - the drop tolerance, try from 1.e-10 to .1 3519b54502bSHong Zhang . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 3529b54502bSHong Zhang - maxrowcount - the max number of nonzeros allowed in a row, best value 3539b54502bSHong Zhang depends on the number of nonzeros in row of original matrix 3549b54502bSHong Zhang 3559b54502bSHong Zhang Options Database Key: 3569b54502bSHong Zhang . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 3579b54502bSHong Zhang 3589b54502bSHong Zhang Level: intermediate 3599b54502bSHong Zhang 3609b54502bSHong Zhang Notes: 3619b54502bSHong Zhang This uses the iludt() code of Saad's SPARSKIT package 3629b54502bSHong Zhang 3639b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 3649b54502bSHong Zhang @*/ 3659b54502bSHong Zhang PetscErrorCode PCILUSetUseDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 3669b54502bSHong Zhang { 3679b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal,PetscReal,PetscInt); 3689b54502bSHong Zhang 3699b54502bSHong Zhang PetscFunctionBegin; 3709b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3719b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseDropTolerance_C",(void (**)(void))&f);CHKERRQ(ierr); 3729b54502bSHong Zhang if (f) { 3739b54502bSHong Zhang ierr = (*f)(pc,dt,dtcol,maxrowcount);CHKERRQ(ierr); 3749b54502bSHong Zhang } 3759b54502bSHong Zhang PetscFunctionReturn(0); 3769b54502bSHong Zhang } 3779b54502bSHong Zhang 3789b54502bSHong Zhang #undef __FUNCT__ 3799b54502bSHong Zhang #define __FUNCT__ "PCILUSetFill" 3809b54502bSHong Zhang /*@ 3819b54502bSHong Zhang PCILUSetFill - Indicate the amount of fill you expect in the factored matrix, 3829b54502bSHong Zhang where fill = number nonzeros in factor/number nonzeros in original matrix. 3839b54502bSHong Zhang 3849b54502bSHong Zhang Collective on PC 3859b54502bSHong Zhang 3869b54502bSHong Zhang Input Parameters: 3879b54502bSHong Zhang + pc - the preconditioner context 3889b54502bSHong Zhang - fill - amount of expected fill 3899b54502bSHong Zhang 3909b54502bSHong Zhang Options Database Key: 3919b54502bSHong Zhang $ -pc_ilu_fill <fill> 3929b54502bSHong Zhang 3939b54502bSHong Zhang Note: 3949b54502bSHong Zhang For sparse matrix factorizations it is difficult to predict how much 3959b54502bSHong Zhang fill to expect. By running with the option -log_info PETSc will print the 3969b54502bSHong Zhang actual amount of fill used; allowing you to set the value accurately for 3979b54502bSHong Zhang future runs. But default PETSc uses a value of 1.0 3989b54502bSHong Zhang 3999b54502bSHong Zhang Level: intermediate 4009b54502bSHong Zhang 4019b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 4029b54502bSHong Zhang 4039b54502bSHong Zhang .seealso: PCLUSetFill() 4049b54502bSHong Zhang @*/ 4059b54502bSHong Zhang PetscErrorCode PCILUSetFill(PC pc,PetscReal fill) 4069b54502bSHong Zhang { 4079b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 4089b54502bSHong Zhang 4099b54502bSHong Zhang PetscFunctionBegin; 4109b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4119b54502bSHong Zhang if (fill < 1.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less than 1.0"); 4129b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetFill_C",(void (**)(void))&f);CHKERRQ(ierr); 4139b54502bSHong Zhang if (f) { 4149b54502bSHong Zhang ierr = (*f)(pc,fill);CHKERRQ(ierr); 4159b54502bSHong Zhang } 4169b54502bSHong Zhang PetscFunctionReturn(0); 4179b54502bSHong Zhang } 4189b54502bSHong Zhang 4199b54502bSHong Zhang #undef __FUNCT__ 4209b54502bSHong Zhang #define __FUNCT__ "PCILUSetMatOrdering" 4219b54502bSHong Zhang /*@C 4229b54502bSHong Zhang PCILUSetMatOrdering - Sets the ordering routine (to reduce fill) to 4239b54502bSHong Zhang be used in the ILU factorization. 4249b54502bSHong Zhang 4259b54502bSHong Zhang Collective on PC 4269b54502bSHong Zhang 4279b54502bSHong Zhang Input Parameters: 4289b54502bSHong Zhang + pc - the preconditioner context 4299b54502bSHong Zhang - ordering - the matrix ordering name, for example, MATORDERING_ND or MATORDERING_RCM 4309b54502bSHong Zhang 4319b54502bSHong Zhang Options Database Key: 4329b54502bSHong Zhang . -pc_ilu_mat_ordering_type <nd,rcm,...> - Sets ordering routine 4339b54502bSHong Zhang 4349b54502bSHong Zhang Level: intermediate 4359b54502bSHong Zhang 4369b54502bSHong Zhang Notes: natural ordering is used by default 4379b54502bSHong Zhang 4389b54502bSHong Zhang .seealso: PCLUSetMatOrdering() 4399b54502bSHong Zhang 4409b54502bSHong Zhang .keywords: PC, ILU, set, matrix, reordering 4419b54502bSHong Zhang 4429b54502bSHong Zhang @*/ 4439b54502bSHong Zhang PetscErrorCode PCILUSetMatOrdering(PC pc,MatOrderingType ordering) 4449b54502bSHong Zhang { 4459b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,MatOrderingType); 4469b54502bSHong Zhang 4479b54502bSHong Zhang PetscFunctionBegin; 4489b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4499b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetMatOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 4509b54502bSHong Zhang if (f) { 4519b54502bSHong Zhang ierr = (*f)(pc,ordering);CHKERRQ(ierr); 4529b54502bSHong Zhang } 4539b54502bSHong Zhang PetscFunctionReturn(0); 4549b54502bSHong Zhang } 4559b54502bSHong Zhang 4569b54502bSHong Zhang #undef __FUNCT__ 4579b54502bSHong Zhang #define __FUNCT__ "PCILUSetReuseOrdering" 4589b54502bSHong Zhang /*@ 4599b54502bSHong Zhang PCILUSetReuseOrdering - When similar matrices are factored, this 4609b54502bSHong Zhang causes the ordering computed in the first factor to be used for all 4619b54502bSHong Zhang following factors; applies to both fill and drop tolerance ILUs. 4629b54502bSHong Zhang 4639b54502bSHong Zhang Collective on PC 4649b54502bSHong Zhang 4659b54502bSHong Zhang Input Parameters: 4669b54502bSHong Zhang + pc - the preconditioner context 4679b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 4689b54502bSHong Zhang 4699b54502bSHong Zhang Options Database Key: 4709b54502bSHong Zhang . -pc_ilu_reuse_ordering - Activate PCILUSetReuseOrdering() 4719b54502bSHong Zhang 4729b54502bSHong Zhang Level: intermediate 4739b54502bSHong Zhang 4749b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 4759b54502bSHong Zhang 4769b54502bSHong Zhang .seealso: PCILUDTSetReuseFill(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 4779b54502bSHong Zhang @*/ 4789b54502bSHong Zhang PetscErrorCode PCILUSetReuseOrdering(PC pc,PetscTruth flag) 4799b54502bSHong Zhang { 4809b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 4819b54502bSHong Zhang 4829b54502bSHong Zhang PetscFunctionBegin; 4839b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4849b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetReuseOrdering_C",(void (**)(void))&f);CHKERRQ(ierr); 4859b54502bSHong Zhang if (f) { 4869b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 4879b54502bSHong Zhang } 4889b54502bSHong Zhang PetscFunctionReturn(0); 4899b54502bSHong Zhang } 4909b54502bSHong Zhang 4919b54502bSHong Zhang #undef __FUNCT__ 4929b54502bSHong Zhang #define __FUNCT__ "PCILUDTSetReuseFill" 4939b54502bSHong Zhang /*@ 4949b54502bSHong Zhang PCILUDTSetReuseFill - When matrices with same nonzero structure are ILUDT factored, 4959b54502bSHong Zhang this causes later ones to use the fill computed in the initial factorization. 4969b54502bSHong Zhang 4979b54502bSHong Zhang Collective on PC 4989b54502bSHong Zhang 4999b54502bSHong Zhang Input Parameters: 5009b54502bSHong Zhang + pc - the preconditioner context 5019b54502bSHong Zhang - flag - PETSC_TRUE to reuse else PETSC_FALSE 5029b54502bSHong Zhang 5039b54502bSHong Zhang Options Database Key: 5049b54502bSHong Zhang . -pc_iludt_reuse_fill - Activates PCILUDTSetReuseFill() 5059b54502bSHong Zhang 5069b54502bSHong Zhang Level: intermediate 5079b54502bSHong Zhang 5089b54502bSHong Zhang .keywords: PC, levels, reordering, factorization, incomplete, ILU 5099b54502bSHong Zhang 5109b54502bSHong Zhang .seealso: PCILUSetReuseOrdering(), PCLUSetReuseOrdering(), PCLUSetReuseFill() 5119b54502bSHong Zhang @*/ 5129b54502bSHong Zhang PetscErrorCode PCILUDTSetReuseFill(PC pc,PetscTruth flag) 5139b54502bSHong Zhang { 5149b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 5159b54502bSHong Zhang 5169b54502bSHong Zhang PetscFunctionBegin; 5179b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5189b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUDTSetReuseFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5199b54502bSHong Zhang if (f) { 5209b54502bSHong Zhang ierr = (*f)(pc,flag);CHKERRQ(ierr); 5219b54502bSHong Zhang } 5229b54502bSHong Zhang PetscFunctionReturn(0); 5239b54502bSHong Zhang } 5249b54502bSHong Zhang 5259b54502bSHong Zhang #undef __FUNCT__ 5269b54502bSHong Zhang #define __FUNCT__ "PCILUSetLevels" 5279b54502bSHong Zhang /*@ 5289b54502bSHong Zhang PCILUSetLevels - Sets the number of levels of fill to use. 5299b54502bSHong Zhang 5309b54502bSHong Zhang Collective on PC 5319b54502bSHong Zhang 5329b54502bSHong Zhang Input Parameters: 5339b54502bSHong Zhang + pc - the preconditioner context 5349b54502bSHong Zhang - levels - number of levels of fill 5359b54502bSHong Zhang 5369b54502bSHong Zhang Options Database Key: 5379b54502bSHong Zhang . -pc_ilu_levels <levels> - Sets fill level 5389b54502bSHong Zhang 5399b54502bSHong Zhang Level: intermediate 5409b54502bSHong Zhang 5419b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 5429b54502bSHong Zhang @*/ 5439b54502bSHong Zhang PetscErrorCode PCILUSetLevels(PC pc,PetscInt levels) 5449b54502bSHong Zhang { 5459b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 5469b54502bSHong Zhang 5479b54502bSHong Zhang PetscFunctionBegin; 5489b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5499b54502bSHong Zhang if (levels < 0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 5509b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetLevels_C",(void (**)(void))&f);CHKERRQ(ierr); 5519b54502bSHong Zhang if (f) { 5529b54502bSHong Zhang ierr = (*f)(pc,levels);CHKERRQ(ierr); 5539b54502bSHong Zhang } 5549b54502bSHong Zhang PetscFunctionReturn(0); 5559b54502bSHong Zhang } 5569b54502bSHong Zhang 5579b54502bSHong Zhang #undef __FUNCT__ 5589b54502bSHong Zhang #define __FUNCT__ "PCILUSetAllowDiagonalFill" 5599b54502bSHong Zhang /*@ 5609b54502bSHong Zhang PCILUSetAllowDiagonalFill - Causes all diagonal matrix entries to be 5619b54502bSHong Zhang treated as level 0 fill even if there is no non-zero location. 5629b54502bSHong Zhang 5639b54502bSHong Zhang Collective on PC 5649b54502bSHong Zhang 5659b54502bSHong Zhang Input Parameters: 5669b54502bSHong Zhang + pc - the preconditioner context 5679b54502bSHong Zhang 5689b54502bSHong Zhang Options Database Key: 5699b54502bSHong Zhang . -pc_ilu_diagonal_fill 5709b54502bSHong Zhang 5719b54502bSHong Zhang Notes: 5729b54502bSHong Zhang Does not apply with 0 fill. 5739b54502bSHong Zhang 5749b54502bSHong Zhang Level: intermediate 5759b54502bSHong Zhang 5769b54502bSHong Zhang .keywords: PC, levels, fill, factorization, incomplete, ILU 5779b54502bSHong Zhang @*/ 5789b54502bSHong Zhang PetscErrorCode PCILUSetAllowDiagonalFill(PC pc) 5799b54502bSHong Zhang { 5809b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 5819b54502bSHong Zhang 5829b54502bSHong Zhang PetscFunctionBegin; 5839b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 5849b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetAllowDiagonalFill_C",(void (**)(void))&f);CHKERRQ(ierr); 5859b54502bSHong Zhang if (f) { 5869b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 5879b54502bSHong Zhang } 5889b54502bSHong Zhang PetscFunctionReturn(0); 5899b54502bSHong Zhang } 5909b54502bSHong Zhang 5919b54502bSHong Zhang #undef __FUNCT__ 5929b54502bSHong Zhang #define __FUNCT__ "PCILUSetUseInPlace" 5939b54502bSHong Zhang /*@ 5949b54502bSHong Zhang PCILUSetUseInPlace - Tells the system to do an in-place incomplete factorization. 5959b54502bSHong Zhang Collective on PC 5969b54502bSHong Zhang 5979b54502bSHong Zhang Input Parameters: 5989b54502bSHong Zhang . pc - the preconditioner context 5999b54502bSHong Zhang 6009b54502bSHong Zhang Options Database Key: 6019b54502bSHong Zhang . -pc_ilu_in_place - Activates in-place factorization 6029b54502bSHong Zhang 6039b54502bSHong Zhang Notes: 6049b54502bSHong Zhang PCILUSetUseInPlace() is intended for use with matrix-free variants of 6059b54502bSHong Zhang Krylov methods, or when a different matrices are employed for the linear 6069b54502bSHong Zhang system and preconditioner, or with ASM preconditioning. Do NOT use 6079b54502bSHong Zhang this option if the linear system 6089b54502bSHong Zhang matrix also serves as the preconditioning matrix, since the factored 6099b54502bSHong Zhang matrix would then overwrite the original matrix. 6109b54502bSHong Zhang 6119b54502bSHong Zhang Only works well with ILU(0). 6129b54502bSHong Zhang 6139b54502bSHong Zhang Level: intermediate 6149b54502bSHong Zhang 6159b54502bSHong Zhang .keywords: PC, set, factorization, inplace, in-place, ILU 6169b54502bSHong Zhang 6179b54502bSHong Zhang .seealso: PCLUSetUseInPlace() 6189b54502bSHong Zhang @*/ 6199b54502bSHong Zhang PetscErrorCode PCILUSetUseInPlace(PC pc) 6209b54502bSHong Zhang { 6219b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC); 6229b54502bSHong Zhang 6239b54502bSHong Zhang PetscFunctionBegin; 6249b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6259b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetUseInPlace_C",(void (**)(void))&f);CHKERRQ(ierr); 6269b54502bSHong Zhang if (f) { 6279b54502bSHong Zhang ierr = (*f)(pc);CHKERRQ(ierr); 6289b54502bSHong Zhang } 6299b54502bSHong Zhang PetscFunctionReturn(0); 6309b54502bSHong Zhang } 6319b54502bSHong Zhang 6329b54502bSHong Zhang #undef __FUNCT__ 6339b54502bSHong Zhang #define __FUNCT__ "PCILUReorderForNonzeroDiagonal" 6349b54502bSHong Zhang /*@ 6359b54502bSHong Zhang PCILUReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 6369b54502bSHong Zhang 6379b54502bSHong Zhang Collective on PC 6389b54502bSHong Zhang 6399b54502bSHong Zhang Input Parameters: 6409b54502bSHong Zhang + pc - the preconditioner context 6419b54502bSHong Zhang - tol - diagonal entries smaller than this in absolute value are considered zero 6429b54502bSHong Zhang 6439b54502bSHong Zhang Options Database Key: 6449b54502bSHong Zhang . -pc_lu_nonzeros_along_diagonal 6459b54502bSHong Zhang 6469b54502bSHong Zhang Level: intermediate 6479b54502bSHong Zhang 6489b54502bSHong Zhang .keywords: PC, set, factorization, direct, fill 6499b54502bSHong Zhang 6509b54502bSHong Zhang .seealso: PCILUSetFill(), PCILUSetDamp(), PCIILUSetZeroPivot(), MatReorderForNonzeroDiagonal() 6519b54502bSHong Zhang @*/ 6529b54502bSHong Zhang PetscErrorCode PCILUReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 6539b54502bSHong Zhang { 6549b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscReal); 6559b54502bSHong Zhang 6569b54502bSHong Zhang PetscFunctionBegin; 6579b54502bSHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6589b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C",(void (**)(void))&f);CHKERRQ(ierr); 6599b54502bSHong Zhang if (f) { 6609b54502bSHong Zhang ierr = (*f)(pc,rtol);CHKERRQ(ierr); 6619b54502bSHong Zhang } 6629b54502bSHong Zhang PetscFunctionReturn(0); 6639b54502bSHong Zhang } 6649b54502bSHong Zhang 6659b54502bSHong Zhang #undef __FUNCT__ 6669b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks" 6679b54502bSHong Zhang /*@ 6689b54502bSHong Zhang PCILUSetPivotInBlocks - Determines if pivoting is done while factoring each block 6699b54502bSHong Zhang with BAIJ or SBAIJ matrices 6709b54502bSHong Zhang 6719b54502bSHong Zhang Collective on PC 6729b54502bSHong Zhang 6739b54502bSHong Zhang Input Parameters: 6749b54502bSHong Zhang + pc - the preconditioner context 6759b54502bSHong Zhang - pivot - PETSC_TRUE or PETSC_FALSE 6769b54502bSHong Zhang 6779b54502bSHong Zhang Options Database Key: 6789b54502bSHong Zhang . -pc_ilu_pivot_in_blocks <true,false> 6799b54502bSHong Zhang 6809b54502bSHong Zhang Level: intermediate 6819b54502bSHong Zhang 6829b54502bSHong Zhang .seealso: PCIILUSetMatOrdering(), PCILUSetPivoting() 6839b54502bSHong Zhang @*/ 6849b54502bSHong Zhang PetscErrorCode PCILUSetPivotInBlocks(PC pc,PetscTruth pivot) 6859b54502bSHong Zhang { 6869b54502bSHong Zhang PetscErrorCode ierr,(*f)(PC,PetscTruth); 6879b54502bSHong Zhang 6889b54502bSHong Zhang PetscFunctionBegin; 6899b54502bSHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCILUSetPivotInBlocks_C",(void (**)(void))&f);CHKERRQ(ierr); 6909b54502bSHong Zhang if (f) { 6919b54502bSHong Zhang ierr = (*f)(pc,pivot);CHKERRQ(ierr); 6929b54502bSHong Zhang } 6939b54502bSHong Zhang PetscFunctionReturn(0); 6949b54502bSHong Zhang } 6959b54502bSHong Zhang 6969b54502bSHong Zhang /* ------------------------------------------------------------------------------------------*/ 6979b54502bSHong Zhang 6989b54502bSHong Zhang EXTERN_C_BEGIN 6999b54502bSHong Zhang #undef __FUNCT__ 7009b54502bSHong Zhang #define __FUNCT__ "PCILUSetPivotInBlocks_ILU" 7019b54502bSHong Zhang PetscErrorCode PCILUSetPivotInBlocks_ILU(PC pc,PetscTruth pivot) 7029b54502bSHong Zhang { 7039b54502bSHong Zhang PC_ILU *dir = (PC_ILU*)pc->data; 7049b54502bSHong Zhang 7059b54502bSHong Zhang PetscFunctionBegin; 7069b54502bSHong Zhang dir->info.pivotinblocks = pivot ? 1.0 : 0.0; 7079b54502bSHong Zhang PetscFunctionReturn(0); 7089b54502bSHong Zhang } 7099b54502bSHong Zhang EXTERN_C_END 7109b54502bSHong Zhang 7119b54502bSHong Zhang #undef __FUNCT__ 7129b54502bSHong Zhang #define __FUNCT__ "PCSetFromOptions_ILU" 7139b54502bSHong Zhang static PetscErrorCode PCSetFromOptions_ILU(PC pc) 7149b54502bSHong Zhang { 7159b54502bSHong Zhang PetscErrorCode ierr; 7169b54502bSHong Zhang PetscInt dtmax = 3,itmp; 7179b54502bSHong Zhang PetscTruth flg,set; 7189b54502bSHong Zhang PetscReal dt[3]; 7199b54502bSHong Zhang char tname[256]; 7209b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 7219b54502bSHong Zhang PetscFList ordlist; 7229b54502bSHong Zhang PetscReal tol; 7239b54502bSHong Zhang 7249b54502bSHong Zhang PetscFunctionBegin; 7259b54502bSHong Zhang ierr = MatOrderingRegisterAll(PETSC_NULL);CHKERRQ(ierr); 7269b54502bSHong Zhang ierr = PetscOptionsHead("ILU Options");CHKERRQ(ierr); 7279b54502bSHong Zhang ierr = PetscOptionsInt("-pc_ilu_levels","levels of fill","PCILUSetLevels",(PetscInt)ilu->info.levels,&itmp,&flg);CHKERRQ(ierr); 7289b54502bSHong Zhang if (flg) ilu->info.levels = itmp; 7299b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_in_place","do factorization in place","PCILUSetUseInPlace",&ilu->inplace);CHKERRQ(ierr); 7309b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_diagonal_fill","Allow fill into empty diagonal entry","PCILUSetAllowDiagonalFill",&flg);CHKERRQ(ierr); 7319b54502bSHong Zhang ilu->info.diagonal_fill = (double) flg; 7329b54502bSHong Zhang ierr = PetscOptionsName("-pc_iludt_reuse_fill","Reuse fill from previous ILUdt","PCILUDTSetReuseFill",&ilu->reusefill);CHKERRQ(ierr); 7339b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_reuse_ordering","Reuse previous reordering","PCILUSetReuseOrdering",&ilu->reuseordering);CHKERRQ(ierr); 734*5e8efad8SHong Zhang 735*5e8efad8SHong Zhang ierr = PetscOptionsName("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",&flg);CHKERRQ(ierr); 7369b54502bSHong Zhang if (flg) { 737*5e8efad8SHong Zhang ierr = PCFactorSetShiftNonzero((PetscReal) PETSC_DECIDE,&ilu->info);CHKERRQ(ierr); 7389b54502bSHong Zhang } 739*5e8efad8SHong Zhang ierr = PetscOptionsReal("-pc_factor_shiftnonzero","Shift added to diagonal","PCFactorSetShiftNonzero",ilu->info.damping,&ilu->info.damping,0);CHKERRQ(ierr); 740*5e8efad8SHong Zhang 7419b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",&flg);CHKERRQ(ierr); 7429b54502bSHong Zhang if (flg) { 7439b54502bSHong Zhang ierr = PetscOptionsInt("-pc_ilu_shift","Manteuffel shift applied to diagonal","PCILUSetShift",(PetscInt)ilu->info.shift,&itmp,&flg); CHKERRQ(ierr); 7449b54502bSHong Zhang if (flg && !itmp) { 7459b54502bSHong Zhang ierr = PCILUSetShift(pc,PETSC_FALSE);CHKERRQ(ierr); 7469b54502bSHong Zhang } else { 7479b54502bSHong Zhang ierr = PCILUSetShift(pc,PETSC_TRUE);CHKERRQ(ierr); 7489b54502bSHong Zhang } 7499b54502bSHong Zhang } 7509b54502bSHong Zhang ierr = PetscOptionsReal("-pc_ilu_zeropivot","Pivot is considered zero if less than","PCILUSetSetZeroPivot",ilu->info.zeropivot,&ilu->info.zeropivot,0);CHKERRQ(ierr); 7519b54502bSHong Zhang 7529b54502bSHong Zhang dt[0] = ilu->info.dt; 7539b54502bSHong Zhang dt[1] = ilu->info.dtcol; 7549b54502bSHong Zhang dt[2] = ilu->info.dtcount; 7559b54502bSHong Zhang ierr = PetscOptionsRealArray("-pc_ilu_use_drop_tolerance","<dt,dtcol,maxrowcount>","PCILUSetUseDropTolerance",dt,&dtmax,&flg);CHKERRQ(ierr); 7569b54502bSHong Zhang if (flg) { 7579b54502bSHong Zhang ierr = PCILUSetUseDropTolerance(pc,dt[0],dt[1],(PetscInt)dt[2]);CHKERRQ(ierr); 7589b54502bSHong Zhang } 7599b54502bSHong Zhang ierr = PetscOptionsReal("-pc_ilu_fill","Expected fill in factorization","PCILUSetFill",ilu->info.fill,&ilu->info.fill,&flg);CHKERRQ(ierr); 7609b54502bSHong Zhang ierr = PetscOptionsName("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",&flg);CHKERRQ(ierr); 7619b54502bSHong Zhang if (flg) { 7629b54502bSHong Zhang tol = PETSC_DECIDE; 7639b54502bSHong Zhang ierr = PetscOptionsReal("-pc_ilu_nonzeros_along_diagonal","Reorder to remove zeros from diagonal","PCILUReorderForNonzeroDiagonal",ilu->nonzerosalongdiagonaltol,&tol,0);CHKERRQ(ierr); 7649b54502bSHong Zhang ierr = PCILUReorderForNonzeroDiagonal(pc,tol);CHKERRQ(ierr); 7659b54502bSHong Zhang } 7669b54502bSHong Zhang 7679b54502bSHong Zhang ierr = MatGetOrderingList(&ordlist);CHKERRQ(ierr); 7689b54502bSHong Zhang ierr = PetscOptionsList("-pc_ilu_mat_ordering_type","Reorder to reduce nonzeros in ILU","PCILUSetMatOrdering",ordlist,ilu->ordering,tname,256,&flg);CHKERRQ(ierr); 7699b54502bSHong Zhang if (flg) { 7709b54502bSHong Zhang ierr = PCILUSetMatOrdering(pc,tname);CHKERRQ(ierr); 7719b54502bSHong Zhang } 7729b54502bSHong Zhang flg = ilu->info.pivotinblocks ? PETSC_TRUE : PETSC_FALSE; 7739b54502bSHong Zhang ierr = PetscOptionsLogical("-pc_ilu_pivot_in_blocks","Pivot inside matrix blocks for BAIJ and SBAIJ","PCILUSetPivotInBlocks",flg,&flg,&set);CHKERRQ(ierr); 7749b54502bSHong Zhang if (set) { 7759b54502bSHong Zhang ierr = PCILUSetPivotInBlocks(pc,flg);CHKERRQ(ierr); 7769b54502bSHong Zhang } 7779b54502bSHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 7789b54502bSHong Zhang PetscFunctionReturn(0); 7799b54502bSHong Zhang } 7809b54502bSHong Zhang 7819b54502bSHong Zhang #undef __FUNCT__ 7829b54502bSHong Zhang #define __FUNCT__ "PCView_ILU" 7839b54502bSHong Zhang static PetscErrorCode PCView_ILU(PC pc,PetscViewer viewer) 7849b54502bSHong Zhang { 7859b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 7869b54502bSHong Zhang PetscErrorCode ierr; 7879b54502bSHong Zhang PetscTruth isstring,iascii; 7889b54502bSHong Zhang 7899b54502bSHong Zhang PetscFunctionBegin; 7909b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 7919b54502bSHong Zhang ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 7929b54502bSHong Zhang if (iascii) { 7939b54502bSHong Zhang if (ilu->usedt) { 7949b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: drop tolerance %g\n",ilu->info.dt);CHKERRQ(ierr); 7959b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: max nonzeros per row %D\n",(PetscInt)ilu->info.dtcount);CHKERRQ(ierr); 7969b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: column permutation tolerance %g\n",ilu->info.dtcol);CHKERRQ(ierr); 7979b54502bSHong Zhang } else if (ilu->info.levels == 1) { 7989b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D level of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 7999b54502bSHong Zhang } else { 8009b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: %D levels of fill\n",(PetscInt)ilu->info.levels);CHKERRQ(ierr); 8019b54502bSHong Zhang } 8029b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: max fill ratio allocated %g\n",ilu->info.fill);CHKERRQ(ierr); 8039b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," ILU: tolerance for zero pivot %g\n",ilu->info.zeropivot);CHKERRQ(ierr); 8049b54502bSHong Zhang if (ilu->info.shift) {ierr = PetscViewerASCIIPrintf(viewer," ILU: using Manteuffel shift\n");CHKERRQ(ierr);} 8059b54502bSHong Zhang if (ilu->inplace) {ierr = PetscViewerASCIIPrintf(viewer," in-place factorization\n");CHKERRQ(ierr);} 8069b54502bSHong Zhang else {ierr = PetscViewerASCIIPrintf(viewer," out-of-place factorization\n");CHKERRQ(ierr);} 8079b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," matrix ordering: %s\n",ilu->ordering);CHKERRQ(ierr); 8089b54502bSHong Zhang if (ilu->reusefill) {ierr = PetscViewerASCIIPrintf(viewer," Reusing fill from past factorization\n");CHKERRQ(ierr);} 8099b54502bSHong Zhang if (ilu->reuseordering) {ierr = PetscViewerASCIIPrintf(viewer," Reusing reordering from past factorization\n");CHKERRQ(ierr);} 8109b54502bSHong Zhang if (ilu->fact) { 8119b54502bSHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Factored matrix follows\n");CHKERRQ(ierr); 8129b54502bSHong Zhang ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 8139b54502bSHong Zhang ierr = PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);CHKERRQ(ierr); 8149b54502bSHong Zhang ierr = MatView(ilu->fact,viewer);CHKERRQ(ierr); 8159b54502bSHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 8169b54502bSHong Zhang ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 8179b54502bSHong Zhang } 8189b54502bSHong Zhang } else if (isstring) { 8199b54502bSHong Zhang ierr = PetscViewerStringSPrintf(viewer," lvls=%D,order=%s",(PetscInt)ilu->info.levels,ilu->ordering);CHKERRQ(ierr);CHKERRQ(ierr); 8209b54502bSHong Zhang } else { 8219b54502bSHong Zhang SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCILU",((PetscObject)viewer)->type_name); 8229b54502bSHong Zhang } 8239b54502bSHong Zhang PetscFunctionReturn(0); 8249b54502bSHong Zhang } 8259b54502bSHong Zhang 8269b54502bSHong Zhang #undef __FUNCT__ 8279b54502bSHong Zhang #define __FUNCT__ "PCSetUp_ILU" 8289b54502bSHong Zhang static PetscErrorCode PCSetUp_ILU(PC pc) 8299b54502bSHong Zhang { 8309b54502bSHong Zhang PetscErrorCode ierr; 8319b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 8329b54502bSHong Zhang 8339b54502bSHong Zhang PetscFunctionBegin; 8349b54502bSHong Zhang if (ilu->inplace) { 8359b54502bSHong Zhang if (!pc->setupcalled) { 8369b54502bSHong Zhang 8379b54502bSHong Zhang /* In-place factorization only makes sense with the natural ordering, 8389b54502bSHong Zhang so we only need to get the ordering once, even if nonzero structure changes */ 8399b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 8409b54502bSHong Zhang if (ilu->row) PetscLogObjectParent(pc,ilu->row); 8419b54502bSHong Zhang if (ilu->col) PetscLogObjectParent(pc,ilu->col); 8429b54502bSHong Zhang } 8439b54502bSHong Zhang 8449b54502bSHong Zhang /* In place ILU only makes sense with fill factor of 1.0 because 8459b54502bSHong Zhang cannot have levels of fill */ 8469b54502bSHong Zhang ilu->info.fill = 1.0; 8479b54502bSHong Zhang ilu->info.diagonal_fill = 0; 8489b54502bSHong Zhang ierr = MatILUFactor(pc->pmat,ilu->row,ilu->col,&ilu->info);CHKERRQ(ierr); 8499b54502bSHong Zhang ilu->fact = pc->pmat; 8509b54502bSHong Zhang } else if (ilu->usedt) { 8519b54502bSHong Zhang if (!pc->setupcalled) { 8529b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 8539b54502bSHong Zhang if (ilu->row) PetscLogObjectParent(pc,ilu->row); 8549b54502bSHong Zhang if (ilu->col) PetscLogObjectParent(pc,ilu->col); 8559b54502bSHong Zhang ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 8569b54502bSHong Zhang PetscLogObjectParent(pc,ilu->fact); 8579b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 8589b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 8599b54502bSHong Zhang if (!ilu->reuseordering) { 8609b54502bSHong Zhang if (ilu->row) {ierr = ISDestroy(ilu->row);CHKERRQ(ierr);} 8619b54502bSHong Zhang if (ilu->col) {ierr = ISDestroy(ilu->col);CHKERRQ(ierr);} 8629b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 8639b54502bSHong Zhang if (ilu->row) PetscLogObjectParent(pc,ilu->row); 8649b54502bSHong Zhang if (ilu->col) PetscLogObjectParent(pc,ilu->col); 8659b54502bSHong Zhang } 8669b54502bSHong Zhang ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 8679b54502bSHong Zhang PetscLogObjectParent(pc,ilu->fact); 8689b54502bSHong Zhang } else if (!ilu->reusefill) { 8699b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 8709b54502bSHong Zhang ierr = MatILUDTFactor(pc->pmat,&ilu->info,ilu->row,ilu->col,&ilu->fact);CHKERRQ(ierr); 8719b54502bSHong Zhang PetscLogObjectParent(pc,ilu->fact); 8729b54502bSHong Zhang } else { 8739b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 8749b54502bSHong Zhang } 8759b54502bSHong Zhang } else { 8769b54502bSHong Zhang if (!pc->setupcalled) { 8779b54502bSHong Zhang /* first time in so compute reordering and symbolic factorization */ 8789b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 8799b54502bSHong Zhang if (ilu->row) PetscLogObjectParent(pc,ilu->row); 8809b54502bSHong Zhang if (ilu->col) PetscLogObjectParent(pc,ilu->col); 8819b54502bSHong Zhang /* Remove zeros along diagonal? */ 8829b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 8839b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 8849b54502bSHong Zhang } 8859b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 8869b54502bSHong Zhang PetscLogObjectParent(pc,ilu->fact); 8879b54502bSHong Zhang } else if (pc->flag != SAME_NONZERO_PATTERN) { 8889b54502bSHong Zhang if (!ilu->reuseordering) { 8899b54502bSHong Zhang /* compute a new ordering for the ILU */ 8909b54502bSHong Zhang ierr = ISDestroy(ilu->row);CHKERRQ(ierr); 8919b54502bSHong Zhang ierr = ISDestroy(ilu->col);CHKERRQ(ierr); 8929b54502bSHong Zhang ierr = MatGetOrdering(pc->pmat,ilu->ordering,&ilu->row,&ilu->col);CHKERRQ(ierr); 8939b54502bSHong Zhang if (ilu->row) PetscLogObjectParent(pc,ilu->row); 8949b54502bSHong Zhang if (ilu->col) PetscLogObjectParent(pc,ilu->col); 8959b54502bSHong Zhang /* Remove zeros along diagonal? */ 8969b54502bSHong Zhang if (ilu->nonzerosalongdiagonal) { 8979b54502bSHong Zhang ierr = MatReorderForNonzeroDiagonal(pc->pmat,ilu->nonzerosalongdiagonaltol,ilu->row,ilu->col);CHKERRQ(ierr); 8989b54502bSHong Zhang } 8999b54502bSHong Zhang } 9009b54502bSHong Zhang ierr = MatDestroy(ilu->fact);CHKERRQ(ierr); 9019b54502bSHong Zhang ierr = MatILUFactorSymbolic(pc->pmat,ilu->row,ilu->col,&ilu->info,&ilu->fact);CHKERRQ(ierr); 9029b54502bSHong Zhang PetscLogObjectParent(pc,ilu->fact); 9039b54502bSHong Zhang } 9049b54502bSHong Zhang ierr = MatLUFactorNumeric(pc->pmat,&ilu->info,&ilu->fact);CHKERRQ(ierr); 9059b54502bSHong Zhang } 9069b54502bSHong Zhang PetscFunctionReturn(0); 9079b54502bSHong Zhang } 9089b54502bSHong Zhang 9099b54502bSHong Zhang #undef __FUNCT__ 9109b54502bSHong Zhang #define __FUNCT__ "PCDestroy_ILU" 9119b54502bSHong Zhang static PetscErrorCode PCDestroy_ILU(PC pc) 9129b54502bSHong Zhang { 9139b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 9149b54502bSHong Zhang PetscErrorCode ierr; 9159b54502bSHong Zhang 9169b54502bSHong Zhang PetscFunctionBegin; 9179b54502bSHong Zhang ierr = PCDestroy_ILU_Internal(pc);CHKERRQ(ierr); 9189b54502bSHong Zhang ierr = PetscStrfree(ilu->ordering);CHKERRQ(ierr); 9199b54502bSHong Zhang ierr = PetscFree(ilu);CHKERRQ(ierr); 9209b54502bSHong Zhang PetscFunctionReturn(0); 9219b54502bSHong Zhang } 9229b54502bSHong Zhang 9239b54502bSHong Zhang #undef __FUNCT__ 9249b54502bSHong Zhang #define __FUNCT__ "PCApply_ILU" 9259b54502bSHong Zhang static PetscErrorCode PCApply_ILU(PC pc,Vec x,Vec y) 9269b54502bSHong Zhang { 9279b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 9289b54502bSHong Zhang PetscErrorCode ierr; 9299b54502bSHong Zhang 9309b54502bSHong Zhang PetscFunctionBegin; 9319b54502bSHong Zhang ierr = MatSolve(ilu->fact,x,y);CHKERRQ(ierr); 9329b54502bSHong Zhang PetscFunctionReturn(0); 9339b54502bSHong Zhang } 9349b54502bSHong Zhang 9359b54502bSHong Zhang #undef __FUNCT__ 9369b54502bSHong Zhang #define __FUNCT__ "PCApplyTranspose_ILU" 9379b54502bSHong Zhang static PetscErrorCode PCApplyTranspose_ILU(PC pc,Vec x,Vec y) 9389b54502bSHong Zhang { 9399b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 9409b54502bSHong Zhang PetscErrorCode ierr; 9419b54502bSHong Zhang 9429b54502bSHong Zhang PetscFunctionBegin; 9439b54502bSHong Zhang ierr = MatSolveTranspose(ilu->fact,x,y);CHKERRQ(ierr); 9449b54502bSHong Zhang PetscFunctionReturn(0); 9459b54502bSHong Zhang } 9469b54502bSHong Zhang 9479b54502bSHong Zhang #undef __FUNCT__ 9489b54502bSHong Zhang #define __FUNCT__ "PCGetFactoredMatrix_ILU" 9499b54502bSHong Zhang static PetscErrorCode PCGetFactoredMatrix_ILU(PC pc,Mat *mat) 9509b54502bSHong Zhang { 9519b54502bSHong Zhang PC_ILU *ilu = (PC_ILU*)pc->data; 9529b54502bSHong Zhang 9539b54502bSHong Zhang PetscFunctionBegin; 9549b54502bSHong Zhang if (!ilu->fact) SETERRQ(PETSC_ERR_ORDER,"Matrix not yet factored; call after KSPSetUp() or PCSetUp()"); 9559b54502bSHong Zhang *mat = ilu->fact; 9569b54502bSHong Zhang PetscFunctionReturn(0); 9579b54502bSHong Zhang } 9589b54502bSHong Zhang 9599b54502bSHong Zhang /*MC 9609b54502bSHong Zhang PCILU - Incomplete factorization preconditioners. 9619b54502bSHong Zhang 9629b54502bSHong Zhang Options Database Keys: 9639b54502bSHong Zhang + -pc_ilu_levels <k> - number of levels of fill for ILU(k) 9649b54502bSHong Zhang . -pc_ilu_in_place - only for ILU(0) with natural ordering, reuses the space of the matrix for 9659b54502bSHong Zhang its factorization (overwrites original matrix) 9669b54502bSHong Zhang . -pc_ilu_diagonal_fill - fill in a zero diagonal even if levels of fill indicate it wouldn't be fill 9679b54502bSHong Zhang . -pc_ilu_reuse_ordering - reuse ordering of factorized matrix from previous factorization 9689b54502bSHong Zhang . -pc_ilu_damping - add damping to diagonal to prevent zero (or very small) pivots 9699b54502bSHong Zhang . -pc_ilu_shift - apply Manteuffel shift to diagonal to force positive definite preconditioner 9709b54502bSHong Zhang . -pc_ilu_zeropivot <tol> - set tolerance for what is considered a zero pivot 9719b54502bSHong Zhang . -pc_ilu_use_drop_tolerance <dt,dtcol,maxrowcount> - use Saad's drop tolerance ILUdt 9729b54502bSHong Zhang . -pc_ilu_fill <nfill> - expected amount of fill in factored matrix compared to original matrix, nfill > 1 9739b54502bSHong Zhang . -pc_ilu_nonzeros_along_diagonal - reorder the matrix before factorization to remove zeros from the diagonal, 9749b54502bSHong Zhang this decreases the chance of getting a zero pivot 9759b54502bSHong Zhang . -pc_ilu_mat_ordering_type <natural,nd,1wd,rcm,qmd> - set the row/column ordering of the factored matrix 9769b54502bSHong Zhang - -pc_ilu_pivot_in_blocks - for block ILU(k) factorization, i.e. with BAIJ matrices with block size larger 9779b54502bSHong Zhang than 1 the diagonal blocks are factored with partial pivoting (this increases the 9789b54502bSHong Zhang stability of the ILU factorization 9799b54502bSHong Zhang 9809b54502bSHong Zhang Level: beginner 9819b54502bSHong Zhang 9829b54502bSHong Zhang Concepts: incomplete factorization 9839b54502bSHong Zhang 9849b54502bSHong Zhang Notes: Only implemented for some matrix formats. Not implemented in parallel 9859b54502bSHong Zhang 9869b54502bSHong Zhang For BAIJ matrices this implements a point block ILU 9879b54502bSHong Zhang 9889b54502bSHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCSOR, MatOrderingType, 9899b54502bSHong Zhang PCILUSetSetZeroPivot(), PCILUSetDamping(), PCILUSetShift(), PCILUSetUseDropTolerance(), 9909b54502bSHong Zhang PCILUSetFill(), PCILUSetMatOrdering(), PCILUSetReuseOrdering(), PCILUDTSetReuseFill(), 9919b54502bSHong Zhang PCILUSetLevels(), PCILUSetUseInPlace(), PCILUSetAllowDiagonalFill(), PCILUSetPivotInBlocks(), 9929b54502bSHong Zhang 9939b54502bSHong Zhang M*/ 9949b54502bSHong Zhang 9959b54502bSHong Zhang EXTERN_C_BEGIN 9969b54502bSHong Zhang #undef __FUNCT__ 9979b54502bSHong Zhang #define __FUNCT__ "PCCreate_ILU" 9989b54502bSHong Zhang PetscErrorCode PCCreate_ILU(PC pc) 9999b54502bSHong Zhang { 10009b54502bSHong Zhang PetscErrorCode ierr; 10019b54502bSHong Zhang PC_ILU *ilu; 10029b54502bSHong Zhang 10039b54502bSHong Zhang PetscFunctionBegin; 10049b54502bSHong Zhang ierr = PetscNew(PC_ILU,&ilu);CHKERRQ(ierr); 10059b54502bSHong Zhang PetscLogObjectMemory(pc,sizeof(PC_ILU)); 10069b54502bSHong Zhang 10079b54502bSHong Zhang ilu->fact = 0; 10089b54502bSHong Zhang ierr = MatFactorInfoInitialize(&ilu->info);CHKERRQ(ierr); 10099b54502bSHong Zhang ilu->info.levels = 0; 10109b54502bSHong Zhang ilu->info.fill = 1.0; 10119b54502bSHong Zhang ilu->col = 0; 10129b54502bSHong Zhang ilu->row = 0; 10139b54502bSHong Zhang ilu->inplace = PETSC_FALSE; 10149b54502bSHong Zhang ierr = PetscStrallocpy(MATORDERING_NATURAL,&ilu->ordering);CHKERRQ(ierr); 10159b54502bSHong Zhang ilu->reuseordering = PETSC_FALSE; 10169b54502bSHong Zhang ilu->usedt = PETSC_FALSE; 10179b54502bSHong Zhang ilu->info.dt = PETSC_DEFAULT; 10189b54502bSHong Zhang ilu->info.dtcount = PETSC_DEFAULT; 10199b54502bSHong Zhang ilu->info.dtcol = PETSC_DEFAULT; 10209b54502bSHong Zhang ilu->info.damping = 0.0; 10219b54502bSHong Zhang ilu->info.shift = PETSC_FALSE; 10229b54502bSHong Zhang ilu->info.shift_fraction = 0.0; 10239b54502bSHong Zhang ilu->info.zeropivot = 1.e-12; 10249b54502bSHong Zhang ilu->info.pivotinblocks = 1.0; 10259b54502bSHong Zhang ilu->reusefill = PETSC_FALSE; 10269b54502bSHong Zhang ilu->info.diagonal_fill = 0; 10279b54502bSHong Zhang pc->data = (void*)ilu; 10289b54502bSHong Zhang 10299b54502bSHong Zhang pc->ops->destroy = PCDestroy_ILU; 10309b54502bSHong Zhang pc->ops->apply = PCApply_ILU; 10319b54502bSHong Zhang pc->ops->applytranspose = PCApplyTranspose_ILU; 10329b54502bSHong Zhang pc->ops->setup = PCSetUp_ILU; 10339b54502bSHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ILU; 10349b54502bSHong Zhang pc->ops->getfactoredmatrix = PCGetFactoredMatrix_ILU; 10359b54502bSHong Zhang pc->ops->view = PCView_ILU; 10369b54502bSHong Zhang pc->ops->applyrichardson = 0; 10379b54502bSHong Zhang 10389b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseDropTolerance_C","PCILUSetUseDropTolerance_ILU", 10399b54502bSHong Zhang PCILUSetUseDropTolerance_ILU);CHKERRQ(ierr); 10409b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetFill_C","PCILUSetFill_ILU", 10419b54502bSHong Zhang PCILUSetFill_ILU);CHKERRQ(ierr); 10429b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetDamping_C","PCILUSetDamping_ILU", 10439b54502bSHong Zhang PCILUSetDamping_ILU);CHKERRQ(ierr); 10449b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetShift_C","PCILUSetShift_ILU", 10459b54502bSHong Zhang PCILUSetShift_ILU);CHKERRQ(ierr); 10469b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetMatOrdering_C","PCILUSetMatOrdering_ILU", 10479b54502bSHong Zhang PCILUSetMatOrdering_ILU);CHKERRQ(ierr); 10489b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetReuseOrdering_C","PCILUSetReuseOrdering_ILU", 10499b54502bSHong Zhang PCILUSetReuseOrdering_ILU);CHKERRQ(ierr); 10509b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUDTSetReuseFill_C","PCILUDTSetReuseFill_ILUDT", 10519b54502bSHong Zhang PCILUDTSetReuseFill_ILUDT);CHKERRQ(ierr); 10529b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetLevels_C","PCILUSetLevels_ILU", 10539b54502bSHong Zhang PCILUSetLevels_ILU);CHKERRQ(ierr); 10549b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetUseInPlace_C","PCILUSetUseInPlace_ILU", 10559b54502bSHong Zhang PCILUSetUseInPlace_ILU);CHKERRQ(ierr); 10569b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetAllowDiagonalFill_C","PCILUSetAllowDiagonalFill_ILU", 10579b54502bSHong Zhang PCILUSetAllowDiagonalFill_ILU);CHKERRQ(ierr); 10589b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetPivotInBlocks_C","PCILUSetPivotInBlocks_ILU", 10599b54502bSHong Zhang PCILUSetPivotInBlocks_ILU);CHKERRQ(ierr); 10609b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUSetZeroPivot_C","PCILUSetZeroPivot_ILU", 10619b54502bSHong Zhang PCILUSetZeroPivot_ILU);CHKERRQ(ierr); 10629b54502bSHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCILUReorderForNonzeroDiagonal_C","PCILUReorderForNonzeroDiagonal_ILU", 10639b54502bSHong Zhang PCILUReorderForNonzeroDiagonal_ILU);CHKERRQ(ierr); 10649b54502bSHong Zhang PetscFunctionReturn(0); 10659b54502bSHong Zhang } 10669b54502bSHong Zhang EXTERN_C_END 1067