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