1c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 2978beb7fSPierre Jolivet #include <petsc/private/matimpl.h> 35e8efad8SHong Zhang 44ac6704cSBarry Smith /* 54ac6704cSBarry Smith If an ordering is not yet set and the matrix is available determine a default ordering 64ac6704cSBarry Smith */ 7d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc) 8d71ae5a4SJacob Faibussowitsch { 97def7855SStefano Zampini PetscBool foundmtype, flg, destroy = PETSC_FALSE; 1026cc229bSBarry Smith const char *prefix; 114ac6704cSBarry Smith 124ac6704cSBarry Smith PetscFunctionBegin; 134ac6704cSBarry Smith if (pc->pmat) { 1426cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1526cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 164ac6704cSBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 179566063dSJacob Faibussowitsch PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL)); 18978beb7fSPierre Jolivet if (foundmtype) { 197def7855SStefano Zampini if (!fact->fact) { 207def7855SStefano Zampini /* factored matrix is not present at this point, we want to create it during PCSetUp. 217def7855SStefano Zampini Since this can be called from setfromoptions, we destroy it when we are done with it */ 227def7855SStefano Zampini PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); 237def7855SStefano Zampini destroy = PETSC_TRUE; 247def7855SStefano Zampini } 2503e5aca4SStefano Zampini if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS); 2603e5aca4SStefano Zampini if (!fact->fact->assembled) { 279566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg)); 28978beb7fSPierre Jolivet if (!flg) { 2926cc229bSBarry Smith Mat B; 309566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B)); 319566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(fact->fact, &B)); 32978beb7fSPierre Jolivet } 334ac6704cSBarry Smith } 344ac6704cSBarry Smith if (!fact->ordering) { 354ac6704cSBarry Smith PetscBool canuseordering; 364ac6704cSBarry Smith MatOrderingType otype; 374ac6704cSBarry Smith 389566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering)); 39*ac530a7eSPierre Jolivet if (canuseordering) PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype)); 40*ac530a7eSPierre Jolivet else otype = MATORDERINGEXTERNAL; 419566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering)); 424ac6704cSBarry Smith } 437def7855SStefano Zampini if (destroy) PetscCall(MatDestroy(&fact->fact)); 444ac6704cSBarry Smith } 45978beb7fSPierre Jolivet } 463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 474ac6704cSBarry Smith } 484ac6704cSBarry Smith 49d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag) 50d71ae5a4SJacob Faibussowitsch { 513d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 523d1c1ea0SBarry Smith 533d1c1ea0SBarry Smith PetscFunctionBegin; 543d1c1ea0SBarry Smith lu->reuseordering = flag; 553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 563d1c1ea0SBarry Smith } 573d1c1ea0SBarry Smith 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag) 59d71ae5a4SJacob Faibussowitsch { 603d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 613d1c1ea0SBarry Smith 623d1c1ea0SBarry Smith PetscFunctionBegin; 633d1c1ea0SBarry Smith lu->reusefill = flag; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 653d1c1ea0SBarry Smith } 663d1c1ea0SBarry Smith 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg) 68d71ae5a4SJacob Faibussowitsch { 693d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 703d1c1ea0SBarry Smith 713d1c1ea0SBarry Smith PetscFunctionBegin; 723d1c1ea0SBarry Smith dir->inplace = flg; 733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 743d1c1ea0SBarry Smith } 753d1c1ea0SBarry Smith 76d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg) 77d71ae5a4SJacob Faibussowitsch { 783d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 793d1c1ea0SBarry Smith 803d1c1ea0SBarry Smith PetscFunctionBegin; 813d1c1ea0SBarry Smith *flg = dir->inplace; 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 833d1c1ea0SBarry Smith } 843d1c1ea0SBarry Smith 85f8260c8fSBarry Smith /*@ 86f1580f4eSBarry Smith PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may 87f8260c8fSBarry Smith set the options for that particular factorization object. 88f8260c8fSBarry Smith 89f8260c8fSBarry Smith Input Parameter: 90f8260c8fSBarry Smith . pc - the preconditioner context 91f8260c8fSBarry Smith 92f1580f4eSBarry Smith Note: 93f1580f4eSBarry Smith After you have called this function (which has to be after the `KSPSetOperators()` or `PCSetOperators()`) you can call `PCFactorGetMatrix()` and then set factor options on that matrix. 9403e5aca4SStefano Zampini This function raises an error if the requested combination of solver package and matrix type is not supported. 95f8260c8fSBarry Smith 962bd2b0e6SSatish Balay Level: intermediate 972bd2b0e6SSatish Balay 98562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 99f8260c8fSBarry Smith @*/ 100d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 101d71ae5a4SJacob Faibussowitsch { 102f8260c8fSBarry Smith PetscFunctionBegin; 103f8260c8fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 104cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106f8260c8fSBarry Smith } 107f8260c8fSBarry Smith 108ee45ca4aSHong Zhang /*@ 109ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 110ee45ca4aSHong Zhang 111c3339decSBarry Smith Logically Collective 112ee45ca4aSHong Zhang 113ee45ca4aSHong Zhang Input Parameters: 114afaefe49SHong Zhang + pc - the preconditioner context 115afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 116ee45ca4aSHong Zhang 117ee45ca4aSHong Zhang Options Database Key: 118ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 119ee45ca4aSHong Zhang 120ee45ca4aSHong Zhang Level: intermediate 121ee45ca4aSHong Zhang 122562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 123ee45ca4aSHong Zhang @*/ 124d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) 125d71ae5a4SJacob Faibussowitsch { 126ee45ca4aSHong Zhang PetscFunctionBegin; 1270700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 128c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, zero, 2); 129cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 1303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 131ee45ca4aSHong Zhang } 132ee45ca4aSHong Zhang 133915743fcSHong Zhang /*@ 134915743fcSHong Zhang PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 135915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 136915743fcSHong Zhang 137c3339decSBarry Smith Logically Collective 138915743fcSHong Zhang 139915743fcSHong Zhang Input Parameters: 140915743fcSHong Zhang + pc - the preconditioner context 141f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS` 142915743fcSHong Zhang 143915743fcSHong Zhang Options Database Key: 14428d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 145915743fcSHong Zhang 146915743fcSHong Zhang Level: intermediate 147915743fcSHong Zhang 148562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 149915743fcSHong Zhang @*/ 150d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) 151d71ae5a4SJacob Faibussowitsch { 152d90ac83dSHong Zhang PetscFunctionBegin; 1530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 154c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 155cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 1563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 157d90ac83dSHong Zhang } 158d90ac83dSHong Zhang 159915743fcSHong Zhang /*@ 160915743fcSHong Zhang PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 161915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 162915743fcSHong Zhang 163c3339decSBarry Smith Logically Collective 164915743fcSHong Zhang 165915743fcSHong Zhang Input Parameters: 166915743fcSHong Zhang + pc - the preconditioner context 167f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default 168915743fcSHong Zhang 169915743fcSHong Zhang Options Database Key: 170f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 171915743fcSHong Zhang 172915743fcSHong Zhang Level: intermediate 173915743fcSHong Zhang 174a94f484eSPierre Jolivet .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 175915743fcSHong Zhang @*/ 176d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) 177d71ae5a4SJacob Faibussowitsch { 178d90ac83dSHong Zhang PetscFunctionBegin; 1790700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 180c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 181cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183d90ac83dSHong Zhang } 184d90ac83dSHong Zhang 1856d33885cSprj- /*@ 186f1580f4eSBarry Smith PCFactorSetDropTolerance - The preconditioner will use an `PCILU` 187f1580f4eSBarry Smith based on a drop tolerance. 18885317021SBarry Smith 189c3339decSBarry Smith Logically Collective 19085317021SBarry Smith 19185317021SBarry Smith Input Parameters: 19285317021SBarry Smith + pc - the preconditioner context 19385317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 19485317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 19585317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 19685317021SBarry Smith depends on the number of nonzeros in row of original matrix 19785317021SBarry Smith 19885317021SBarry Smith Options Database Key: 199b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 20085317021SBarry Smith 20185317021SBarry Smith Level: intermediate 20285317021SBarry Smith 203f1580f4eSBarry Smith Note: 20485317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 20585317021SBarry Smith matrix. We don't know how to compute reasonable values. 20685317021SBarry Smith 207562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU` 2086d33885cSprj- @*/ 209d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) 210d71ae5a4SJacob Faibussowitsch { 21185317021SBarry Smith PetscFunctionBegin; 2120700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 213064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveReal(pc, dtcol, 3); 214064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 215cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21785317021SBarry Smith } 21885317021SBarry Smith 219c7f610a1SBarry Smith /*@ 220c7f610a1SBarry Smith PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 221c7f610a1SBarry Smith 222c7f610a1SBarry Smith Not Collective 223c7f610a1SBarry Smith 2242fe279fdSBarry Smith Input Parameter: 225c7f610a1SBarry Smith . pc - the preconditioner context 226c7f610a1SBarry Smith 227c7f610a1SBarry Smith Output Parameter: 228c7f610a1SBarry Smith . pivot - the tolerance 229c7f610a1SBarry Smith 230c7f610a1SBarry Smith Level: intermediate 231c7f610a1SBarry Smith 232562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 233c7f610a1SBarry Smith @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) 235d71ae5a4SJacob Faibussowitsch { 236c7f610a1SBarry Smith PetscFunctionBegin; 237c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 238cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 240c7f610a1SBarry Smith } 241c7f610a1SBarry Smith 242c7f610a1SBarry Smith /*@ 243c7f610a1SBarry Smith PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 244c7f610a1SBarry Smith 245c7f610a1SBarry Smith Not Collective 246c7f610a1SBarry Smith 2472fe279fdSBarry Smith Input Parameter: 248c7f610a1SBarry Smith . pc - the preconditioner context 249c7f610a1SBarry Smith 250c7f610a1SBarry Smith Output Parameter: 251c7f610a1SBarry Smith . shift - how much to shift the diagonal entry 252c7f610a1SBarry Smith 253c7f610a1SBarry Smith Level: intermediate 254c7f610a1SBarry Smith 255562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 256c7f610a1SBarry Smith @*/ 257d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) 258d71ae5a4SJacob Faibussowitsch { 259c7f610a1SBarry Smith PetscFunctionBegin; 260c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 261cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 2623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 263c7f610a1SBarry Smith } 264c7f610a1SBarry Smith 265c7f610a1SBarry Smith /*@ 266c7f610a1SBarry Smith PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 267c7f610a1SBarry Smith 268c7f610a1SBarry Smith Not Collective 269c7f610a1SBarry Smith 270f1580f4eSBarry Smith Input Parameter: 271c7f610a1SBarry Smith . pc - the preconditioner context 272c7f610a1SBarry Smith 273c7f610a1SBarry Smith Output Parameter: 274f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 275c7f610a1SBarry Smith 276c7f610a1SBarry Smith Level: intermediate 277c7f610a1SBarry Smith 278562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 279c7f610a1SBarry Smith @*/ 280d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) 281d71ae5a4SJacob Faibussowitsch { 282c7f610a1SBarry Smith PetscFunctionBegin; 283c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 284cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 2853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 286c7f610a1SBarry Smith } 287c7f610a1SBarry Smith 2882591b318SToby Isaac /*@ 2892591b318SToby Isaac PCFactorGetLevels - Gets the number of levels of fill to use. 2902591b318SToby Isaac 291c3339decSBarry Smith Logically Collective 2922591b318SToby Isaac 2932fe279fdSBarry Smith Input Parameter: 2942591b318SToby Isaac . pc - the preconditioner context 2952591b318SToby Isaac 2962591b318SToby Isaac Output Parameter: 2972591b318SToby Isaac . levels - number of levels of fill 2982591b318SToby Isaac 2992591b318SToby Isaac Level: intermediate 3002591b318SToby Isaac 301562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetLevels()` 3022591b318SToby Isaac @*/ 303d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) 304d71ae5a4SJacob Faibussowitsch { 3052591b318SToby Isaac PetscFunctionBegin; 3062591b318SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 307cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 3083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3092591b318SToby Isaac } 3102591b318SToby Isaac 31185317021SBarry Smith /*@ 31285317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 31385317021SBarry Smith 314c3339decSBarry Smith Logically Collective 31585317021SBarry Smith 31685317021SBarry Smith Input Parameters: 31785317021SBarry Smith + pc - the preconditioner context 31885317021SBarry Smith - levels - number of levels of fill 31985317021SBarry Smith 32085317021SBarry Smith Options Database Key: 32185317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 32285317021SBarry Smith 32385317021SBarry Smith Level: intermediate 32485317021SBarry Smith 325562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetLevels()` 32685317021SBarry Smith @*/ 327d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) 328d71ae5a4SJacob Faibussowitsch { 32985317021SBarry Smith PetscFunctionBegin; 3300700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3315f80ce2aSJacob Faibussowitsch PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 332c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, levels, 2); 333cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 3343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33585317021SBarry Smith } 33685317021SBarry Smith 33785317021SBarry Smith /*@ 33885317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 33985317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 34085317021SBarry Smith 341c3339decSBarry Smith Logically Collective 34285317021SBarry Smith 34385317021SBarry Smith Input Parameters: 34485317021SBarry Smith + pc - the preconditioner context 345f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 34685317021SBarry Smith 34785317021SBarry Smith Options Database Key: 34867b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 34985317021SBarry Smith 350f1580f4eSBarry Smith Note: 35185317021SBarry Smith Does not apply with 0 fill. 35285317021SBarry Smith 35385317021SBarry Smith Level: intermediate 35485317021SBarry Smith 355562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 35685317021SBarry Smith @*/ 357d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) 358d71ae5a4SJacob Faibussowitsch { 35985317021SBarry Smith PetscFunctionBegin; 3600700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 361cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36392e9c092SBarry Smith } 36492e9c092SBarry Smith 36592e9c092SBarry Smith /*@ 36692e9c092SBarry Smith PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 36792e9c092SBarry Smith treated as level 0 fill even if there is no non-zero location. 36892e9c092SBarry Smith 369c3339decSBarry Smith Logically Collective 37092e9c092SBarry Smith 37192e9c092SBarry Smith Input Parameter: 37292e9c092SBarry Smith . pc - the preconditioner context 37392e9c092SBarry Smith 37492e9c092SBarry Smith Output Parameter: 375f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 37692e9c092SBarry Smith 377f1580f4eSBarry Smith Note: 37892e9c092SBarry Smith Does not apply with 0 fill. 37992e9c092SBarry Smith 38092e9c092SBarry Smith Level: intermediate 38192e9c092SBarry Smith 382562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 38392e9c092SBarry Smith @*/ 384d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) 385d71ae5a4SJacob Faibussowitsch { 38692e9c092SBarry Smith PetscFunctionBegin; 38792e9c092SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 388cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39085317021SBarry Smith } 39185317021SBarry Smith 39285317021SBarry Smith /*@ 39385317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 39485317021SBarry Smith 395c3339decSBarry Smith Logically Collective 39685317021SBarry Smith 39785317021SBarry Smith Input Parameters: 39885317021SBarry Smith + pc - the preconditioner context 399feefa0e1SJacob Faibussowitsch - rtol - diagonal entries smaller than this in absolute value are considered zero 40085317021SBarry Smith 40185317021SBarry Smith Options Database Key: 402147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 40385317021SBarry Smith 40485317021SBarry Smith Level: intermediate 40585317021SBarry Smith 406ef959800SPierre Jolivet .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftAmount()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 40785317021SBarry Smith @*/ 408d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) 409d71ae5a4SJacob Faibussowitsch { 41085317021SBarry Smith PetscFunctionBegin; 4110700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 412c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, rtol, 2); 413cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 4143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41585317021SBarry Smith } 41685317021SBarry Smith 417cc4c1da9SBarry Smith /*@ 418f1580f4eSBarry Smith PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 41985317021SBarry Smith 420c3339decSBarry Smith Logically Collective 42185317021SBarry Smith 42285317021SBarry Smith Input Parameters: 42385317021SBarry Smith + pc - the preconditioner context 424f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 42585317021SBarry Smith 42685317021SBarry Smith Options Database Key: 4273ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 42885317021SBarry Smith 42985317021SBarry Smith Level: intermediate 43085317021SBarry Smith 43185317021SBarry Smith Note: 432998e4596SBarry Smith The default type is set by searching for available types based on the order of the calls to `MatSolverTypeRegister()` in `MatInitializePackage()`. 433998e4596SBarry Smith Since different PETSc configurations may have different external solvers, seemingly identical runs with different PETSc configurations may use a different solver. 434998e4596SBarry Smith For example if one configuration had --download-mumps while a different one had --download-superlu_dist. 43585317021SBarry Smith 436998e4596SBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverTypeRegister()`, 437998e4596SBarry Smith `MatInitializePackage()`, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS`, `MatSolverTypeGet()` 43885317021SBarry Smith @*/ 439d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) 440d71ae5a4SJacob Faibussowitsch { 44185317021SBarry Smith PetscFunctionBegin; 4420700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 443cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44585317021SBarry Smith } 44685317021SBarry Smith 447cc4c1da9SBarry Smith /*@ 448f1580f4eSBarry Smith PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization 4497112b564SBarry Smith 450c5eb9154SBarry Smith Not Collective 4517112b564SBarry Smith 4527112b564SBarry Smith Input Parameter: 4537112b564SBarry Smith . pc - the preconditioner context 4547112b564SBarry Smith 4557112b564SBarry Smith Output Parameter: 456f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4577112b564SBarry Smith 4587112b564SBarry Smith Level: intermediate 4597112b564SBarry Smith 460562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `MATSOLVERSUPERLU`, 46142747ad1SJacob Faibussowitsch `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4627112b564SBarry Smith @*/ 463d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) 464d71ae5a4SJacob Faibussowitsch { 4655f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(PC, MatSolverType *); 4667112b564SBarry Smith 4677112b564SBarry Smith PetscFunctionBegin; 4680700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4694f572ea9SToby Isaac PetscAssertPointer(stype, 2); 4709566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 4719566063dSJacob Faibussowitsch if (f) PetscCall((*f)(pc, stype)); 4725f80ce2aSJacob Faibussowitsch else *stype = NULL; 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4747112b564SBarry Smith } 4757112b564SBarry Smith 47685317021SBarry Smith /*@ 47785317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 47885317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 47985317021SBarry Smith 480c5eb9154SBarry Smith Not Collective, each process can expect a different amount of fill 48185317021SBarry Smith 48285317021SBarry Smith Input Parameters: 48385317021SBarry Smith + pc - the preconditioner context 48485317021SBarry Smith - fill - amount of expected fill 48585317021SBarry Smith 48685317021SBarry Smith Options Database Key: 48785317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 48885317021SBarry Smith 48985317021SBarry Smith Level: intermediate 49085317021SBarry Smith 491f1580f4eSBarry Smith Notes: 49285317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 49385317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 49485317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 49585317021SBarry Smith future runs. Default PETSc uses a value of 5.0 49685317021SBarry Smith 497f1580f4eSBarry Smith This is ignored for most solver packages 49801a79839SBarry Smith 499f1580f4eSBarry Smith This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels. 500f1580f4eSBarry Smith 501562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()` 50285317021SBarry Smith @*/ 503d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) 504d71ae5a4SJacob Faibussowitsch { 50585317021SBarry Smith PetscFunctionBegin; 5060700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 507467446fbSPierre Jolivet PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less than 1.0"); 508cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51085317021SBarry Smith } 51185317021SBarry Smith 51285317021SBarry Smith /*@ 51304c3f3b8SBarry Smith PCFactorSetUseInPlace - Tells the preconditioner to do an in-place factorization. 51485317021SBarry Smith 515c3339decSBarry Smith Logically Collective 51685317021SBarry Smith 51785317021SBarry Smith Input Parameters: 5188e37d05fSBarry Smith + pc - the preconditioner context 519f1580f4eSBarry Smith - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 52085317021SBarry Smith 52185317021SBarry Smith Options Database Key: 5228e37d05fSBarry Smith . -pc_factor_in_place <true,false> - Activate/deactivate in-place factorization 52385317021SBarry Smith 524f1580f4eSBarry Smith Note: 52504c3f3b8SBarry Smith For dense matrices, this enables the solution of much larger problems. 52604c3f3b8SBarry Smith For sparse matrices the factorization cannot be done truly in-place 52704c3f3b8SBarry Smith so this does not save memory during the factorization, but after the matrix 52804c3f3b8SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 52904c3f3b8SBarry Smith space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 53004c3f3b8SBarry Smith 531f1580f4eSBarry Smith `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when 53285317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 533f1580f4eSBarry Smith a call to `KSPSetOperators()`. 53485317021SBarry Smith This is because the Krylov space methods require an application of the 53585317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 53685317021SBarry Smith been factored in-place, replacing the original matrix. 53785317021SBarry Smith 53885317021SBarry Smith Level: intermediate 53985317021SBarry Smith 540562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `Mat`, `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()` 54185317021SBarry Smith @*/ 542d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) 543d71ae5a4SJacob Faibussowitsch { 54485317021SBarry Smith PetscFunctionBegin; 5450700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 546cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 5473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5488e37d05fSBarry Smith } 5498e37d05fSBarry Smith 5508e37d05fSBarry Smith /*@ 5518e37d05fSBarry Smith PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 5528e37d05fSBarry Smith 553c3339decSBarry Smith Logically Collective 5548e37d05fSBarry Smith 5558e37d05fSBarry Smith Input Parameter: 5568e37d05fSBarry Smith . pc - the preconditioner context 5578e37d05fSBarry Smith 5588e37d05fSBarry Smith Output Parameter: 559f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 5608e37d05fSBarry Smith 5618e37d05fSBarry Smith Level: intermediate 5628e37d05fSBarry Smith 563562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()` 5648e37d05fSBarry Smith @*/ 565d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) 566d71ae5a4SJacob Faibussowitsch { 5678e37d05fSBarry Smith PetscFunctionBegin; 5688e37d05fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 569cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 5703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57185317021SBarry Smith } 57285317021SBarry Smith 5735d83a8b1SBarry Smith /*@ 57485317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 575f1580f4eSBarry Smith be used in the `PCLU`, `PCCHOLESKY`, `PCILU`, or `PCICC` preconditioners 57685317021SBarry Smith 577c3339decSBarry Smith Logically Collective 57885317021SBarry Smith 57985317021SBarry Smith Input Parameters: 58085317021SBarry Smith + pc - the preconditioner context 581f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM` 58285317021SBarry Smith 58385317021SBarry Smith Options Database Key: 5842c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 58585317021SBarry Smith 58685317021SBarry Smith Level: intermediate 58785317021SBarry Smith 58895452b02SPatrick Sanan Notes: 5894ac6704cSBarry Smith Nested dissection is used by default for some of PETSc's sparse matrix formats 59085317021SBarry Smith 591f1580f4eSBarry Smith For `PCCHOLESKY` and `PCICC` and the `MATSBAIJ` format the only reordering available is natural since only the upper half of the matrix is stored 5929bd791bbSBarry Smith and reordering this matrix is very expensive. 5939bd791bbSBarry Smith 594f1580f4eSBarry Smith You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering. 5959bd791bbSBarry Smith 596f1580f4eSBarry Smith `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others. 5972c7c0729SBarry Smith 598562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM` 59985317021SBarry Smith @*/ 600d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) 601d71ae5a4SJacob Faibussowitsch { 60285317021SBarry Smith PetscFunctionBegin; 603c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 604cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 6053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60685317021SBarry Smith } 60785317021SBarry Smith 60885317021SBarry Smith /*@ 6098ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 61085317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 611f1580f4eSBarry Smith it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used. 61285317021SBarry Smith 613c3339decSBarry Smith Logically Collective 61485317021SBarry Smith 61585317021SBarry Smith Input Parameters: 61685317021SBarry Smith + pc - the preconditioner context 61785317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 61885317021SBarry Smith 61985317021SBarry Smith Options Database Key: 620147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 62185317021SBarry Smith 62285317021SBarry Smith Level: intermediate 62385317021SBarry Smith 624562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 62585317021SBarry Smith @*/ 626d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) 627d71ae5a4SJacob Faibussowitsch { 62885317021SBarry Smith PetscFunctionBegin; 629c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 630c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, dtcol, 2); 631cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 6323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63385317021SBarry Smith } 63485317021SBarry Smith 63585317021SBarry Smith /*@ 63685317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 637f1580f4eSBarry Smith with `MATBAIJ` or `MATSBAIJ` matrices 63885317021SBarry Smith 639c3339decSBarry Smith Logically Collective 64085317021SBarry Smith 64185317021SBarry Smith Input Parameters: 64285317021SBarry Smith + pc - the preconditioner context 643f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE` 64485317021SBarry Smith 64585317021SBarry Smith Options Database Key: 646f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ` 64785317021SBarry Smith 64885317021SBarry Smith Level: intermediate 64985317021SBarry Smith 650562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 65185317021SBarry Smith @*/ 652d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) 653d71ae5a4SJacob Faibussowitsch { 65485317021SBarry Smith PetscFunctionBegin; 655c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 656acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, pivot, 2); 657cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 6583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65985317021SBarry Smith } 66085317021SBarry Smith 66185317021SBarry Smith /*@ 662288e7d53SBarry Smith PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 66385317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 66485317021SBarry Smith 665c3339decSBarry Smith Logically Collective 66685317021SBarry Smith 66785317021SBarry Smith Input Parameters: 66885317021SBarry Smith + pc - the preconditioner context 669f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 67085317021SBarry Smith 67185317021SBarry Smith Options Database Key: 672f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 67385317021SBarry Smith 67485317021SBarry Smith Level: intermediate 67585317021SBarry Smith 676562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()` 67785317021SBarry Smith @*/ 678d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) 679d71ae5a4SJacob Faibussowitsch { 68085317021SBarry Smith PetscFunctionBegin; 681064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 682acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 683cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 6843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68585317021SBarry Smith } 6863d1c1ea0SBarry Smith 687d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) 688d71ae5a4SJacob Faibussowitsch { 6893d1c1ea0SBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 6903d1c1ea0SBarry Smith 6913d1c1ea0SBarry Smith PetscFunctionBegin; 6929566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&fact->info)); 6934ac6704cSBarry Smith fact->factortype = ftype; 6943d1c1ea0SBarry Smith fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 6953d1c1ea0SBarry Smith fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 6963d1c1ea0SBarry Smith fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 6973d1c1ea0SBarry Smith fact->info.pivotinblocks = 1.0; 6983d1c1ea0SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 6993d1c1ea0SBarry Smith 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 7049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 7069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 7129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 7139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 7149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 7159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 7169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 7179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 7199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 7203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7213d1c1ea0SBarry Smith } 7222e956fe4SStefano Zampini 723d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc) 724d71ae5a4SJacob Faibussowitsch { 7252e956fe4SStefano Zampini PetscFunctionBegin; 7262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 7272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 7282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 7292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 7302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 7312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 7322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 7332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 7342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 7352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 7362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 7372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 7382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 7392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 7402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 7412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 7422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 7432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 7442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 7452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 7462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 7472e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 7483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7492e956fe4SStefano Zampini } 750