15e8efad8SHong Zhang 2c6db04a5SJed Brown #include <../src/ksp/pc/impls/factor/factor.h> /*I "petscpc.h" I*/ 3978beb7fSPierre Jolivet #include <petsc/private/matimpl.h> 45e8efad8SHong Zhang 54ac6704cSBarry Smith /* 64ac6704cSBarry Smith If an ordering is not yet set and the matrix is available determine a default ordering 74ac6704cSBarry Smith */ 8d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc) 9d71ae5a4SJacob Faibussowitsch { 10978beb7fSPierre Jolivet PetscBool foundmtype, flg; 1126cc229bSBarry Smith const char *prefix; 124ac6704cSBarry Smith 134ac6704cSBarry Smith PetscFunctionBegin; 144ac6704cSBarry Smith if (pc->pmat) { 1526cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc, &prefix)); 1626cc229bSBarry Smith PetscCall(MatSetOptionsPrefixFactor(pc->pmat, prefix)); 174ac6704cSBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 189566063dSJacob Faibussowitsch PetscCall(MatSolverTypeGet(fact->solvertype, ((PetscObject)pc->pmat)->type_name, fact->factortype, NULL, &foundmtype, NULL)); 19978beb7fSPierre Jolivet if (foundmtype) { 204ac6704cSBarry Smith if (!fact->fact) { 219566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); 22978beb7fSPierre Jolivet } else if (!fact->fact->assembled) { 239566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(fact->solvertype, fact->fact->solvertype, &flg)); 24978beb7fSPierre Jolivet if (!flg) { 2526cc229bSBarry Smith Mat B; 269566063dSJacob Faibussowitsch PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &B)); 279566063dSJacob Faibussowitsch PetscCall(MatHeaderReplace(fact->fact, &B)); 28978beb7fSPierre Jolivet } 294ac6704cSBarry Smith } 304ac6704cSBarry Smith if (!fact->ordering) { 314ac6704cSBarry Smith PetscBool canuseordering; 324ac6704cSBarry Smith MatOrderingType otype; 334ac6704cSBarry Smith 349566063dSJacob Faibussowitsch PetscCall(MatFactorGetCanUseOrdering(fact->fact, &canuseordering)); 354ac6704cSBarry Smith if (canuseordering) { 369566063dSJacob Faibussowitsch PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype)); 374ac6704cSBarry Smith } else otype = MATORDERINGEXTERNAL; 389566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering)); 394ac6704cSBarry Smith } 404ac6704cSBarry Smith } 41978beb7fSPierre Jolivet } 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 434ac6704cSBarry Smith } 444ac6704cSBarry Smith 45d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag) 46d71ae5a4SJacob Faibussowitsch { 473d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 483d1c1ea0SBarry Smith 493d1c1ea0SBarry Smith PetscFunctionBegin; 503d1c1ea0SBarry Smith lu->reuseordering = flag; 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 523d1c1ea0SBarry Smith } 533d1c1ea0SBarry Smith 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag) 55d71ae5a4SJacob Faibussowitsch { 563d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 573d1c1ea0SBarry Smith 583d1c1ea0SBarry Smith PetscFunctionBegin; 593d1c1ea0SBarry Smith lu->reusefill = flag; 603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 613d1c1ea0SBarry Smith } 623d1c1ea0SBarry Smith 63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg) 64d71ae5a4SJacob Faibussowitsch { 653d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 663d1c1ea0SBarry Smith 673d1c1ea0SBarry Smith PetscFunctionBegin; 683d1c1ea0SBarry Smith dir->inplace = flg; 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 703d1c1ea0SBarry Smith } 713d1c1ea0SBarry Smith 72d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg) 73d71ae5a4SJacob Faibussowitsch { 743d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 753d1c1ea0SBarry Smith 763d1c1ea0SBarry Smith PetscFunctionBegin; 773d1c1ea0SBarry Smith *flg = dir->inplace; 783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 793d1c1ea0SBarry Smith } 803d1c1ea0SBarry Smith 81f8260c8fSBarry Smith /*@ 82f1580f4eSBarry Smith PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may 83f8260c8fSBarry Smith set the options for that particular factorization object. 84f8260c8fSBarry Smith 85f8260c8fSBarry Smith Input Parameter: 86f8260c8fSBarry Smith . pc - the preconditioner context 87f8260c8fSBarry Smith 88f1580f4eSBarry Smith Note: 89f1580f4eSBarry 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. 90f8260c8fSBarry Smith 912bd2b0e6SSatish Balay Level: intermediate 922bd2b0e6SSatish Balay 93f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 94f8260c8fSBarry Smith @*/ 95d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 96d71ae5a4SJacob Faibussowitsch { 97f8260c8fSBarry Smith PetscFunctionBegin; 98f8260c8fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 99cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101f8260c8fSBarry Smith } 102f8260c8fSBarry Smith 103ee45ca4aSHong Zhang /*@ 104ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 105ee45ca4aSHong Zhang 106c3339decSBarry Smith Logically Collective 107ee45ca4aSHong Zhang 108ee45ca4aSHong Zhang Input Parameters: 109afaefe49SHong Zhang + pc - the preconditioner context 110afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 111ee45ca4aSHong Zhang 112ee45ca4aSHong Zhang Options Database Key: 113ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 114ee45ca4aSHong Zhang 115ee45ca4aSHong Zhang Level: intermediate 116ee45ca4aSHong Zhang 117f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 118ee45ca4aSHong Zhang @*/ 119d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) 120d71ae5a4SJacob Faibussowitsch { 121ee45ca4aSHong Zhang PetscFunctionBegin; 1220700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 123c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, zero, 2); 124cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126ee45ca4aSHong Zhang } 127ee45ca4aSHong Zhang 128915743fcSHong Zhang /*@ 129915743fcSHong Zhang PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 130915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 131915743fcSHong Zhang 132c3339decSBarry Smith Logically Collective 133915743fcSHong Zhang 134915743fcSHong Zhang Input Parameters: 135915743fcSHong Zhang + pc - the preconditioner context 136f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS` 137915743fcSHong Zhang 138915743fcSHong Zhang Options Database Key: 13928d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 140915743fcSHong Zhang 141915743fcSHong Zhang Level: intermediate 142915743fcSHong Zhang 143f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 144915743fcSHong Zhang @*/ 145d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) 146d71ae5a4SJacob Faibussowitsch { 147d90ac83dSHong Zhang PetscFunctionBegin; 1480700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 149c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 150cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152d90ac83dSHong Zhang } 153d90ac83dSHong Zhang 154915743fcSHong Zhang /*@ 155915743fcSHong Zhang PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 156915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 157915743fcSHong Zhang 158c3339decSBarry Smith Logically Collective 159915743fcSHong Zhang 160915743fcSHong Zhang Input Parameters: 161915743fcSHong Zhang + pc - the preconditioner context 162f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default 163915743fcSHong Zhang 164915743fcSHong Zhang Options Database Key: 165f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 166915743fcSHong Zhang 167915743fcSHong Zhang Level: intermediate 168915743fcSHong Zhang 169f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, ``PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 170915743fcSHong Zhang @*/ 171d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) 172d71ae5a4SJacob Faibussowitsch { 173d90ac83dSHong Zhang PetscFunctionBegin; 1740700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 175c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 176cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178d90ac83dSHong Zhang } 179d90ac83dSHong Zhang 1806d33885cSprj- /*@ 181f1580f4eSBarry Smith PCFactorSetDropTolerance - The preconditioner will use an `PCILU` 182f1580f4eSBarry Smith based on a drop tolerance. 18385317021SBarry Smith 184c3339decSBarry Smith Logically Collective 18585317021SBarry Smith 18685317021SBarry Smith Input Parameters: 18785317021SBarry Smith + pc - the preconditioner context 18885317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 18985317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 19085317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 19185317021SBarry Smith depends on the number of nonzeros in row of original matrix 19285317021SBarry Smith 19385317021SBarry Smith Options Database Key: 194b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 19585317021SBarry Smith 19685317021SBarry Smith Level: intermediate 19785317021SBarry Smith 198f1580f4eSBarry Smith Note: 19985317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 20085317021SBarry Smith matrix. We don't know how to compute reasonable values. 20185317021SBarry Smith 202f1580f4eSBarry Smith .seealso: `PCILU` 2036d33885cSprj- @*/ 204d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) 205d71ae5a4SJacob Faibussowitsch { 20685317021SBarry Smith PetscFunctionBegin; 2070700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 208064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveReal(pc, dtcol, 3); 209064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 210cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21285317021SBarry Smith } 21385317021SBarry Smith 214c7f610a1SBarry Smith /*@ 215c7f610a1SBarry Smith PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 216c7f610a1SBarry Smith 217c7f610a1SBarry Smith Not Collective 218c7f610a1SBarry Smith 219*2fe279fdSBarry Smith Input Parameter: 220c7f610a1SBarry Smith . pc - the preconditioner context 221c7f610a1SBarry Smith 222c7f610a1SBarry Smith Output Parameter: 223c7f610a1SBarry Smith . pivot - the tolerance 224c7f610a1SBarry Smith 225c7f610a1SBarry Smith Level: intermediate 226c7f610a1SBarry Smith 227f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 228c7f610a1SBarry Smith @*/ 229d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) 230d71ae5a4SJacob Faibussowitsch { 231c7f610a1SBarry Smith PetscFunctionBegin; 232c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 233cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 2343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 235c7f610a1SBarry Smith } 236c7f610a1SBarry Smith 237c7f610a1SBarry Smith /*@ 238c7f610a1SBarry Smith PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 239c7f610a1SBarry Smith 240c7f610a1SBarry Smith Not Collective 241c7f610a1SBarry Smith 242*2fe279fdSBarry Smith Input Parameter: 243c7f610a1SBarry Smith . pc - the preconditioner context 244c7f610a1SBarry Smith 245c7f610a1SBarry Smith Output Parameter: 246c7f610a1SBarry Smith . shift - how much to shift the diagonal entry 247c7f610a1SBarry Smith 248c7f610a1SBarry Smith Level: intermediate 249c7f610a1SBarry Smith 250f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 251c7f610a1SBarry Smith @*/ 252d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) 253d71ae5a4SJacob Faibussowitsch { 254c7f610a1SBarry Smith PetscFunctionBegin; 255c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 256cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 258c7f610a1SBarry Smith } 259c7f610a1SBarry Smith 260c7f610a1SBarry Smith /*@ 261c7f610a1SBarry Smith PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 262c7f610a1SBarry Smith 263c7f610a1SBarry Smith Not Collective 264c7f610a1SBarry Smith 265f1580f4eSBarry Smith Input Parameter: 266c7f610a1SBarry Smith . pc - the preconditioner context 267c7f610a1SBarry Smith 268c7f610a1SBarry Smith Output Parameter: 269f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 270c7f610a1SBarry Smith 271c7f610a1SBarry Smith Level: intermediate 272c7f610a1SBarry Smith 273f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 274c7f610a1SBarry Smith @*/ 275d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) 276d71ae5a4SJacob Faibussowitsch { 277c7f610a1SBarry Smith PetscFunctionBegin; 278c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 279cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 2803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 281c7f610a1SBarry Smith } 282c7f610a1SBarry Smith 2832591b318SToby Isaac /*@ 2842591b318SToby Isaac PCFactorGetLevels - Gets the number of levels of fill to use. 2852591b318SToby Isaac 286c3339decSBarry Smith Logically Collective 2872591b318SToby Isaac 288*2fe279fdSBarry Smith Input Parameter: 2892591b318SToby Isaac . pc - the preconditioner context 2902591b318SToby Isaac 2912591b318SToby Isaac Output Parameter: 2922591b318SToby Isaac . levels - number of levels of fill 2932591b318SToby Isaac 2942591b318SToby Isaac Level: intermediate 2952591b318SToby Isaac 296f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetLevels()` 2972591b318SToby Isaac @*/ 298d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) 299d71ae5a4SJacob Faibussowitsch { 3002591b318SToby Isaac PetscFunctionBegin; 3012591b318SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 302cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3042591b318SToby Isaac } 3052591b318SToby Isaac 30685317021SBarry Smith /*@ 30785317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 30885317021SBarry Smith 309c3339decSBarry Smith Logically Collective 31085317021SBarry Smith 31185317021SBarry Smith Input Parameters: 31285317021SBarry Smith + pc - the preconditioner context 31385317021SBarry Smith - levels - number of levels of fill 31485317021SBarry Smith 31585317021SBarry Smith Options Database Key: 31685317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 31785317021SBarry Smith 31885317021SBarry Smith Level: intermediate 31985317021SBarry Smith 320f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetLevels()` 32185317021SBarry Smith @*/ 322d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) 323d71ae5a4SJacob Faibussowitsch { 32485317021SBarry Smith PetscFunctionBegin; 3250700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3265f80ce2aSJacob Faibussowitsch PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 327c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, levels, 2); 328cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 3293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33085317021SBarry Smith } 33185317021SBarry Smith 33285317021SBarry Smith /*@ 33385317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 33485317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 33585317021SBarry Smith 336c3339decSBarry Smith Logically Collective 33785317021SBarry Smith 33885317021SBarry Smith Input Parameters: 33985317021SBarry Smith + pc - the preconditioner context 340f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 34185317021SBarry Smith 34285317021SBarry Smith Options Database Key: 34367b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 34485317021SBarry Smith 345f1580f4eSBarry Smith Note: 34685317021SBarry Smith Does not apply with 0 fill. 34785317021SBarry Smith 34885317021SBarry Smith Level: intermediate 34985317021SBarry Smith 350f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 35185317021SBarry Smith @*/ 352d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) 353d71ae5a4SJacob Faibussowitsch { 35485317021SBarry Smith PetscFunctionBegin; 3550700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 356cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35892e9c092SBarry Smith } 35992e9c092SBarry Smith 36092e9c092SBarry Smith /*@ 36192e9c092SBarry Smith PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 36292e9c092SBarry Smith treated as level 0 fill even if there is no non-zero location. 36392e9c092SBarry Smith 364c3339decSBarry Smith Logically Collective 36592e9c092SBarry Smith 36692e9c092SBarry Smith Input Parameter: 36792e9c092SBarry Smith . pc - the preconditioner context 36892e9c092SBarry Smith 36992e9c092SBarry Smith Output Parameter: 370f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 37192e9c092SBarry Smith 372f1580f4eSBarry Smith Note: 37392e9c092SBarry Smith Does not apply with 0 fill. 37492e9c092SBarry Smith 37592e9c092SBarry Smith Level: intermediate 37692e9c092SBarry Smith 377f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 37892e9c092SBarry Smith @*/ 379d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) 380d71ae5a4SJacob Faibussowitsch { 38192e9c092SBarry Smith PetscFunctionBegin; 38292e9c092SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 383cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 3843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38585317021SBarry Smith } 38685317021SBarry Smith 38785317021SBarry Smith /*@ 38885317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 38985317021SBarry Smith 390c3339decSBarry Smith Logically Collective 39185317021SBarry Smith 39285317021SBarry Smith Input Parameters: 39385317021SBarry Smith + pc - the preconditioner context 39485317021SBarry Smith - tol - diagonal entries smaller than this in absolute value are considered zero 39585317021SBarry Smith 39685317021SBarry Smith Options Database Key: 397147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 39885317021SBarry Smith 39985317021SBarry Smith Level: intermediate 40085317021SBarry Smith 401f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 40285317021SBarry Smith @*/ 403d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) 404d71ae5a4SJacob Faibussowitsch { 40585317021SBarry Smith PetscFunctionBegin; 4060700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 407c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, rtol, 2); 408cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 4093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41085317021SBarry Smith } 41185317021SBarry Smith 412bf6011e8SBarry Smith /*@C 413f1580f4eSBarry Smith PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 41485317021SBarry Smith 415c3339decSBarry Smith Logically Collective 41685317021SBarry Smith 41785317021SBarry Smith Input Parameters: 41885317021SBarry Smith + pc - the preconditioner context 419f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 42085317021SBarry Smith 42185317021SBarry Smith Options Database Key: 4223ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 42385317021SBarry Smith 42485317021SBarry Smith Level: intermediate 42585317021SBarry Smith 42685317021SBarry Smith Note: 42785317021SBarry Smith By default this will use the PETSc factorization if it exists 42885317021SBarry Smith 429f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 430f1580f4eSBarry Smith `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 43185317021SBarry Smith @*/ 432d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) 433d71ae5a4SJacob Faibussowitsch { 43485317021SBarry Smith PetscFunctionBegin; 4350700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 436cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 4373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43885317021SBarry Smith } 43985317021SBarry Smith 440bf6011e8SBarry Smith /*@C 441f1580f4eSBarry Smith PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization 4427112b564SBarry Smith 443c5eb9154SBarry Smith Not Collective 4447112b564SBarry Smith 4457112b564SBarry Smith Input Parameter: 4467112b564SBarry Smith . pc - the preconditioner context 4477112b564SBarry Smith 4487112b564SBarry Smith Output Parameter: 449f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4507112b564SBarry Smith 4517112b564SBarry Smith Level: intermediate 4527112b564SBarry Smith 453f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 454f1580f4eSBarry Smith `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4557112b564SBarry Smith @*/ 456d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) 457d71ae5a4SJacob Faibussowitsch { 4585f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(PC, MatSolverType *); 4597112b564SBarry Smith 4607112b564SBarry Smith PetscFunctionBegin; 4610700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4625f80ce2aSJacob Faibussowitsch PetscValidPointer(stype, 2); 4639566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 4649566063dSJacob Faibussowitsch if (f) PetscCall((*f)(pc, stype)); 4655f80ce2aSJacob Faibussowitsch else *stype = NULL; 4663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4677112b564SBarry Smith } 4687112b564SBarry Smith 46985317021SBarry Smith /*@ 47085317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 47185317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 47285317021SBarry Smith 473c5eb9154SBarry Smith Not Collective, each process can expect a different amount of fill 47485317021SBarry Smith 47585317021SBarry Smith Input Parameters: 47685317021SBarry Smith + pc - the preconditioner context 47785317021SBarry Smith - fill - amount of expected fill 47885317021SBarry Smith 47985317021SBarry Smith Options Database Key: 48085317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 48185317021SBarry Smith 48285317021SBarry Smith Level: intermediate 48385317021SBarry Smith 484f1580f4eSBarry Smith Notes: 48585317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 48685317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 48785317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 48885317021SBarry Smith future runs. Default PETSc uses a value of 5.0 48985317021SBarry Smith 490f1580f4eSBarry Smith This is ignored for most solver packages 49101a79839SBarry Smith 492f1580f4eSBarry Smith This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels. 493f1580f4eSBarry Smith 494f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()` 49585317021SBarry Smith @*/ 496d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) 497d71ae5a4SJacob Faibussowitsch { 49885317021SBarry Smith PetscFunctionBegin; 4990700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5005f80ce2aSJacob Faibussowitsch PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0"); 501cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 5023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50385317021SBarry Smith } 50485317021SBarry Smith 50585317021SBarry Smith /*@ 50685317021SBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 50785317021SBarry Smith For dense matrices, this enables the solution of much larger problems. 50885317021SBarry Smith For sparse matrices the factorization cannot be done truly in-place 50985317021SBarry Smith so this does not save memory during the factorization, but after the matrix 51085317021SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 511ec5066bdSBarry Smith space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 51285317021SBarry Smith 513c3339decSBarry Smith Logically Collective 51485317021SBarry Smith 51585317021SBarry Smith Input Parameters: 5168e37d05fSBarry Smith + pc - the preconditioner context 517f1580f4eSBarry Smith - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 51885317021SBarry Smith 51985317021SBarry Smith Options Database Key: 5208e37d05fSBarry Smith . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization 52185317021SBarry Smith 522f1580f4eSBarry Smith Note: 523f1580f4eSBarry Smith `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when 52485317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 525f1580f4eSBarry Smith a call to `KSPSetOperators()`. 52685317021SBarry Smith This is because the Krylov space methods require an application of the 52785317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 52885317021SBarry Smith been factored in-place, replacing the original matrix. 52985317021SBarry Smith 53085317021SBarry Smith Level: intermediate 53185317021SBarry Smith 532f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()` 53385317021SBarry Smith @*/ 534d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) 535d71ae5a4SJacob Faibussowitsch { 53685317021SBarry Smith PetscFunctionBegin; 5370700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 538cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 5393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5408e37d05fSBarry Smith } 5418e37d05fSBarry Smith 5428e37d05fSBarry Smith /*@ 5438e37d05fSBarry Smith PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 5448e37d05fSBarry Smith 545c3339decSBarry Smith Logically Collective 5468e37d05fSBarry Smith 5478e37d05fSBarry Smith Input Parameter: 5488e37d05fSBarry Smith . pc - the preconditioner context 5498e37d05fSBarry Smith 5508e37d05fSBarry Smith Output Parameter: 551f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 5528e37d05fSBarry Smith 5538e37d05fSBarry Smith Level: intermediate 5548e37d05fSBarry Smith 555f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()` 5568e37d05fSBarry Smith @*/ 557d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) 558d71ae5a4SJacob Faibussowitsch { 5598e37d05fSBarry Smith PetscFunctionBegin; 5608e37d05fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 561cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56385317021SBarry Smith } 56485317021SBarry Smith 56585317021SBarry Smith /*@C 56685317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 567f1580f4eSBarry Smith be used in the `PCLU`, `PCCHOLESKY`, `PCILU`, or `PCICC` preconditioners 56885317021SBarry Smith 569c3339decSBarry Smith Logically Collective 57085317021SBarry Smith 57185317021SBarry Smith Input Parameters: 57285317021SBarry Smith + pc - the preconditioner context 573f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM` 57485317021SBarry Smith 57585317021SBarry Smith Options Database Key: 5762c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 57785317021SBarry Smith 57885317021SBarry Smith Level: intermediate 57985317021SBarry Smith 58095452b02SPatrick Sanan Notes: 5814ac6704cSBarry Smith Nested dissection is used by default for some of PETSc's sparse matrix formats 58285317021SBarry Smith 583f1580f4eSBarry 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 5849bd791bbSBarry Smith and reordering this matrix is very expensive. 5859bd791bbSBarry Smith 586f1580f4eSBarry Smith You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering. 5879bd791bbSBarry Smith 588f1580f4eSBarry Smith `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others. 5892c7c0729SBarry Smith 590f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM` 59185317021SBarry Smith @*/ 592d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) 593d71ae5a4SJacob Faibussowitsch { 59485317021SBarry Smith PetscFunctionBegin; 595c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 596cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 5973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59885317021SBarry Smith } 59985317021SBarry Smith 60085317021SBarry Smith /*@ 6018ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 60285317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 603f1580f4eSBarry Smith it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used. 60485317021SBarry Smith 605c3339decSBarry Smith Logically Collective 60685317021SBarry Smith 60785317021SBarry Smith Input Parameters: 60885317021SBarry Smith + pc - the preconditioner context 60985317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 61085317021SBarry Smith 61185317021SBarry Smith Options Database Key: 612147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 61385317021SBarry Smith 61485317021SBarry Smith Level: intermediate 61585317021SBarry Smith 616f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 61785317021SBarry Smith @*/ 618d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) 619d71ae5a4SJacob Faibussowitsch { 62085317021SBarry Smith PetscFunctionBegin; 621c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 622c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, dtcol, 2); 623cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 6243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62585317021SBarry Smith } 62685317021SBarry Smith 62785317021SBarry Smith /*@ 62885317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 629f1580f4eSBarry Smith with `MATBAIJ` or `MATSBAIJ` matrices 63085317021SBarry Smith 631c3339decSBarry Smith Logically Collective 63285317021SBarry Smith 63385317021SBarry Smith Input Parameters: 63485317021SBarry Smith + pc - the preconditioner context 635f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE` 63685317021SBarry Smith 63785317021SBarry Smith Options Database Key: 638f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ` 63985317021SBarry Smith 64085317021SBarry Smith Level: intermediate 64185317021SBarry Smith 642f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 64385317021SBarry Smith @*/ 644d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) 645d71ae5a4SJacob Faibussowitsch { 64685317021SBarry Smith PetscFunctionBegin; 647c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 648acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, pivot, 2); 649cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 6503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65185317021SBarry Smith } 65285317021SBarry Smith 65385317021SBarry Smith /*@ 654288e7d53SBarry Smith PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 65585317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 65685317021SBarry Smith 657c3339decSBarry Smith Logically Collective 65885317021SBarry Smith 65985317021SBarry Smith Input Parameters: 66085317021SBarry Smith + pc - the preconditioner context 661f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 66285317021SBarry Smith 66385317021SBarry Smith Options Database Key: 664f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 66585317021SBarry Smith 66685317021SBarry Smith Level: intermediate 66785317021SBarry Smith 668f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()` 66985317021SBarry Smith @*/ 670d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) 671d71ae5a4SJacob Faibussowitsch { 67285317021SBarry Smith PetscFunctionBegin; 673064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 674acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 675cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 6763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67785317021SBarry Smith } 6783d1c1ea0SBarry Smith 679d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) 680d71ae5a4SJacob Faibussowitsch { 6813d1c1ea0SBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 6823d1c1ea0SBarry Smith 6833d1c1ea0SBarry Smith PetscFunctionBegin; 6849566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&fact->info)); 6854ac6704cSBarry Smith fact->factortype = ftype; 6863d1c1ea0SBarry Smith fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 6873d1c1ea0SBarry Smith fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 6883d1c1ea0SBarry Smith fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 6893d1c1ea0SBarry Smith fact->info.pivotinblocks = 1.0; 6903d1c1ea0SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 6913d1c1ea0SBarry Smith 6929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 6949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 6959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 6969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 6979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 6989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 6999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 7049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 7069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 7123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7133d1c1ea0SBarry Smith } 7142e956fe4SStefano Zampini 715d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc) 716d71ae5a4SJacob Faibussowitsch { 7172e956fe4SStefano Zampini PetscFunctionBegin; 7182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 7192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 7202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 7212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 7222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 7232e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 7242e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 7252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 7262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 7272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 7282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 7292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 7302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 7312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 7322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 7332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 7342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 7352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 7362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 7372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 7382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 7392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 7403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7412e956fe4SStefano Zampini } 742