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)); 394ac6704cSBarry Smith if (canuseordering) { 409566063dSJacob Faibussowitsch PetscCall(MatFactorGetPreferredOrdering(fact->fact, fact->factortype, &otype)); 414ac6704cSBarry Smith } else otype = MATORDERINGEXTERNAL; 429566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(otype, (char **)&fact->ordering)); 434ac6704cSBarry Smith } 447def7855SStefano Zampini if (destroy) PetscCall(MatDestroy(&fact->fact)); 454ac6704cSBarry Smith } 46978beb7fSPierre Jolivet } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 484ac6704cSBarry Smith } 494ac6704cSBarry Smith 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc, PetscBool flag) 51d71ae5a4SJacob Faibussowitsch { 523d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 533d1c1ea0SBarry Smith 543d1c1ea0SBarry Smith PetscFunctionBegin; 553d1c1ea0SBarry Smith lu->reuseordering = flag; 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 573d1c1ea0SBarry Smith } 583d1c1ea0SBarry Smith 59d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc, PetscBool flag) 60d71ae5a4SJacob Faibussowitsch { 613d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor *)pc->data; 623d1c1ea0SBarry Smith 633d1c1ea0SBarry Smith PetscFunctionBegin; 643d1c1ea0SBarry Smith lu->reusefill = flag; 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 663d1c1ea0SBarry Smith } 673d1c1ea0SBarry Smith 68d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc, PetscBool flg) 69d71ae5a4SJacob Faibussowitsch { 703d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 713d1c1ea0SBarry Smith 723d1c1ea0SBarry Smith PetscFunctionBegin; 733d1c1ea0SBarry Smith dir->inplace = flg; 743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 753d1c1ea0SBarry Smith } 763d1c1ea0SBarry Smith 77d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc, PetscBool *flg) 78d71ae5a4SJacob Faibussowitsch { 793d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor *)pc->data; 803d1c1ea0SBarry Smith 813d1c1ea0SBarry Smith PetscFunctionBegin; 823d1c1ea0SBarry Smith *flg = dir->inplace; 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 843d1c1ea0SBarry Smith } 853d1c1ea0SBarry Smith 86f8260c8fSBarry Smith /*@ 87f1580f4eSBarry Smith PCFactorSetUpMatSolverType - Can be called after `KSPSetOperators()` or `PCSetOperators()`, causes `MatGetFactor()` to be called so then one may 88f8260c8fSBarry Smith set the options for that particular factorization object. 89f8260c8fSBarry Smith 90f8260c8fSBarry Smith Input Parameter: 91f8260c8fSBarry Smith . pc - the preconditioner context 92f8260c8fSBarry Smith 93f1580f4eSBarry Smith Note: 94f1580f4eSBarry 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. 9503e5aca4SStefano Zampini This function raises an error if the requested combination of solver package and matrix type is not supported. 96f8260c8fSBarry Smith 972bd2b0e6SSatish Balay Level: intermediate 982bd2b0e6SSatish Balay 99*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 100f8260c8fSBarry Smith @*/ 101d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 102d71ae5a4SJacob Faibussowitsch { 103f8260c8fSBarry Smith PetscFunctionBegin; 104f8260c8fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 105cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107f8260c8fSBarry Smith } 108f8260c8fSBarry Smith 109ee45ca4aSHong Zhang /*@ 110ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 111ee45ca4aSHong Zhang 112c3339decSBarry Smith Logically Collective 113ee45ca4aSHong Zhang 114ee45ca4aSHong Zhang Input Parameters: 115afaefe49SHong Zhang + pc - the preconditioner context 116afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 117ee45ca4aSHong Zhang 118ee45ca4aSHong Zhang Options Database Key: 119ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 120ee45ca4aSHong Zhang 121ee45ca4aSHong Zhang Level: intermediate 122ee45ca4aSHong Zhang 123*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 124ee45ca4aSHong Zhang @*/ 125d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) 126d71ae5a4SJacob Faibussowitsch { 127ee45ca4aSHong Zhang PetscFunctionBegin; 1280700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 129c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, zero, 2); 130cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 132ee45ca4aSHong Zhang } 133ee45ca4aSHong Zhang 134915743fcSHong Zhang /*@ 135915743fcSHong Zhang PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 136915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 137915743fcSHong Zhang 138c3339decSBarry Smith Logically Collective 139915743fcSHong Zhang 140915743fcSHong Zhang Input Parameters: 141915743fcSHong Zhang + pc - the preconditioner context 142f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS` 143915743fcSHong Zhang 144915743fcSHong Zhang Options Database Key: 14528d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 146915743fcSHong Zhang 147915743fcSHong Zhang Level: intermediate 148915743fcSHong Zhang 149*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 150915743fcSHong Zhang @*/ 151d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) 152d71ae5a4SJacob Faibussowitsch { 153d90ac83dSHong Zhang PetscFunctionBegin; 1540700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 155c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 156cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158d90ac83dSHong Zhang } 159d90ac83dSHong Zhang 160915743fcSHong Zhang /*@ 161915743fcSHong Zhang PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 162915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 163915743fcSHong Zhang 164c3339decSBarry Smith Logically Collective 165915743fcSHong Zhang 166915743fcSHong Zhang Input Parameters: 167915743fcSHong Zhang + pc - the preconditioner context 168f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default 169915743fcSHong Zhang 170915743fcSHong Zhang Options Database Key: 171f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 172915743fcSHong Zhang 173915743fcSHong Zhang Level: intermediate 174915743fcSHong Zhang 175*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCHOLESKY`, `PCLU`, ``PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 176915743fcSHong Zhang @*/ 177d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) 178d71ae5a4SJacob Faibussowitsch { 179d90ac83dSHong Zhang PetscFunctionBegin; 1800700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 181c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 182cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 184d90ac83dSHong Zhang } 185d90ac83dSHong Zhang 1866d33885cSprj- /*@ 187f1580f4eSBarry Smith PCFactorSetDropTolerance - The preconditioner will use an `PCILU` 188f1580f4eSBarry Smith based on a drop tolerance. 18985317021SBarry Smith 190c3339decSBarry Smith Logically Collective 19185317021SBarry Smith 19285317021SBarry Smith Input Parameters: 19385317021SBarry Smith + pc - the preconditioner context 19485317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 19585317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 19685317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 19785317021SBarry Smith depends on the number of nonzeros in row of original matrix 19885317021SBarry Smith 19985317021SBarry Smith Options Database Key: 200b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 20185317021SBarry Smith 20285317021SBarry Smith Level: intermediate 20385317021SBarry Smith 204f1580f4eSBarry Smith Note: 20585317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 20685317021SBarry Smith matrix. We don't know how to compute reasonable values. 20785317021SBarry Smith 208*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU` 2096d33885cSprj- @*/ 210d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) 211d71ae5a4SJacob Faibussowitsch { 21285317021SBarry Smith PetscFunctionBegin; 2130700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 214064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveReal(pc, dtcol, 3); 215064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 216cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21885317021SBarry Smith } 21985317021SBarry Smith 220c7f610a1SBarry Smith /*@ 221c7f610a1SBarry Smith PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 222c7f610a1SBarry Smith 223c7f610a1SBarry Smith Not Collective 224c7f610a1SBarry Smith 2252fe279fdSBarry Smith Input Parameter: 226c7f610a1SBarry Smith . pc - the preconditioner context 227c7f610a1SBarry Smith 228c7f610a1SBarry Smith Output Parameter: 229c7f610a1SBarry Smith . pivot - the tolerance 230c7f610a1SBarry Smith 231c7f610a1SBarry Smith Level: intermediate 232c7f610a1SBarry Smith 233*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 234c7f610a1SBarry Smith @*/ 235d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) 236d71ae5a4SJacob Faibussowitsch { 237c7f610a1SBarry Smith PetscFunctionBegin; 238c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 239cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 241c7f610a1SBarry Smith } 242c7f610a1SBarry Smith 243c7f610a1SBarry Smith /*@ 244c7f610a1SBarry Smith PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 245c7f610a1SBarry Smith 246c7f610a1SBarry Smith Not Collective 247c7f610a1SBarry Smith 2482fe279fdSBarry Smith Input Parameter: 249c7f610a1SBarry Smith . pc - the preconditioner context 250c7f610a1SBarry Smith 251c7f610a1SBarry Smith Output Parameter: 252c7f610a1SBarry Smith . shift - how much to shift the diagonal entry 253c7f610a1SBarry Smith 254c7f610a1SBarry Smith Level: intermediate 255c7f610a1SBarry Smith 256*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 257c7f610a1SBarry Smith @*/ 258d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) 259d71ae5a4SJacob Faibussowitsch { 260c7f610a1SBarry Smith PetscFunctionBegin; 261c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 262cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 264c7f610a1SBarry Smith } 265c7f610a1SBarry Smith 266c7f610a1SBarry Smith /*@ 267c7f610a1SBarry Smith PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 268c7f610a1SBarry Smith 269c7f610a1SBarry Smith Not Collective 270c7f610a1SBarry Smith 271f1580f4eSBarry Smith Input Parameter: 272c7f610a1SBarry Smith . pc - the preconditioner context 273c7f610a1SBarry Smith 274c7f610a1SBarry Smith Output Parameter: 275f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 276c7f610a1SBarry Smith 277c7f610a1SBarry Smith Level: intermediate 278c7f610a1SBarry Smith 279*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 280c7f610a1SBarry Smith @*/ 281d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) 282d71ae5a4SJacob Faibussowitsch { 283c7f610a1SBarry Smith PetscFunctionBegin; 284c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 285cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 287c7f610a1SBarry Smith } 288c7f610a1SBarry Smith 2892591b318SToby Isaac /*@ 2902591b318SToby Isaac PCFactorGetLevels - Gets the number of levels of fill to use. 2912591b318SToby Isaac 292c3339decSBarry Smith Logically Collective 2932591b318SToby Isaac 2942fe279fdSBarry Smith Input Parameter: 2952591b318SToby Isaac . pc - the preconditioner context 2962591b318SToby Isaac 2972591b318SToby Isaac Output Parameter: 2982591b318SToby Isaac . levels - number of levels of fill 2992591b318SToby Isaac 3002591b318SToby Isaac Level: intermediate 3012591b318SToby Isaac 302*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetLevels()` 3032591b318SToby Isaac @*/ 304d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) 305d71ae5a4SJacob Faibussowitsch { 3062591b318SToby Isaac PetscFunctionBegin; 3072591b318SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 308cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3102591b318SToby Isaac } 3112591b318SToby Isaac 31285317021SBarry Smith /*@ 31385317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 31485317021SBarry Smith 315c3339decSBarry Smith Logically Collective 31685317021SBarry Smith 31785317021SBarry Smith Input Parameters: 31885317021SBarry Smith + pc - the preconditioner context 31985317021SBarry Smith - levels - number of levels of fill 32085317021SBarry Smith 32185317021SBarry Smith Options Database Key: 32285317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 32385317021SBarry Smith 32485317021SBarry Smith Level: intermediate 32585317021SBarry Smith 326*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetLevels()` 32785317021SBarry Smith @*/ 328d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) 329d71ae5a4SJacob Faibussowitsch { 33085317021SBarry Smith PetscFunctionBegin; 3310700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3325f80ce2aSJacob Faibussowitsch PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 333c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, levels, 2); 334cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 3353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33685317021SBarry Smith } 33785317021SBarry Smith 33885317021SBarry Smith /*@ 33985317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 34085317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 34185317021SBarry Smith 342c3339decSBarry Smith Logically Collective 34385317021SBarry Smith 34485317021SBarry Smith Input Parameters: 34585317021SBarry Smith + pc - the preconditioner context 346f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 34785317021SBarry Smith 34885317021SBarry Smith Options Database Key: 34967b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 35085317021SBarry Smith 351f1580f4eSBarry Smith Note: 35285317021SBarry Smith Does not apply with 0 fill. 35385317021SBarry Smith 35485317021SBarry Smith Level: intermediate 35585317021SBarry Smith 356*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 35785317021SBarry Smith @*/ 358d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) 359d71ae5a4SJacob Faibussowitsch { 36085317021SBarry Smith PetscFunctionBegin; 3610700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 362cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 36492e9c092SBarry Smith } 36592e9c092SBarry Smith 36692e9c092SBarry Smith /*@ 36792e9c092SBarry Smith PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 36892e9c092SBarry Smith treated as level 0 fill even if there is no non-zero location. 36992e9c092SBarry Smith 370c3339decSBarry Smith Logically Collective 37192e9c092SBarry Smith 37292e9c092SBarry Smith Input Parameter: 37392e9c092SBarry Smith . pc - the preconditioner context 37492e9c092SBarry Smith 37592e9c092SBarry Smith Output Parameter: 376f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 37792e9c092SBarry Smith 378f1580f4eSBarry Smith Note: 37992e9c092SBarry Smith Does not apply with 0 fill. 38092e9c092SBarry Smith 38192e9c092SBarry Smith Level: intermediate 38292e9c092SBarry Smith 383*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 38492e9c092SBarry Smith @*/ 385d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) 386d71ae5a4SJacob Faibussowitsch { 38792e9c092SBarry Smith PetscFunctionBegin; 38892e9c092SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 389cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 3903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 39185317021SBarry Smith } 39285317021SBarry Smith 39385317021SBarry Smith /*@ 39485317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 39585317021SBarry Smith 396c3339decSBarry Smith Logically Collective 39785317021SBarry Smith 39885317021SBarry Smith Input Parameters: 39985317021SBarry Smith + pc - the preconditioner context 400feefa0e1SJacob Faibussowitsch - rtol - diagonal entries smaller than this in absolute value are considered zero 40185317021SBarry Smith 40285317021SBarry Smith Options Database Key: 403147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 40485317021SBarry Smith 40585317021SBarry Smith Level: intermediate 40685317021SBarry Smith 407*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 40885317021SBarry Smith @*/ 409d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) 410d71ae5a4SJacob Faibussowitsch { 41185317021SBarry Smith PetscFunctionBegin; 4120700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 413c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, rtol, 2); 414cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 4153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41685317021SBarry Smith } 41785317021SBarry Smith 418bf6011e8SBarry Smith /*@C 419f1580f4eSBarry Smith PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 42085317021SBarry Smith 421c3339decSBarry Smith Logically Collective 42285317021SBarry Smith 42385317021SBarry Smith Input Parameters: 42485317021SBarry Smith + pc - the preconditioner context 425f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 42685317021SBarry Smith 42785317021SBarry Smith Options Database Key: 4283ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 42985317021SBarry Smith 43085317021SBarry Smith Level: intermediate 43185317021SBarry Smith 43285317021SBarry Smith Note: 43385317021SBarry Smith By default this will use the PETSc factorization if it exists 43485317021SBarry Smith 435*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, 436f1580f4eSBarry Smith `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 43785317021SBarry Smith @*/ 438d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) 439d71ae5a4SJacob Faibussowitsch { 44085317021SBarry Smith PetscFunctionBegin; 4410700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 442cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 4433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 44485317021SBarry Smith } 44585317021SBarry Smith 446bf6011e8SBarry Smith /*@C 447f1580f4eSBarry Smith PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization 4487112b564SBarry Smith 449c5eb9154SBarry Smith Not Collective 4507112b564SBarry Smith 4517112b564SBarry Smith Input Parameter: 4527112b564SBarry Smith . pc - the preconditioner context 4537112b564SBarry Smith 4547112b564SBarry Smith Output Parameter: 455f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4567112b564SBarry Smith 4577112b564SBarry Smith Level: intermediate 4587112b564SBarry Smith 459*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `MATSOLVERSUPERLU`, 46042747ad1SJacob Faibussowitsch `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4617112b564SBarry Smith @*/ 462d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) 463d71ae5a4SJacob Faibussowitsch { 4645f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(PC, MatSolverType *); 4657112b564SBarry Smith 4667112b564SBarry Smith PetscFunctionBegin; 4670700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4684f572ea9SToby Isaac PetscAssertPointer(stype, 2); 4699566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 4709566063dSJacob Faibussowitsch if (f) PetscCall((*f)(pc, stype)); 4715f80ce2aSJacob Faibussowitsch else *stype = NULL; 4723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4737112b564SBarry Smith } 4747112b564SBarry Smith 47585317021SBarry Smith /*@ 47685317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 47785317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 47885317021SBarry Smith 479c5eb9154SBarry Smith Not Collective, each process can expect a different amount of fill 48085317021SBarry Smith 48185317021SBarry Smith Input Parameters: 48285317021SBarry Smith + pc - the preconditioner context 48385317021SBarry Smith - fill - amount of expected fill 48485317021SBarry Smith 48585317021SBarry Smith Options Database Key: 48685317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 48785317021SBarry Smith 48885317021SBarry Smith Level: intermediate 48985317021SBarry Smith 490f1580f4eSBarry Smith Notes: 49185317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 49285317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 49385317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 49485317021SBarry Smith future runs. Default PETSc uses a value of 5.0 49585317021SBarry Smith 496f1580f4eSBarry Smith This is ignored for most solver packages 49701a79839SBarry Smith 498f1580f4eSBarry Smith This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels. 499f1580f4eSBarry Smith 500*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()` 50185317021SBarry Smith @*/ 502d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) 503d71ae5a4SJacob Faibussowitsch { 50485317021SBarry Smith PetscFunctionBegin; 5050700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 506467446fbSPierre Jolivet PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less than 1.0"); 507cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50985317021SBarry Smith } 51085317021SBarry Smith 51185317021SBarry Smith /*@ 51204c3f3b8SBarry Smith PCFactorSetUseInPlace - Tells the preconditioner to do an in-place factorization. 51385317021SBarry Smith 514c3339decSBarry Smith Logically Collective 51585317021SBarry Smith 51685317021SBarry Smith Input Parameters: 5178e37d05fSBarry Smith + pc - the preconditioner context 518f1580f4eSBarry Smith - flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 51985317021SBarry Smith 52085317021SBarry Smith Options Database Key: 5218e37d05fSBarry Smith . -pc_factor_in_place <true,false> - Activate/deactivate in-place factorization 52285317021SBarry Smith 523f1580f4eSBarry Smith Note: 52404c3f3b8SBarry Smith For dense matrices, this enables the solution of much larger problems. 52504c3f3b8SBarry Smith For sparse matrices the factorization cannot be done truly in-place 52604c3f3b8SBarry Smith so this does not save memory during the factorization, but after the matrix 52704c3f3b8SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 52804c3f3b8SBarry Smith space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 52904c3f3b8SBarry Smith 530f1580f4eSBarry Smith `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when 53185317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 532f1580f4eSBarry Smith a call to `KSPSetOperators()`. 53385317021SBarry Smith This is because the Krylov space methods require an application of the 53485317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 53585317021SBarry Smith been factored in-place, replacing the original matrix. 53685317021SBarry Smith 53785317021SBarry Smith Level: intermediate 53885317021SBarry Smith 539*562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `Mat`, `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()` 54085317021SBarry Smith @*/ 541d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) 542d71ae5a4SJacob Faibussowitsch { 54385317021SBarry Smith PetscFunctionBegin; 5440700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 545cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 5463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5478e37d05fSBarry Smith } 5488e37d05fSBarry Smith 5498e37d05fSBarry Smith /*@ 5508e37d05fSBarry Smith PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 5518e37d05fSBarry Smith 552c3339decSBarry Smith Logically Collective 5538e37d05fSBarry Smith 5548e37d05fSBarry Smith Input Parameter: 5558e37d05fSBarry Smith . pc - the preconditioner context 5568e37d05fSBarry Smith 5578e37d05fSBarry Smith Output Parameter: 558f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 5598e37d05fSBarry Smith 5608e37d05fSBarry Smith Level: intermediate 5618e37d05fSBarry Smith 562*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()` 5638e37d05fSBarry Smith @*/ 564d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) 565d71ae5a4SJacob Faibussowitsch { 5668e37d05fSBarry Smith PetscFunctionBegin; 5678e37d05fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 568cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 5693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 57085317021SBarry Smith } 57185317021SBarry Smith 57285317021SBarry Smith /*@C 57385317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 574f1580f4eSBarry Smith be used in the `PCLU`, `PCCHOLESKY`, `PCILU`, or `PCICC` preconditioners 57585317021SBarry Smith 576c3339decSBarry Smith Logically Collective 57785317021SBarry Smith 57885317021SBarry Smith Input Parameters: 57985317021SBarry Smith + pc - the preconditioner context 580f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM` 58185317021SBarry Smith 58285317021SBarry Smith Options Database Key: 5832c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 58485317021SBarry Smith 58585317021SBarry Smith Level: intermediate 58685317021SBarry Smith 58795452b02SPatrick Sanan Notes: 5884ac6704cSBarry Smith Nested dissection is used by default for some of PETSc's sparse matrix formats 58985317021SBarry Smith 590f1580f4eSBarry 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 5919bd791bbSBarry Smith and reordering this matrix is very expensive. 5929bd791bbSBarry Smith 593f1580f4eSBarry Smith You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering. 5949bd791bbSBarry Smith 595f1580f4eSBarry Smith `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others. 5962c7c0729SBarry Smith 597*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM` 59885317021SBarry Smith @*/ 599d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) 600d71ae5a4SJacob Faibussowitsch { 60185317021SBarry Smith PetscFunctionBegin; 602c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 603cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 6043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 60585317021SBarry Smith } 60685317021SBarry Smith 60785317021SBarry Smith /*@ 6088ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 60985317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 610f1580f4eSBarry Smith it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used. 61185317021SBarry Smith 612c3339decSBarry Smith Logically Collective 61385317021SBarry Smith 61485317021SBarry Smith Input Parameters: 61585317021SBarry Smith + pc - the preconditioner context 61685317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 61785317021SBarry Smith 61885317021SBarry Smith Options Database Key: 619147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 62085317021SBarry Smith 62185317021SBarry Smith Level: intermediate 62285317021SBarry Smith 623*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 62485317021SBarry Smith @*/ 625d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) 626d71ae5a4SJacob Faibussowitsch { 62785317021SBarry Smith PetscFunctionBegin; 628c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 629c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, dtcol, 2); 630cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63285317021SBarry Smith } 63385317021SBarry Smith 63485317021SBarry Smith /*@ 63585317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 636f1580f4eSBarry Smith with `MATBAIJ` or `MATSBAIJ` matrices 63785317021SBarry Smith 638c3339decSBarry Smith Logically Collective 63985317021SBarry Smith 64085317021SBarry Smith Input Parameters: 64185317021SBarry Smith + pc - the preconditioner context 642f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE` 64385317021SBarry Smith 64485317021SBarry Smith Options Database Key: 645f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ` 64685317021SBarry Smith 64785317021SBarry Smith Level: intermediate 64885317021SBarry Smith 649*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 65085317021SBarry Smith @*/ 651d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) 652d71ae5a4SJacob Faibussowitsch { 65385317021SBarry Smith PetscFunctionBegin; 654c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 655acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, pivot, 2); 656cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 6573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65885317021SBarry Smith } 65985317021SBarry Smith 66085317021SBarry Smith /*@ 661288e7d53SBarry Smith PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 66285317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 66385317021SBarry Smith 664c3339decSBarry Smith Logically Collective 66585317021SBarry Smith 66685317021SBarry Smith Input Parameters: 66785317021SBarry Smith + pc - the preconditioner context 668f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 66985317021SBarry Smith 67085317021SBarry Smith Options Database Key: 671f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 67285317021SBarry Smith 67385317021SBarry Smith Level: intermediate 67485317021SBarry Smith 675*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()` 67685317021SBarry Smith @*/ 677d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) 678d71ae5a4SJacob Faibussowitsch { 67985317021SBarry Smith PetscFunctionBegin; 680064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 681acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 682cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 68485317021SBarry Smith } 6853d1c1ea0SBarry Smith 686d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) 687d71ae5a4SJacob Faibussowitsch { 6883d1c1ea0SBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 6893d1c1ea0SBarry Smith 6903d1c1ea0SBarry Smith PetscFunctionBegin; 6919566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&fact->info)); 6924ac6704cSBarry Smith fact->factortype = ftype; 6933d1c1ea0SBarry Smith fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 6943d1c1ea0SBarry Smith fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 6953d1c1ea0SBarry Smith fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 6963d1c1ea0SBarry Smith fact->info.pivotinblocks = 1.0; 6973d1c1ea0SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 6983d1c1ea0SBarry Smith 6999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 7049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 7069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 7129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 7139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 7149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 7159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 7169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 7179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 7193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7203d1c1ea0SBarry Smith } 7212e956fe4SStefano Zampini 722d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc) 723d71ae5a4SJacob Faibussowitsch { 7242e956fe4SStefano Zampini PetscFunctionBegin; 7252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 7262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 7272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 7282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 7292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 7302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 7312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 7322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 7332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 7342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 7352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 7362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 7372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 7382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 7392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 7402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 7412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 7422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 7432e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 7442e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 7452e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 7462e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 7473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7482e956fe4SStefano Zampini } 749