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) { 20*03e5aca4SStefano Zampini if (!fact->fact) { PetscCall(MatGetFactor(pc->pmat, fact->solvertype, fact->factortype, &fact->fact)); } 21*03e5aca4SStefano Zampini if (!fact->fact) PetscFunctionReturn(PETSC_SUCCESS); 22*03e5aca4SStefano Zampini 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. 90*03e5aca4SStefano Zampini This function raises an error if the requested combination of solver package and matrix type is not supported. 91f8260c8fSBarry Smith 922bd2b0e6SSatish Balay Level: intermediate 932bd2b0e6SSatish Balay 94f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 95f8260c8fSBarry Smith @*/ 96d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 97d71ae5a4SJacob Faibussowitsch { 98f8260c8fSBarry Smith PetscFunctionBegin; 99f8260c8fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 100cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUpMatSolverType_C", (PC), (pc)); 1013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 102f8260c8fSBarry Smith } 103f8260c8fSBarry Smith 104ee45ca4aSHong Zhang /*@ 105ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 106ee45ca4aSHong Zhang 107c3339decSBarry Smith Logically Collective 108ee45ca4aSHong Zhang 109ee45ca4aSHong Zhang Input Parameters: 110afaefe49SHong Zhang + pc - the preconditioner context 111afaefe49SHong Zhang - zero - all pivots smaller than this will be considered zero 112ee45ca4aSHong Zhang 113ee45ca4aSHong Zhang Options Database Key: 114ee45ca4aSHong Zhang . -pc_factor_zeropivot <zero> - Sets tolerance for what is considered a zero pivot 115ee45ca4aSHong Zhang 116ee45ca4aSHong Zhang Level: intermediate 117ee45ca4aSHong Zhang 118f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 119ee45ca4aSHong Zhang @*/ 120d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetZeroPivot(PC pc, PetscReal zero) 121d71ae5a4SJacob Faibussowitsch { 122ee45ca4aSHong Zhang PetscFunctionBegin; 1230700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 124c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, zero, 2); 125cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetZeroPivot_C", (PC, PetscReal), (pc, zero)); 1263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 127ee45ca4aSHong Zhang } 128ee45ca4aSHong Zhang 129915743fcSHong Zhang /*@ 130915743fcSHong Zhang PCFactorSetShiftType - adds a particular type of quantity to the diagonal of the matrix during 131915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 132915743fcSHong Zhang 133c3339decSBarry Smith Logically Collective 134915743fcSHong Zhang 135915743fcSHong Zhang Input Parameters: 136915743fcSHong Zhang + pc - the preconditioner context 137f1580f4eSBarry Smith - shifttype - type of shift; one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, `MAT_SHIFT_INBLOCKS` 138915743fcSHong Zhang 139915743fcSHong Zhang Options Database Key: 14028d58a37SPierre Jolivet . -pc_factor_shift_type <shifttype> - Sets shift type; use '-help' for a list of available types 141915743fcSHong Zhang 142915743fcSHong Zhang Level: intermediate 143915743fcSHong Zhang 144f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 145915743fcSHong Zhang @*/ 146d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftType(PC pc, MatFactorShiftType shifttype) 147d71ae5a4SJacob Faibussowitsch { 148d90ac83dSHong Zhang PetscFunctionBegin; 1490700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 150c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, shifttype, 2); 151cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftType_C", (PC, MatFactorShiftType), (pc, shifttype)); 1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153d90ac83dSHong Zhang } 154d90ac83dSHong Zhang 155915743fcSHong Zhang /*@ 156915743fcSHong Zhang PCFactorSetShiftAmount - adds a quantity to the diagonal of the matrix during 157915743fcSHong Zhang numerical factorization, thus the matrix has nonzero pivots 158915743fcSHong Zhang 159c3339decSBarry Smith Logically Collective 160915743fcSHong Zhang 161915743fcSHong Zhang Input Parameters: 162915743fcSHong Zhang + pc - the preconditioner context 163f1580f4eSBarry Smith - shiftamount - amount of shift or `PETSC_DECIDE` for the default 164915743fcSHong Zhang 165915743fcSHong Zhang Options Database Key: 166f1580f4eSBarry Smith . -pc_factor_shift_amount <shiftamount> - Sets shift amount or -1 for the default 167915743fcSHong Zhang 168915743fcSHong Zhang Level: intermediate 169915743fcSHong Zhang 170f1580f4eSBarry Smith .seealso: `PCCHOLESKY`, `PCLU`, ``PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 171915743fcSHong Zhang @*/ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetShiftAmount(PC pc, PetscReal shiftamount) 173d71ae5a4SJacob Faibussowitsch { 174d90ac83dSHong Zhang PetscFunctionBegin; 1750700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 176c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, shiftamount, 2); 177cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetShiftAmount_C", (PC, PetscReal), (pc, shiftamount)); 1783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 179d90ac83dSHong Zhang } 180d90ac83dSHong Zhang 1816d33885cSprj- /*@ 182f1580f4eSBarry Smith PCFactorSetDropTolerance - The preconditioner will use an `PCILU` 183f1580f4eSBarry Smith based on a drop tolerance. 18485317021SBarry Smith 185c3339decSBarry Smith Logically Collective 18685317021SBarry Smith 18785317021SBarry Smith Input Parameters: 18885317021SBarry Smith + pc - the preconditioner context 18985317021SBarry Smith . dt - the drop tolerance, try from 1.e-10 to .1 19085317021SBarry Smith . dtcol - tolerance for column pivot, good values [0.1 to 0.01] 19185317021SBarry Smith - maxrowcount - the max number of nonzeros allowed in a row, best value 19285317021SBarry Smith depends on the number of nonzeros in row of original matrix 19385317021SBarry Smith 19485317021SBarry Smith Options Database Key: 195b7c853c4SBarry Smith . -pc_factor_drop_tolerance <dt,dtcol,maxrowcount> - Sets drop tolerance 19685317021SBarry Smith 19785317021SBarry Smith Level: intermediate 19885317021SBarry Smith 199f1580f4eSBarry Smith Note: 20085317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 20185317021SBarry Smith matrix. We don't know how to compute reasonable values. 20285317021SBarry Smith 203f1580f4eSBarry Smith .seealso: `PCILU` 2046d33885cSprj- @*/ 205d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetDropTolerance(PC pc, PetscReal dt, PetscReal dtcol, PetscInt maxrowcount) 206d71ae5a4SJacob Faibussowitsch { 20785317021SBarry Smith PetscFunctionBegin; 2080700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 209064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveReal(pc, dtcol, 3); 210064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(pc, maxrowcount, 4); 211cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetDropTolerance_C", (PC, PetscReal, PetscReal, PetscInt), (pc, dt, dtcol, maxrowcount)); 2123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21385317021SBarry Smith } 21485317021SBarry Smith 215c7f610a1SBarry Smith /*@ 216c7f610a1SBarry Smith PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 217c7f610a1SBarry Smith 218c7f610a1SBarry Smith Not Collective 219c7f610a1SBarry Smith 2202fe279fdSBarry Smith Input Parameter: 221c7f610a1SBarry Smith . pc - the preconditioner context 222c7f610a1SBarry Smith 223c7f610a1SBarry Smith Output Parameter: 224c7f610a1SBarry Smith . pivot - the tolerance 225c7f610a1SBarry Smith 226c7f610a1SBarry Smith Level: intermediate 227c7f610a1SBarry Smith 228f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetZeroPivot()` 229c7f610a1SBarry Smith @*/ 230d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetZeroPivot(PC pc, PetscReal *pivot) 231d71ae5a4SJacob Faibussowitsch { 232c7f610a1SBarry Smith PetscFunctionBegin; 233c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 234cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetZeroPivot_C", (PC, PetscReal *), (pc, pivot)); 2353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 236c7f610a1SBarry Smith } 237c7f610a1SBarry Smith 238c7f610a1SBarry Smith /*@ 239c7f610a1SBarry Smith PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 240c7f610a1SBarry Smith 241c7f610a1SBarry Smith Not Collective 242c7f610a1SBarry Smith 2432fe279fdSBarry Smith Input Parameter: 244c7f610a1SBarry Smith . pc - the preconditioner context 245c7f610a1SBarry Smith 246c7f610a1SBarry Smith Output Parameter: 247c7f610a1SBarry Smith . shift - how much to shift the diagonal entry 248c7f610a1SBarry Smith 249c7f610a1SBarry Smith Level: intermediate 250c7f610a1SBarry Smith 251f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 252c7f610a1SBarry Smith @*/ 253d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftAmount(PC pc, PetscReal *shift) 254d71ae5a4SJacob Faibussowitsch { 255c7f610a1SBarry Smith PetscFunctionBegin; 256c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 257cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftAmount_C", (PC, PetscReal *), (pc, shift)); 2583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 259c7f610a1SBarry Smith } 260c7f610a1SBarry Smith 261c7f610a1SBarry Smith /*@ 262c7f610a1SBarry Smith PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 263c7f610a1SBarry Smith 264c7f610a1SBarry Smith Not Collective 265c7f610a1SBarry Smith 266f1580f4eSBarry Smith Input Parameter: 267c7f610a1SBarry Smith . pc - the preconditioner context 268c7f610a1SBarry Smith 269c7f610a1SBarry Smith Output Parameter: 270f1580f4eSBarry Smith . type - one of `MAT_SHIFT_NONE`, `MAT_SHIFT_NONZERO`, `MAT_SHIFT_POSITIVE_DEFINITE`, or `MAT_SHIFT_INBLOCKS` 271c7f610a1SBarry Smith 272c7f610a1SBarry Smith Level: intermediate 273c7f610a1SBarry Smith 274f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 275c7f610a1SBarry Smith @*/ 276d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetShiftType(PC pc, MatFactorShiftType *type) 277d71ae5a4SJacob Faibussowitsch { 278c7f610a1SBarry Smith PetscFunctionBegin; 279c7f610a1SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 280cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetShiftType_C", (PC, MatFactorShiftType *), (pc, type)); 2813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 282c7f610a1SBarry Smith } 283c7f610a1SBarry Smith 2842591b318SToby Isaac /*@ 2852591b318SToby Isaac PCFactorGetLevels - Gets the number of levels of fill to use. 2862591b318SToby Isaac 287c3339decSBarry Smith Logically Collective 2882591b318SToby Isaac 2892fe279fdSBarry Smith Input Parameter: 2902591b318SToby Isaac . pc - the preconditioner context 2912591b318SToby Isaac 2922591b318SToby Isaac Output Parameter: 2932591b318SToby Isaac . levels - number of levels of fill 2942591b318SToby Isaac 2952591b318SToby Isaac Level: intermediate 2962591b318SToby Isaac 297f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetLevels()` 2982591b318SToby Isaac @*/ 299d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetLevels(PC pc, PetscInt *levels) 300d71ae5a4SJacob Faibussowitsch { 3012591b318SToby Isaac PetscFunctionBegin; 3022591b318SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 303cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetLevels_C", (PC, PetscInt *), (pc, levels)); 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3052591b318SToby Isaac } 3062591b318SToby Isaac 30785317021SBarry Smith /*@ 30885317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 30985317021SBarry Smith 310c3339decSBarry Smith Logically Collective 31185317021SBarry Smith 31285317021SBarry Smith Input Parameters: 31385317021SBarry Smith + pc - the preconditioner context 31485317021SBarry Smith - levels - number of levels of fill 31585317021SBarry Smith 31685317021SBarry Smith Options Database Key: 31785317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 31885317021SBarry Smith 31985317021SBarry Smith Level: intermediate 32085317021SBarry Smith 321f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetLevels()` 32285317021SBarry Smith @*/ 323d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetLevels(PC pc, PetscInt levels) 324d71ae5a4SJacob Faibussowitsch { 32585317021SBarry Smith PetscFunctionBegin; 3260700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3275f80ce2aSJacob Faibussowitsch PetscCheck(levels >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "negative levels"); 328c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, levels, 2); 329cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetLevels_C", (PC, PetscInt), (pc, levels)); 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33185317021SBarry Smith } 33285317021SBarry Smith 33385317021SBarry Smith /*@ 33485317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 33585317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 33685317021SBarry Smith 337c3339decSBarry Smith Logically Collective 33885317021SBarry Smith 33985317021SBarry Smith Input Parameters: 34085317021SBarry Smith + pc - the preconditioner context 341f1580f4eSBarry Smith - flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 34285317021SBarry Smith 34385317021SBarry Smith Options Database Key: 34467b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 34585317021SBarry Smith 346f1580f4eSBarry Smith Note: 34785317021SBarry Smith Does not apply with 0 fill. 34885317021SBarry Smith 34985317021SBarry Smith Level: intermediate 35085317021SBarry Smith 351f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorGetAllowDiagonalFill()` 35285317021SBarry Smith @*/ 353d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc, PetscBool flg) 354d71ae5a4SJacob Faibussowitsch { 35585317021SBarry Smith PetscFunctionBegin; 3560700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 357cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetAllowDiagonalFill_C", (PC, PetscBool), (pc, flg)); 3583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35992e9c092SBarry Smith } 36092e9c092SBarry Smith 36192e9c092SBarry Smith /*@ 36292e9c092SBarry Smith PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 36392e9c092SBarry Smith treated as level 0 fill even if there is no non-zero location. 36492e9c092SBarry Smith 365c3339decSBarry Smith Logically Collective 36692e9c092SBarry Smith 36792e9c092SBarry Smith Input Parameter: 36892e9c092SBarry Smith . pc - the preconditioner context 36992e9c092SBarry Smith 37092e9c092SBarry Smith Output Parameter: 371f1580f4eSBarry Smith . flg - `PETSC_TRUE` to turn on, `PETSC_FALSE` to turn off 37292e9c092SBarry Smith 373f1580f4eSBarry Smith Note: 37492e9c092SBarry Smith Does not apply with 0 fill. 37592e9c092SBarry Smith 37692e9c092SBarry Smith Level: intermediate 37792e9c092SBarry Smith 378f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetAllowDiagonalFill()` 37992e9c092SBarry Smith @*/ 380d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc, PetscBool *flg) 381d71ae5a4SJacob Faibussowitsch { 38292e9c092SBarry Smith PetscFunctionBegin; 38392e9c092SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 384cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetAllowDiagonalFill_C", (PC, PetscBool *), (pc, flg)); 3853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38685317021SBarry Smith } 38785317021SBarry Smith 38885317021SBarry Smith /*@ 38985317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 39085317021SBarry Smith 391c3339decSBarry Smith Logically Collective 39285317021SBarry Smith 39385317021SBarry Smith Input Parameters: 39485317021SBarry Smith + pc - the preconditioner context 39585317021SBarry Smith - tol - diagonal entries smaller than this in absolute value are considered zero 39685317021SBarry Smith 39785317021SBarry Smith Options Database Key: 398147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 39985317021SBarry Smith 40085317021SBarry Smith Level: intermediate 40185317021SBarry Smith 402f1580f4eSBarry Smith .seealso: `PCILU`, `PCICC`, `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 40385317021SBarry Smith @*/ 404d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc, PetscReal rtol) 405d71ae5a4SJacob Faibussowitsch { 40685317021SBarry Smith PetscFunctionBegin; 4070700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 408c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, rtol, 2); 409cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorReorderForNonzeroDiagonal_C", (PC, PetscReal), (pc, rtol)); 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41185317021SBarry Smith } 41285317021SBarry Smith 413bf6011e8SBarry Smith /*@C 414f1580f4eSBarry Smith PCFactorSetMatSolverType - sets the solver package that is used to perform the factorization 41585317021SBarry Smith 416c3339decSBarry Smith Logically Collective 41785317021SBarry Smith 41885317021SBarry Smith Input Parameters: 41985317021SBarry Smith + pc - the preconditioner context 420f1580f4eSBarry Smith - stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 42185317021SBarry Smith 42285317021SBarry Smith Options Database Key: 4233ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 42485317021SBarry Smith 42585317021SBarry Smith Level: intermediate 42685317021SBarry Smith 42785317021SBarry Smith Note: 42885317021SBarry Smith By default this will use the PETSc factorization if it exists 42985317021SBarry Smith 430f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 431f1580f4eSBarry Smith `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 43285317021SBarry Smith @*/ 433d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatSolverType(PC pc, MatSolverType stype) 434d71ae5a4SJacob Faibussowitsch { 43585317021SBarry Smith PetscFunctionBegin; 4360700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 437cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatSolverType_C", (PC, MatSolverType), (pc, stype)); 4383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43985317021SBarry Smith } 44085317021SBarry Smith 441bf6011e8SBarry Smith /*@C 442f1580f4eSBarry Smith PCFactorGetMatSolverType - gets the solver package that is used to perform the factorization 4437112b564SBarry Smith 444c5eb9154SBarry Smith Not Collective 4457112b564SBarry Smith 4467112b564SBarry Smith Input Parameter: 4477112b564SBarry Smith . pc - the preconditioner context 4487112b564SBarry Smith 4497112b564SBarry Smith Output Parameter: 450f1580f4eSBarry Smith . stype - for example, `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4517112b564SBarry Smith 4527112b564SBarry Smith Level: intermediate 4537112b564SBarry Smith 454f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()`, `MatSolverType`, 455f1580f4eSBarry Smith `MATSOLVERSUPERLU`, `MATSOLVERSUPERLU_DIST`, `MATSOLVERMUMPS` 4567112b564SBarry Smith @*/ 457d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatSolverType(PC pc, MatSolverType *stype) 458d71ae5a4SJacob Faibussowitsch { 4595f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(PC, MatSolverType *); 4607112b564SBarry Smith 4617112b564SBarry Smith PetscFunctionBegin; 4620700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4635f80ce2aSJacob Faibussowitsch PetscValidPointer(stype, 2); 4649566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", &f)); 4659566063dSJacob Faibussowitsch if (f) PetscCall((*f)(pc, stype)); 4665f80ce2aSJacob Faibussowitsch else *stype = NULL; 4673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4687112b564SBarry Smith } 4697112b564SBarry Smith 47085317021SBarry Smith /*@ 47185317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 47285317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 47385317021SBarry Smith 474c5eb9154SBarry Smith Not Collective, each process can expect a different amount of fill 47585317021SBarry Smith 47685317021SBarry Smith Input Parameters: 47785317021SBarry Smith + pc - the preconditioner context 47885317021SBarry Smith - fill - amount of expected fill 47985317021SBarry Smith 48085317021SBarry Smith Options Database Key: 48185317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 48285317021SBarry Smith 48385317021SBarry Smith Level: intermediate 48485317021SBarry Smith 485f1580f4eSBarry Smith Notes: 48685317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 48785317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 48885317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 48985317021SBarry Smith future runs. Default PETSc uses a value of 5.0 49085317021SBarry Smith 491f1580f4eSBarry Smith This is ignored for most solver packages 49201a79839SBarry Smith 493f1580f4eSBarry Smith This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with `PCFactorSetLevels()` or -pc_factor_levels. 494f1580f4eSBarry Smith 495f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseFill()` 49685317021SBarry Smith @*/ 497d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetFill(PC pc, PetscReal fill) 498d71ae5a4SJacob Faibussowitsch { 49985317021SBarry Smith PetscFunctionBegin; 5000700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5015f80ce2aSJacob Faibussowitsch PetscCheck(fill >= 1.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Fill factor cannot be less then 1.0"); 502cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetFill_C", (PC, PetscReal), (pc, fill)); 5033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50485317021SBarry Smith } 50585317021SBarry Smith 50685317021SBarry Smith /*@ 50785317021SBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 50885317021SBarry Smith For dense matrices, this enables the solution of much larger problems. 50985317021SBarry Smith For sparse matrices the factorization cannot be done truly in-place 51085317021SBarry Smith so this does not save memory during the factorization, but after the matrix 51185317021SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 512ec5066bdSBarry Smith space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 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: 524f1580f4eSBarry Smith `PCFactorSetUseInplace()` can only be used with the `KSP` method `KSPPREONLY` or when 52585317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 526f1580f4eSBarry Smith a call to `KSPSetOperators()`. 52785317021SBarry Smith This is because the Krylov space methods require an application of the 52885317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 52985317021SBarry Smith been factored in-place, replacing the original matrix. 53085317021SBarry Smith 53185317021SBarry Smith Level: intermediate 53285317021SBarry Smith 533f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorGetUseInPlace()` 53485317021SBarry Smith @*/ 535d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetUseInPlace(PC pc, PetscBool flg) 536d71ae5a4SJacob Faibussowitsch { 53785317021SBarry Smith PetscFunctionBegin; 5380700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 539cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetUseInPlace_C", (PC, PetscBool), (pc, flg)); 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5418e37d05fSBarry Smith } 5428e37d05fSBarry Smith 5438e37d05fSBarry Smith /*@ 5448e37d05fSBarry Smith PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 5458e37d05fSBarry Smith 546c3339decSBarry Smith Logically Collective 5478e37d05fSBarry Smith 5488e37d05fSBarry Smith Input Parameter: 5498e37d05fSBarry Smith . pc - the preconditioner context 5508e37d05fSBarry Smith 5518e37d05fSBarry Smith Output Parameter: 552f1580f4eSBarry Smith . flg - `PETSC_TRUE` to enable, `PETSC_FALSE` to disable 5538e37d05fSBarry Smith 5548e37d05fSBarry Smith Level: intermediate 5558e37d05fSBarry Smith 556f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetUseInPlace()` 5578e37d05fSBarry Smith @*/ 558d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetUseInPlace(PC pc, PetscBool *flg) 559d71ae5a4SJacob Faibussowitsch { 5608e37d05fSBarry Smith PetscFunctionBegin; 5618e37d05fSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 562cac4c232SBarry Smith PetscUseMethod(pc, "PCFactorGetUseInPlace_C", (PC, PetscBool *), (pc, flg)); 5633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56485317021SBarry Smith } 56585317021SBarry Smith 56685317021SBarry Smith /*@C 56785317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 568f1580f4eSBarry Smith be used in the `PCLU`, `PCCHOLESKY`, `PCILU`, or `PCICC` preconditioners 56985317021SBarry Smith 570c3339decSBarry Smith Logically Collective 57185317021SBarry Smith 57285317021SBarry Smith Input Parameters: 57385317021SBarry Smith + pc - the preconditioner context 574f1580f4eSBarry Smith - ordering - the matrix ordering name, for example, `MATORDERINGND` or `MATORDERINGRCM` 57585317021SBarry Smith 57685317021SBarry Smith Options Database Key: 5772c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 57885317021SBarry Smith 57985317021SBarry Smith Level: intermediate 58085317021SBarry Smith 58195452b02SPatrick Sanan Notes: 5824ac6704cSBarry Smith Nested dissection is used by default for some of PETSc's sparse matrix formats 58385317021SBarry Smith 584f1580f4eSBarry 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 5859bd791bbSBarry Smith and reordering this matrix is very expensive. 5869bd791bbSBarry Smith 587f1580f4eSBarry Smith You can use a `MATSEQAIJ` matrix with Cholesky and ICC and use any ordering. 5889bd791bbSBarry Smith 589f1580f4eSBarry Smith `MATORDERINGEXTERNAL` means PETSc will not compute an ordering and the package will use its own ordering, usable with `MATSOLVERCHOLMOD`, `MATSOLVERUMFPACK`, and others. 5902c7c0729SBarry Smith 591f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `MatOrderingType`, `MATORDERINGEXTERNAL`, `MATORDERINGND`, `MATORDERINGRCM` 59285317021SBarry Smith @*/ 593d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetMatOrderingType(PC pc, MatOrderingType ordering) 594d71ae5a4SJacob Faibussowitsch { 59585317021SBarry Smith PetscFunctionBegin; 596c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 597cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetMatOrderingType_C", (PC, MatOrderingType), (pc, ordering)); 5983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 59985317021SBarry Smith } 60085317021SBarry Smith 60185317021SBarry Smith /*@ 6028ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 60385317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 604f1580f4eSBarry Smith it is never done. For the MATLAB and `MATSOLVERSUPERLU` factorization this is used. 60585317021SBarry Smith 606c3339decSBarry Smith Logically Collective 60785317021SBarry Smith 60885317021SBarry Smith Input Parameters: 60985317021SBarry Smith + pc - the preconditioner context 61085317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 61185317021SBarry Smith 61285317021SBarry Smith Options Database Key: 613147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 61485317021SBarry Smith 61585317021SBarry Smith Level: intermediate 61685317021SBarry Smith 617f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 61885317021SBarry Smith @*/ 619d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetColumnPivot(PC pc, PetscReal dtcol) 620d71ae5a4SJacob Faibussowitsch { 62185317021SBarry Smith PetscFunctionBegin; 622c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 623c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, dtcol, 2); 624cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetColumnPivot_C", (PC, PetscReal), (pc, dtcol)); 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62685317021SBarry Smith } 62785317021SBarry Smith 62885317021SBarry Smith /*@ 62985317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 630f1580f4eSBarry Smith with `MATBAIJ` or `MATSBAIJ` matrices 63185317021SBarry Smith 632c3339decSBarry Smith Logically Collective 63385317021SBarry Smith 63485317021SBarry Smith Input Parameters: 63585317021SBarry Smith + pc - the preconditioner context 636f1580f4eSBarry Smith - pivot - `PETSC_TRUE` or `PETSC_FALSE` 63785317021SBarry Smith 63885317021SBarry Smith Options Database Key: 639f1580f4eSBarry Smith . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for `MATBAIJ` and `MATSBAIJ` 64085317021SBarry Smith 64185317021SBarry Smith Level: intermediate 64285317021SBarry Smith 643f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 64485317021SBarry Smith @*/ 645d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetPivotInBlocks(PC pc, PetscBool pivot) 646d71ae5a4SJacob Faibussowitsch { 64785317021SBarry Smith PetscFunctionBegin; 648c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 649acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, pivot, 2); 650cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetPivotInBlocks_C", (PC, PetscBool), (pc, pivot)); 6513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65285317021SBarry Smith } 65385317021SBarry Smith 65485317021SBarry Smith /*@ 655288e7d53SBarry Smith PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 65685317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 65785317021SBarry Smith 658c3339decSBarry Smith Logically Collective 65985317021SBarry Smith 66085317021SBarry Smith Input Parameters: 66185317021SBarry Smith + pc - the preconditioner context 662f1580f4eSBarry Smith - flag - `PETSC_TRUE` to reuse else `PETSC_FALSE` 66385317021SBarry Smith 66485317021SBarry Smith Options Database Key: 665f1580f4eSBarry Smith . -pc_factor_reuse_fill - Activates `PCFactorSetReuseFill()` 66685317021SBarry Smith 66785317021SBarry Smith Level: intermediate 66885317021SBarry Smith 669f1580f4eSBarry Smith .seealso: `PCLU`, `PCCHOLESKY`, `PCILU`, `PCICC`, `PCFactorSetReuseOrdering()`, `PCFactorSetFill()` 67085317021SBarry Smith @*/ 671d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorSetReuseFill(PC pc, PetscBool flag) 672d71ae5a4SJacob Faibussowitsch { 67385317021SBarry Smith PetscFunctionBegin; 674064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 675acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, flag, 2); 676cac4c232SBarry Smith PetscTryMethod(pc, "PCFactorSetReuseFill_C", (PC, PetscBool), (pc, flag)); 6773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 67885317021SBarry Smith } 6793d1c1ea0SBarry Smith 680d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorInitialize(PC pc, MatFactorType ftype) 681d71ae5a4SJacob Faibussowitsch { 6823d1c1ea0SBarry Smith PC_Factor *fact = (PC_Factor *)pc->data; 6833d1c1ea0SBarry Smith 6843d1c1ea0SBarry Smith PetscFunctionBegin; 6859566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&fact->info)); 6864ac6704cSBarry Smith fact->factortype = ftype; 6873d1c1ea0SBarry Smith fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 6883d1c1ea0SBarry Smith fact->info.shiftamount = 100.0 * PETSC_MACHINE_EPSILON; 6893d1c1ea0SBarry Smith fact->info.zeropivot = 100.0 * PETSC_MACHINE_EPSILON; 6903d1c1ea0SBarry Smith fact->info.pivotinblocks = 1.0; 6913d1c1ea0SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 6923d1c1ea0SBarry Smith 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", PCFactorSetZeroPivot_Factor)); 6949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", PCFactorGetZeroPivot_Factor)); 6959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", PCFactorSetShiftType_Factor)); 6969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", PCFactorGetShiftType_Factor)); 6979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", PCFactorSetShiftAmount_Factor)); 6989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", PCFactorGetShiftAmount_Factor)); 6999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", PCFactorGetMatSolverType_Factor)); 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", PCFactorSetMatSolverType_Factor)); 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", PCFactorSetUpMatSolverType_Factor)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", PCFactorSetFill_Factor)); 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", PCFactorSetMatOrderingType_Factor)); 7049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", PCFactorSetLevels_Factor)); 7059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", PCFactorGetLevels_Factor)); 7069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", PCFactorSetAllowDiagonalFill_Factor)); 7079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", PCFactorGetAllowDiagonalFill_Factor)); 7089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", PCFactorSetPivotInBlocks_Factor)); 7099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", PCFactorSetUseInPlace_Factor)); 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", PCFactorGetUseInPlace_Factor)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", PCFactorSetReuseOrdering_Factor)); 7129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", PCFactorSetReuseFill_Factor)); 7133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7143d1c1ea0SBarry Smith } 7152e956fe4SStefano Zampini 716d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorClearComposedFunctions(PC pc) 717d71ae5a4SJacob Faibussowitsch { 7182e956fe4SStefano Zampini PetscFunctionBegin; 7192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetZeroPivot_C", NULL)); 7202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetZeroPivot_C", NULL)); 7212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftType_C", NULL)); 7222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftType_C", NULL)); 7232e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetShiftAmount_C", NULL)); 7242e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetShiftAmount_C", NULL)); 7252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetMatSolverType_C", NULL)); 7262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatSolverType_C", NULL)); 7272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUpMatSolverType_C", NULL)); 7282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetFill_C", NULL)); 7292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetMatOrderingType_C", NULL)); 7302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetLevels_C", NULL)); 7312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetLevels_C", NULL)); 7322e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetAllowDiagonalFill_C", NULL)); 7332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetAllowDiagonalFill_C", NULL)); 7342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetPivotInBlocks_C", NULL)); 7352e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetUseInPlace_C", NULL)); 7362e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorGetUseInPlace_C", NULL)); 7372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseOrdering_C", NULL)); 7382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetReuseFill_C", NULL)); 7392e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorReorderForNonzeroDiagonal_C", NULL)); 7402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCFactorSetDropTolerance_C", NULL)); 7413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7422e956fe4SStefano Zampini } 743