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 */ 84ac6704cSBarry Smith PetscErrorCode PCFactorSetDefaultOrdering_Factor(PC pc) 94ac6704cSBarry Smith { 10978beb7fSPierre Jolivet PetscBool foundmtype,flg; 11*26cc229bSBarry Smith const char *prefix; 124ac6704cSBarry Smith 134ac6704cSBarry Smith PetscFunctionBegin; 144ac6704cSBarry Smith if (pc->pmat) { 15*26cc229bSBarry Smith PetscCall(PCGetOptionsPrefix(pc,&prefix)); 16*26cc229bSBarry 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) { 25*26cc229bSBarry 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 } 424ac6704cSBarry Smith PetscFunctionReturn(0); 434ac6704cSBarry Smith } 444ac6704cSBarry Smith 453d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseOrdering_Factor(PC pc,PetscBool flag) 463d1c1ea0SBarry Smith { 473d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 483d1c1ea0SBarry Smith 493d1c1ea0SBarry Smith PetscFunctionBegin; 503d1c1ea0SBarry Smith lu->reuseordering = flag; 513d1c1ea0SBarry Smith PetscFunctionReturn(0); 523d1c1ea0SBarry Smith } 533d1c1ea0SBarry Smith 543d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetReuseFill_Factor(PC pc,PetscBool flag) 553d1c1ea0SBarry Smith { 563d1c1ea0SBarry Smith PC_Factor *lu = (PC_Factor*)pc->data; 573d1c1ea0SBarry Smith 583d1c1ea0SBarry Smith PetscFunctionBegin; 593d1c1ea0SBarry Smith lu->reusefill = flag; 603d1c1ea0SBarry Smith PetscFunctionReturn(0); 613d1c1ea0SBarry Smith } 623d1c1ea0SBarry Smith 633d1c1ea0SBarry Smith static PetscErrorCode PCFactorSetUseInPlace_Factor(PC pc,PetscBool flg) 643d1c1ea0SBarry Smith { 653d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 663d1c1ea0SBarry Smith 673d1c1ea0SBarry Smith PetscFunctionBegin; 683d1c1ea0SBarry Smith dir->inplace = flg; 693d1c1ea0SBarry Smith PetscFunctionReturn(0); 703d1c1ea0SBarry Smith } 713d1c1ea0SBarry Smith 723d1c1ea0SBarry Smith static PetscErrorCode PCFactorGetUseInPlace_Factor(PC pc,PetscBool *flg) 733d1c1ea0SBarry Smith { 743d1c1ea0SBarry Smith PC_Factor *dir = (PC_Factor*)pc->data; 753d1c1ea0SBarry Smith 763d1c1ea0SBarry Smith PetscFunctionBegin; 773d1c1ea0SBarry Smith *flg = dir->inplace; 783d1c1ea0SBarry Smith PetscFunctionReturn(0); 793d1c1ea0SBarry Smith } 803d1c1ea0SBarry Smith 81f8260c8fSBarry Smith /*@ 823ca39a21SBarry 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 8895452b02SPatrick Sanan Notes: 8995452b02SPatrick Sanan 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 93db781477SPatrick Sanan .seealso: `PCFactorSetMatSolverType()`, `PCFactorGetMatrix()` 94f8260c8fSBarry Smith @*/ 953ca39a21SBarry Smith PetscErrorCode PCFactorSetUpMatSolverType(PC pc) 96f8260c8fSBarry Smith { 97f8260c8fSBarry Smith PetscFunctionBegin; 98f8260c8fSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 99cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetUpMatSolverType_C",(PC),(pc)); 100b3a44c85SBarry Smith PetscFunctionReturn(0); 101f8260c8fSBarry Smith } 102f8260c8fSBarry Smith 103ee45ca4aSHong Zhang /*@ 104ee45ca4aSHong Zhang PCFactorSetZeroPivot - Sets the size at which smaller pivots are declared to be zero 105ee45ca4aSHong Zhang 106ad4df100SBarry Smith Logically Collective on PC 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 117db781477SPatrick Sanan .seealso: `PCFactorSetShiftType()`, `PCFactorSetShiftAmount()` 118ee45ca4aSHong Zhang @*/ 1197087cfbeSBarry Smith PetscErrorCode PCFactorSetZeroPivot(PC pc,PetscReal zero) 120ee45ca4aSHong Zhang { 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)); 125ee45ca4aSHong Zhang PetscFunctionReturn(0); 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 132ad4df100SBarry Smith Logically Collective on PC 133915743fcSHong Zhang 134915743fcSHong Zhang Input Parameters: 135915743fcSHong Zhang + pc - the preconditioner context 136915743fcSHong Zhang - 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 143db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftAmount()` 144915743fcSHong Zhang @*/ 1457087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftType(PC pc,MatFactorShiftType shifttype) 146d90ac83dSHong Zhang { 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)); 151d90ac83dSHong Zhang PetscFunctionReturn(0); 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 158ad4df100SBarry Smith Logically Collective on PC 159915743fcSHong Zhang 160915743fcSHong Zhang Input Parameters: 161915743fcSHong Zhang + pc - the preconditioner context 162915743fcSHong Zhang - shiftamount - amount of shift 163915743fcSHong Zhang 164915743fcSHong Zhang Options Database Key: 165915743fcSHong Zhang . -pc_factor_shift_amount <shiftamount> - Sets shift amount or PETSC_DECIDE for the default 166915743fcSHong Zhang 167915743fcSHong Zhang Level: intermediate 168915743fcSHong Zhang 169db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()`, `PCFactorSetShiftType()` 170915743fcSHong Zhang @*/ 1717087cfbeSBarry Smith PetscErrorCode PCFactorSetShiftAmount(PC pc,PetscReal shiftamount) 172d90ac83dSHong Zhang { 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)); 177d90ac83dSHong Zhang PetscFunctionReturn(0); 178d90ac83dSHong Zhang } 179d90ac83dSHong Zhang 1806d33885cSprj- /*@ 181b7c853c4SBarry Smith PCFactorSetDropTolerance - The preconditioner will use an ILU 18278fc6b22SHong Zhang based on a drop tolerance. (Under development) 18385317021SBarry Smith 184ad4df100SBarry Smith Logically Collective on PC 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 19885317021SBarry Smith There are NO default values for the 3 parameters, you must set them with reasonable values for your 19985317021SBarry Smith matrix. We don't know how to compute reasonable values. 20085317021SBarry Smith 2016d33885cSprj- @*/ 2027087cfbeSBarry Smith PetscErrorCode PCFactorSetDropTolerance(PC pc,PetscReal dt,PetscReal dtcol,PetscInt maxrowcount) 20385317021SBarry Smith { 20485317021SBarry Smith PetscFunctionBegin; 2050700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 206064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveReal(pc,dtcol,3); 207064a246eSJacob Faibussowitsch PetscValidLogicalCollectiveInt(pc,maxrowcount,4); 208cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetDropTolerance_C",(PC,PetscReal,PetscReal,PetscInt),(pc,dt,dtcol,maxrowcount)); 20985317021SBarry Smith PetscFunctionReturn(0); 21085317021SBarry Smith } 21185317021SBarry Smith 212c7f610a1SBarry Smith /*@ 213c7f610a1SBarry Smith PCFactorGetZeroPivot - Gets the tolerance used to define a zero privot 214c7f610a1SBarry Smith 215c7f610a1SBarry Smith Not Collective 216c7f610a1SBarry Smith 217c7f610a1SBarry Smith Input Parameters: 218c7f610a1SBarry Smith . pc - the preconditioner context 219c7f610a1SBarry Smith 220c7f610a1SBarry Smith Output Parameter: 221c7f610a1SBarry Smith . pivot - the tolerance 222c7f610a1SBarry Smith 223c7f610a1SBarry Smith Level: intermediate 224c7f610a1SBarry Smith 225db781477SPatrick Sanan .seealso: `PCFactorSetZeroPivot()` 226c7f610a1SBarry Smith @*/ 227c7f610a1SBarry Smith PetscErrorCode PCFactorGetZeroPivot(PC pc,PetscReal *pivot) 228c7f610a1SBarry Smith { 229c7f610a1SBarry Smith PetscFunctionBegin; 230c7f610a1SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 231cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetZeroPivot_C",(PC,PetscReal*),(pc,pivot)); 232c7f610a1SBarry Smith PetscFunctionReturn(0); 233c7f610a1SBarry Smith } 234c7f610a1SBarry Smith 235c7f610a1SBarry Smith /*@ 236c7f610a1SBarry Smith PCFactorGetShiftAmount - Gets the tolerance used to define a zero privot 237c7f610a1SBarry Smith 238c7f610a1SBarry Smith Not Collective 239c7f610a1SBarry Smith 240c7f610a1SBarry Smith Input Parameters: 241c7f610a1SBarry Smith . pc - the preconditioner context 242c7f610a1SBarry Smith 243c7f610a1SBarry Smith Output Parameter: 244c7f610a1SBarry Smith . shift - how much to shift the diagonal entry 245c7f610a1SBarry Smith 246c7f610a1SBarry Smith Level: intermediate 247c7f610a1SBarry Smith 248db781477SPatrick Sanan .seealso: `PCFactorSetShiftAmount()`, `PCFactorSetShiftType()`, `PCFactorGetShiftType()` 249c7f610a1SBarry Smith @*/ 250c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftAmount(PC pc,PetscReal *shift) 251c7f610a1SBarry Smith { 252c7f610a1SBarry Smith PetscFunctionBegin; 253c7f610a1SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 254cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetShiftAmount_C",(PC,PetscReal*),(pc,shift)); 255c7f610a1SBarry Smith PetscFunctionReturn(0); 256c7f610a1SBarry Smith } 257c7f610a1SBarry Smith 258c7f610a1SBarry Smith /*@ 259c7f610a1SBarry Smith PCFactorGetShiftType - Gets the type of shift, if any, done when a zero pivot is detected 260c7f610a1SBarry Smith 261c7f610a1SBarry Smith Not Collective 262c7f610a1SBarry Smith 263c7f610a1SBarry Smith Input Parameters: 264c7f610a1SBarry Smith . pc - the preconditioner context 265c7f610a1SBarry Smith 266c7f610a1SBarry Smith Output Parameter: 267c7f610a1SBarry Smith . type - one of MAT_SHIFT_NONE, MAT_SHIFT_NONZERO, MAT_SHIFT_POSITIVE_DEFINITE, or MAT_SHIFT_INBLOCKS 268c7f610a1SBarry Smith 269c7f610a1SBarry Smith Level: intermediate 270c7f610a1SBarry Smith 271db781477SPatrick Sanan .seealso: `PCFactorSetShiftType()`, `MatFactorShiftType`, `PCFactorSetShiftAmount()`, `PCFactorGetShiftAmount()` 272c7f610a1SBarry Smith @*/ 273c7f610a1SBarry Smith PetscErrorCode PCFactorGetShiftType(PC pc,MatFactorShiftType *type) 274c7f610a1SBarry Smith { 275c7f610a1SBarry Smith PetscFunctionBegin; 276c7f610a1SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 277cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetShiftType_C",(PC,MatFactorShiftType*),(pc,type)); 278c7f610a1SBarry Smith PetscFunctionReturn(0); 279c7f610a1SBarry Smith } 280c7f610a1SBarry Smith 2812591b318SToby Isaac /*@ 2822591b318SToby Isaac PCFactorGetLevels - Gets the number of levels of fill to use. 2832591b318SToby Isaac 2842591b318SToby Isaac Logically Collective on PC 2852591b318SToby Isaac 2862591b318SToby Isaac Input Parameters: 2872591b318SToby Isaac . pc - the preconditioner context 2882591b318SToby Isaac 2892591b318SToby Isaac Output Parameter: 2902591b318SToby Isaac . levels - number of levels of fill 2912591b318SToby Isaac 2922591b318SToby Isaac Level: intermediate 2932591b318SToby Isaac 2942591b318SToby Isaac @*/ 2952591b318SToby Isaac PetscErrorCode PCFactorGetLevels(PC pc,PetscInt *levels) 2962591b318SToby Isaac { 2972591b318SToby Isaac PetscFunctionBegin; 2982591b318SToby Isaac PetscValidHeaderSpecific(pc,PC_CLASSID,1); 299cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetLevels_C",(PC,PetscInt*),(pc,levels)); 3002591b318SToby Isaac PetscFunctionReturn(0); 3012591b318SToby Isaac } 3022591b318SToby Isaac 30385317021SBarry Smith /*@ 30485317021SBarry Smith PCFactorSetLevels - Sets the number of levels of fill to use. 30585317021SBarry Smith 306ad4df100SBarry Smith Logically Collective on PC 30785317021SBarry Smith 30885317021SBarry Smith Input Parameters: 30985317021SBarry Smith + pc - the preconditioner context 31085317021SBarry Smith - levels - number of levels of fill 31185317021SBarry Smith 31285317021SBarry Smith Options Database Key: 31385317021SBarry Smith . -pc_factor_levels <levels> - Sets fill level 31485317021SBarry Smith 31585317021SBarry Smith Level: intermediate 31685317021SBarry Smith 31785317021SBarry Smith @*/ 3187087cfbeSBarry Smith PetscErrorCode PCFactorSetLevels(PC pc,PetscInt levels) 31985317021SBarry Smith { 32085317021SBarry Smith PetscFunctionBegin; 3210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3225f80ce2aSJacob Faibussowitsch PetscCheck(levels >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"negative levels"); 323c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,levels,2); 324cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetLevels_C",(PC,PetscInt),(pc,levels)); 32585317021SBarry Smith PetscFunctionReturn(0); 32685317021SBarry Smith } 32785317021SBarry Smith 32885317021SBarry Smith /*@ 32985317021SBarry Smith PCFactorSetAllowDiagonalFill - Causes all diagonal matrix entries to be 33085317021SBarry Smith treated as level 0 fill even if there is no non-zero location. 33185317021SBarry Smith 332ad4df100SBarry Smith Logically Collective on PC 33385317021SBarry Smith 33485317021SBarry Smith Input Parameters: 33585317021SBarry Smith + pc - the preconditioner context 33692e9c092SBarry Smith - flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 33785317021SBarry Smith 33885317021SBarry Smith Options Database Key: 33967b8a455SSatish Balay . -pc_factor_diagonal_fill <bool> - allow the diagonal fill 34085317021SBarry Smith 34185317021SBarry Smith Notes: 34285317021SBarry Smith Does not apply with 0 fill. 34385317021SBarry Smith 34485317021SBarry Smith Level: intermediate 34585317021SBarry Smith 346db781477SPatrick Sanan .seealso: `PCFactorGetAllowDiagonalFill()` 34785317021SBarry Smith @*/ 34892e9c092SBarry Smith PetscErrorCode PCFactorSetAllowDiagonalFill(PC pc,PetscBool flg) 34985317021SBarry Smith { 35085317021SBarry Smith PetscFunctionBegin; 3510700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 352cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetAllowDiagonalFill_C",(PC,PetscBool),(pc,flg)); 35392e9c092SBarry Smith PetscFunctionReturn(0); 35492e9c092SBarry Smith } 35592e9c092SBarry Smith 35692e9c092SBarry Smith /*@ 35792e9c092SBarry Smith PCFactorGetAllowDiagonalFill - Determines if all diagonal matrix entries are 35892e9c092SBarry Smith treated as level 0 fill even if there is no non-zero location. 35992e9c092SBarry Smith 36092e9c092SBarry Smith Logically Collective on PC 36192e9c092SBarry Smith 36292e9c092SBarry Smith Input Parameter: 36392e9c092SBarry Smith . pc - the preconditioner context 36492e9c092SBarry Smith 36592e9c092SBarry Smith Output Parameter: 36692e9c092SBarry Smith . flg - PETSC_TRUE to turn on, PETSC_FALSE to turn off 36792e9c092SBarry Smith 36892e9c092SBarry Smith Notes: 36992e9c092SBarry Smith Does not apply with 0 fill. 37092e9c092SBarry Smith 37192e9c092SBarry Smith Level: intermediate 37292e9c092SBarry Smith 373db781477SPatrick Sanan .seealso: `PCFactorSetAllowDiagonalFill()` 37492e9c092SBarry Smith @*/ 37592e9c092SBarry Smith PetscErrorCode PCFactorGetAllowDiagonalFill(PC pc,PetscBool *flg) 37692e9c092SBarry Smith { 37792e9c092SBarry Smith PetscFunctionBegin; 37892e9c092SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 379cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetAllowDiagonalFill_C",(PC,PetscBool*),(pc,flg)); 38085317021SBarry Smith PetscFunctionReturn(0); 38185317021SBarry Smith } 38285317021SBarry Smith 38385317021SBarry Smith /*@ 38485317021SBarry Smith PCFactorReorderForNonzeroDiagonal - reorders rows/columns of matrix to remove zeros from diagonal 38585317021SBarry Smith 386ad4df100SBarry Smith Logically Collective on PC 38785317021SBarry Smith 38885317021SBarry Smith Input Parameters: 38985317021SBarry Smith + pc - the preconditioner context 39085317021SBarry Smith - tol - diagonal entries smaller than this in absolute value are considered zero 39185317021SBarry Smith 39285317021SBarry Smith Options Database Key: 393147403d9SBarry Smith . -pc_factor_nonzeros_along_diagonal <tol> - perform the reordering with the given tolerance 39485317021SBarry Smith 39585317021SBarry Smith Level: intermediate 39685317021SBarry Smith 397db781477SPatrick Sanan .seealso: `PCFactorSetFill()`, `PCFactorSetShiftNonzero()`, `PCFactorSetZeroPivot()`, `MatReorderForNonzeroDiagonal()` 39885317021SBarry Smith @*/ 3997087cfbeSBarry Smith PetscErrorCode PCFactorReorderForNonzeroDiagonal(PC pc,PetscReal rtol) 40085317021SBarry Smith { 40185317021SBarry Smith PetscFunctionBegin; 4020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 403c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,rtol,2); 404cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorReorderForNonzeroDiagonal_C",(PC,PetscReal),(pc,rtol)); 40585317021SBarry Smith PetscFunctionReturn(0); 40685317021SBarry Smith } 40785317021SBarry Smith 408bf6011e8SBarry Smith /*@C 4093ca39a21SBarry Smith PCFactorSetMatSolverType - sets the software that is used to perform the factorization 41085317021SBarry Smith 411ad4df100SBarry Smith Logically Collective on PC 41285317021SBarry Smith 41385317021SBarry Smith Input Parameters: 41485317021SBarry Smith + pc - the preconditioner context 415f60c3dc2SHong Zhang - stype - for example, superlu, superlu_dist 41685317021SBarry Smith 41785317021SBarry Smith Options Database Key: 4183ca39a21SBarry Smith . -pc_factor_mat_solver_type <stype> - petsc, superlu, superlu_dist, mumps, cusparse 41985317021SBarry Smith 42085317021SBarry Smith Level: intermediate 42185317021SBarry Smith 42285317021SBarry Smith Note: 42385317021SBarry Smith By default this will use the PETSc factorization if it exists 42485317021SBarry Smith 425db781477SPatrick Sanan .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()` 42685317021SBarry Smith @*/ 427ea799195SBarry Smith PetscErrorCode PCFactorSetMatSolverType(PC pc,MatSolverType stype) 42885317021SBarry Smith { 42985317021SBarry Smith PetscFunctionBegin; 4300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 431cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetMatSolverType_C",(PC,MatSolverType),(pc,stype)); 43285317021SBarry Smith PetscFunctionReturn(0); 43385317021SBarry Smith } 43485317021SBarry Smith 435bf6011e8SBarry Smith /*@C 4363ca39a21SBarry Smith PCFactorGetMatSolverType - gets the software that is used to perform the factorization 4377112b564SBarry Smith 438c5eb9154SBarry Smith Not Collective 4397112b564SBarry Smith 4407112b564SBarry Smith Input Parameter: 4417112b564SBarry Smith . pc - the preconditioner context 4427112b564SBarry Smith 4437112b564SBarry Smith Output Parameter: 4440298fd71SBarry Smith . stype - for example, superlu, superlu_dist (NULL if the PC does not have a solver package) 4457112b564SBarry Smith 4467112b564SBarry Smith Level: intermediate 4477112b564SBarry Smith 448db781477SPatrick Sanan .seealso: `MatGetFactor()`, `MatSolverType`, `PCFactorGetMatSolverType()` 4497112b564SBarry Smith @*/ 450ea799195SBarry Smith PetscErrorCode PCFactorGetMatSolverType(PC pc,MatSolverType *stype) 4517112b564SBarry Smith { 4525f80ce2aSJacob Faibussowitsch PetscErrorCode (*f)(PC,MatSolverType*); 4537112b564SBarry Smith 4547112b564SBarry Smith PetscFunctionBegin; 4550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4565f80ce2aSJacob Faibussowitsch PetscValidPointer(stype,2); 4579566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",&f)); 4589566063dSJacob Faibussowitsch if (f) PetscCall((*f)(pc,stype)); 4595f80ce2aSJacob Faibussowitsch else *stype = NULL; 4607112b564SBarry Smith PetscFunctionReturn(0); 4617112b564SBarry Smith } 4627112b564SBarry Smith 46385317021SBarry Smith /*@ 46485317021SBarry Smith PCFactorSetFill - Indicate the amount of fill you expect in the factored matrix, 46585317021SBarry Smith fill = number nonzeros in factor/number nonzeros in original matrix. 46685317021SBarry Smith 467c5eb9154SBarry Smith Not Collective, each process can expect a different amount of fill 46885317021SBarry Smith 46985317021SBarry Smith Input Parameters: 47085317021SBarry Smith + pc - the preconditioner context 47185317021SBarry Smith - fill - amount of expected fill 47285317021SBarry Smith 47385317021SBarry Smith Options Database Key: 47485317021SBarry Smith . -pc_factor_fill <fill> - Sets fill amount 47585317021SBarry Smith 47685317021SBarry Smith Level: intermediate 47785317021SBarry Smith 47885317021SBarry Smith Note: 47985317021SBarry Smith For sparse matrix factorizations it is difficult to predict how much 48085317021SBarry Smith fill to expect. By running with the option -info PETSc will print the 48185317021SBarry Smith actual amount of fill used; allowing you to set the value accurately for 48285317021SBarry Smith future runs. Default PETSc uses a value of 5.0 48385317021SBarry Smith 48401a79839SBarry Smith This parameter has NOTHING to do with the levels-of-fill of ILU(). That is set with PCFactorSetLevels() or -pc_factor_levels. 48501a79839SBarry Smith 48685317021SBarry Smith @*/ 4877087cfbeSBarry Smith PetscErrorCode PCFactorSetFill(PC pc,PetscReal fill) 48885317021SBarry Smith { 48985317021SBarry Smith PetscFunctionBegin; 4900700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4915f80ce2aSJacob Faibussowitsch PetscCheck(fill >= 1.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Fill factor cannot be less then 1.0"); 492cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetFill_C",(PC,PetscReal),(pc,fill)); 49385317021SBarry Smith PetscFunctionReturn(0); 49485317021SBarry Smith } 49585317021SBarry Smith 49685317021SBarry Smith /*@ 49785317021SBarry Smith PCFactorSetUseInPlace - Tells the system to do an in-place factorization. 49885317021SBarry Smith For dense matrices, this enables the solution of much larger problems. 49985317021SBarry Smith For sparse matrices the factorization cannot be done truly in-place 50085317021SBarry Smith so this does not save memory during the factorization, but after the matrix 50185317021SBarry Smith is factored, the original unfactored matrix is freed, thus recovering that 502ec5066bdSBarry Smith space. For ICC(0) and ILU(0) with the default natural ordering the factorization is done efficiently in-place. 50385317021SBarry Smith 504ad4df100SBarry Smith Logically Collective on PC 50585317021SBarry Smith 50685317021SBarry Smith Input Parameters: 5078e37d05fSBarry Smith + pc - the preconditioner context 5088e37d05fSBarry Smith - flg - PETSC_TRUE to enable, PETSC_FALSE to disable 50985317021SBarry Smith 51085317021SBarry Smith Options Database Key: 5118e37d05fSBarry Smith . -pc_factor_in_place <true,false>- Activate/deactivate in-place factorization 51285317021SBarry Smith 51385317021SBarry Smith Notes: 51485317021SBarry Smith PCFactorSetUseInplace() can only be used with the KSP method KSPPREONLY or when 51585317021SBarry Smith a different matrix is provided for the multiply and the preconditioner in 51685317021SBarry Smith a call to KSPSetOperators(). 51785317021SBarry Smith This is because the Krylov space methods require an application of the 51885317021SBarry Smith matrix multiplication, which is not possible here because the matrix has 51985317021SBarry Smith been factored in-place, replacing the original matrix. 52085317021SBarry Smith 52185317021SBarry Smith Level: intermediate 52285317021SBarry Smith 523db781477SPatrick Sanan .seealso: `PCFactorGetUseInPlace()` 52485317021SBarry Smith @*/ 5258e37d05fSBarry Smith PetscErrorCode PCFactorSetUseInPlace(PC pc,PetscBool flg) 52685317021SBarry Smith { 52785317021SBarry Smith PetscFunctionBegin; 5280700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 529cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetUseInPlace_C",(PC,PetscBool),(pc,flg)); 5308e37d05fSBarry Smith PetscFunctionReturn(0); 5318e37d05fSBarry Smith } 5328e37d05fSBarry Smith 5338e37d05fSBarry Smith /*@ 5348e37d05fSBarry Smith PCFactorGetUseInPlace - Determines if an in-place factorization is being used. 5358e37d05fSBarry Smith 5368e37d05fSBarry Smith Logically Collective on PC 5378e37d05fSBarry Smith 5388e37d05fSBarry Smith Input Parameter: 5398e37d05fSBarry Smith . pc - the preconditioner context 5408e37d05fSBarry Smith 5418e37d05fSBarry Smith Output Parameter: 5428e37d05fSBarry Smith . flg - PETSC_TRUE to enable, PETSC_FALSE to disable 5438e37d05fSBarry Smith 5448e37d05fSBarry Smith Level: intermediate 5458e37d05fSBarry Smith 546db781477SPatrick Sanan .seealso: `PCFactorSetUseInPlace()` 5478e37d05fSBarry Smith @*/ 5488e37d05fSBarry Smith PetscErrorCode PCFactorGetUseInPlace(PC pc,PetscBool *flg) 5498e37d05fSBarry Smith { 5508e37d05fSBarry Smith PetscFunctionBegin; 5518e37d05fSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 552cac4c232SBarry Smith PetscUseMethod(pc,"PCFactorGetUseInPlace_C",(PC,PetscBool*),(pc,flg)); 55385317021SBarry Smith PetscFunctionReturn(0); 55485317021SBarry Smith } 55585317021SBarry Smith 55685317021SBarry Smith /*@C 55785317021SBarry Smith PCFactorSetMatOrderingType - Sets the ordering routine (to reduce fill) to 5582c7c0729SBarry Smith be used in the LU, ILU, Cholesky, and ICC factorizations. 55985317021SBarry Smith 560ad4df100SBarry Smith Logically Collective on PC 56185317021SBarry Smith 56285317021SBarry Smith Input Parameters: 56385317021SBarry Smith + pc - the preconditioner context 5642692d6eeSBarry Smith - ordering - the matrix ordering name, for example, MATORDERINGND or MATORDERINGRCM 56585317021SBarry Smith 56685317021SBarry Smith Options Database Key: 5672c7c0729SBarry Smith . -pc_factor_mat_ordering_type <nd,rcm,...,external> - Sets ordering routine 56885317021SBarry Smith 56985317021SBarry Smith Level: intermediate 57085317021SBarry Smith 57195452b02SPatrick Sanan Notes: 5724ac6704cSBarry Smith Nested dissection is used by default for some of PETSc's sparse matrix formats 57385317021SBarry Smith 5749bd791bbSBarry Smith For Cholesky and ICC and the SBAIJ format the only reordering available is natural since only the upper half of the matrix is stored 5759bd791bbSBarry Smith and reordering this matrix is very expensive. 5769bd791bbSBarry Smith 5774ac6704cSBarry Smith You can use a SeqAIJ matrix with Cholesky and ICC and use any ordering. 5789bd791bbSBarry Smith 5794ac6704cSBarry Smith MATORDERINGEXTERNAL means PETSc will not compute an ordering and the package will use its own ordering, usable with MATSOLVERCHOLMOD, MATSOLVERUMFPACK, and others. 5802c7c0729SBarry Smith 581db781477SPatrick Sanan .seealso: `MatOrderingType` 58285317021SBarry Smith 58385317021SBarry Smith @*/ 58419fd82e9SBarry Smith PetscErrorCode PCFactorSetMatOrderingType(PC pc,MatOrderingType ordering) 58585317021SBarry Smith { 58685317021SBarry Smith PetscFunctionBegin; 587c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 588cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetMatOrderingType_C",(PC,MatOrderingType),(pc,ordering)); 58985317021SBarry Smith PetscFunctionReturn(0); 59085317021SBarry Smith } 59185317021SBarry Smith 59285317021SBarry Smith /*@ 5938ff23777SHong Zhang PCFactorSetColumnPivot - Determines when column pivoting is done during matrix factorization. 59485317021SBarry Smith For PETSc dense matrices column pivoting is always done, for PETSc sparse matrices 595e3c5b3baSBarry Smith it is never done. For the MATLAB and SuperLU factorization this is used. 59685317021SBarry Smith 597ad4df100SBarry Smith Logically Collective on PC 59885317021SBarry Smith 59985317021SBarry Smith Input Parameters: 60085317021SBarry Smith + pc - the preconditioner context 60185317021SBarry Smith - dtcol - 0.0 implies no pivoting, 1.0 complete pivoting (slower, requires more memory but more stable) 60285317021SBarry Smith 60385317021SBarry Smith Options Database Key: 604147403d9SBarry Smith . -pc_factor_pivoting <dtcol> - perform the pivoting with the given tolerance 60585317021SBarry Smith 60685317021SBarry Smith Level: intermediate 60785317021SBarry Smith 608db781477SPatrick Sanan .seealso: `PCILUSetMatOrdering()`, `PCFactorSetPivotInBlocks()` 60985317021SBarry Smith @*/ 6107087cfbeSBarry Smith PetscErrorCode PCFactorSetColumnPivot(PC pc,PetscReal dtcol) 61185317021SBarry Smith { 61285317021SBarry Smith PetscFunctionBegin; 613c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 614c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,dtcol,2); 615cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetColumnPivot_C",(PC,PetscReal),(pc,dtcol)); 61685317021SBarry Smith PetscFunctionReturn(0); 61785317021SBarry Smith } 61885317021SBarry Smith 61985317021SBarry Smith /*@ 62085317021SBarry Smith PCFactorSetPivotInBlocks - Determines if pivoting is done while factoring each block 62185317021SBarry Smith with BAIJ or SBAIJ matrices 62285317021SBarry Smith 623ad4df100SBarry Smith Logically Collective on PC 62485317021SBarry Smith 62585317021SBarry Smith Input Parameters: 62685317021SBarry Smith + pc - the preconditioner context 62785317021SBarry Smith - pivot - PETSC_TRUE or PETSC_FALSE 62885317021SBarry Smith 62985317021SBarry Smith Options Database Key: 63067b8a455SSatish Balay . -pc_factor_pivot_in_blocks <true,false> - Pivot inside matrix dense blocks for BAIJ and SBAIJ 63185317021SBarry Smith 63285317021SBarry Smith Level: intermediate 63385317021SBarry Smith 634db781477SPatrick Sanan .seealso: `PCILUSetMatOrdering()`, `PCFactorSetColumnPivot()` 63585317021SBarry Smith @*/ 6367087cfbeSBarry Smith PetscErrorCode PCFactorSetPivotInBlocks(PC pc,PetscBool pivot) 63785317021SBarry Smith { 63885317021SBarry Smith PetscFunctionBegin; 639c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 640acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,pivot,2); 641cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetPivotInBlocks_C",(PC,PetscBool),(pc,pivot)); 64285317021SBarry Smith PetscFunctionReturn(0); 64385317021SBarry Smith } 64485317021SBarry Smith 64585317021SBarry Smith /*@ 646288e7d53SBarry Smith PCFactorSetReuseFill - When matrices with different nonzero structure are factored, 64785317021SBarry Smith this causes later ones to use the fill ratio computed in the initial factorization. 64885317021SBarry Smith 649ad4df100SBarry Smith Logically Collective on PC 65085317021SBarry Smith 65185317021SBarry Smith Input Parameters: 65285317021SBarry Smith + pc - the preconditioner context 65385317021SBarry Smith - flag - PETSC_TRUE to reuse else PETSC_FALSE 65485317021SBarry Smith 65585317021SBarry Smith Options Database Key: 65685317021SBarry Smith . -pc_factor_reuse_fill - Activates PCFactorSetReuseFill() 65785317021SBarry Smith 65885317021SBarry Smith Level: intermediate 65985317021SBarry Smith 660db781477SPatrick Sanan .seealso: `PCFactorSetReuseOrdering()` 66185317021SBarry Smith @*/ 6627087cfbeSBarry Smith PetscErrorCode PCFactorSetReuseFill(PC pc,PetscBool flag) 66385317021SBarry Smith { 66485317021SBarry Smith PetscFunctionBegin; 665064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(pc,PC_CLASSID,1); 666acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,flag,2); 667cac4c232SBarry Smith PetscTryMethod(pc,"PCFactorSetReuseFill_C",(PC,PetscBool),(pc,flag)); 66885317021SBarry Smith PetscFunctionReturn(0); 66985317021SBarry Smith } 6703d1c1ea0SBarry Smith 6714ac6704cSBarry Smith PetscErrorCode PCFactorInitialize(PC pc,MatFactorType ftype) 6723d1c1ea0SBarry Smith { 6733d1c1ea0SBarry Smith PC_Factor *fact = (PC_Factor*)pc->data; 6743d1c1ea0SBarry Smith 6753d1c1ea0SBarry Smith PetscFunctionBegin; 6769566063dSJacob Faibussowitsch PetscCall(MatFactorInfoInitialize(&fact->info)); 6774ac6704cSBarry Smith fact->factortype = ftype; 6783d1c1ea0SBarry Smith fact->info.shifttype = (PetscReal)MAT_SHIFT_NONE; 6793d1c1ea0SBarry Smith fact->info.shiftamount = 100.0*PETSC_MACHINE_EPSILON; 6803d1c1ea0SBarry Smith fact->info.zeropivot = 100.0*PETSC_MACHINE_EPSILON; 6813d1c1ea0SBarry Smith fact->info.pivotinblocks = 1.0; 6823d1c1ea0SBarry Smith pc->ops->getfactoredmatrix = PCFactorGetMatrix_Factor; 6833d1c1ea0SBarry Smith 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",PCFactorSetZeroPivot_Factor)); 6859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",PCFactorGetZeroPivot_Factor)); 6869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",PCFactorSetShiftType_Factor)); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",PCFactorGetShiftType_Factor)); 6889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",PCFactorSetShiftAmount_Factor)); 6899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",PCFactorGetShiftAmount_Factor)); 6909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",PCFactorGetMatSolverType_Factor)); 6919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",PCFactorSetMatSolverType_Factor)); 6929566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",PCFactorSetUpMatSolverType_Factor)); 6939566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",PCFactorSetFill_Factor)); 6949566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",PCFactorSetMatOrderingType_Factor)); 6959566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",PCFactorSetLevels_Factor)); 6969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",PCFactorGetLevels_Factor)); 6979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",PCFactorSetAllowDiagonalFill_Factor)); 6989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",PCFactorGetAllowDiagonalFill_Factor)); 6999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",PCFactorSetPivotInBlocks_Factor)); 7009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",PCFactorSetUseInPlace_Factor)); 7019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",PCFactorGetUseInPlace_Factor)); 7029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",PCFactorSetReuseOrdering_Factor)); 7039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",PCFactorSetReuseFill_Factor)); 7043d1c1ea0SBarry Smith PetscFunctionReturn(0); 7053d1c1ea0SBarry Smith } 7062e956fe4SStefano Zampini 7072e956fe4SStefano Zampini PetscErrorCode PCFactorClearComposedFunctions(PC pc) 7082e956fe4SStefano Zampini { 7092e956fe4SStefano Zampini PetscFunctionBegin; 7102e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetZeroPivot_C",NULL)); 7112e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetZeroPivot_C",NULL)); 7122e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftType_C",NULL)); 7132e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftType_C",NULL)); 7142e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetShiftAmount_C",NULL)); 7152e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetShiftAmount_C",NULL)); 7162e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetMatSolverType_C",NULL)); 7172e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatSolverType_C",NULL)); 7182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUpMatSolverType_C",NULL)); 7192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetFill_C",NULL)); 7202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetMatOrderingType_C",NULL)); 7212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetLevels_C",NULL)); 7222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetLevels_C",NULL)); 7232e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetAllowDiagonalFill_C",NULL)); 7242e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetAllowDiagonalFill_C",NULL)); 7252e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetPivotInBlocks_C",NULL)); 7262e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetUseInPlace_C",NULL)); 7272e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorGetUseInPlace_C",NULL)); 7282e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseOrdering_C",NULL)); 7292e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetReuseFill_C",NULL)); 7302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorReorderForNonzeroDiagonal_C",NULL)); 7312e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCFactorSetDropTolerance_C",NULL)); 7322e956fe4SStefano Zampini PetscFunctionReturn(0); 7332e956fe4SStefano Zampini } 734