14b9ad928SBarry Smith /* 24b9ad928SBarry Smith The PC (preconditioner) interface routines, callable by users. 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscksp.h" I*/ 51e25c274SJed Brown #include <petscdm.h> 64b9ad928SBarry Smith 74b9ad928SBarry Smith /* Logging support */ 80700a824SBarry Smith PetscClassId PC_CLASSID; 9c677e75fSPierre Jolivet PetscLogEvent PC_SetUp, PC_SetUpOnBlocks, PC_Apply, PC_MatApply, PC_ApplyCoarse, PC_ApplyMultiple, PC_ApplySymmetricLeft; 10011ea8aeSBarry Smith PetscLogEvent PC_ApplySymmetricRight, PC_ModifySubMatrices, PC_ApplyOnBlocks, PC_ApplyTransposeOnBlocks; 11ab83eea4SMatthew G. Knepley PetscInt PetscMGLevelId; 120316ec64SBarry Smith PetscLogStage PCMPIStage; 134b9ad928SBarry Smith 14d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDefaultType_Private(PC pc, const char *type[]) 15d71ae5a4SJacob Faibussowitsch { 162e70eadcSBarry Smith PetscMPIInt size; 17da74c943SJose E. Roman PetscBool hasop, flg1, flg2, set, flg3, isnormal; 18566e8bf2SBarry Smith 19566e8bf2SBarry Smith PetscFunctionBegin; 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 21566e8bf2SBarry Smith if (pc->pmat) { 229566063dSJacob Faibussowitsch PetscCall(MatHasOperation(pc->pmat, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 23566e8bf2SBarry Smith if (size == 1) { 249566063dSJacob Faibussowitsch PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ICC, &flg1)); 259566063dSJacob Faibussowitsch PetscCall(MatGetFactorAvailable(pc->pmat, "petsc", MAT_FACTOR_ILU, &flg2)); 269566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &flg3)); 27da74c943SJose E. Roman PetscCall(PetscObjectTypeCompareAny((PetscObject)pc->pmat, &isnormal, MATNORMAL, MATNORMALHERMITIAN, NULL)); 28ba8a7149SBarry Smith if (flg1 && (!flg2 || (set && flg3))) { 29566e8bf2SBarry Smith *type = PCICC; 30566e8bf2SBarry Smith } else if (flg2) { 31566e8bf2SBarry Smith *type = PCILU; 32da74c943SJose E. Roman } else if (isnormal) { 33da74c943SJose E. Roman *type = PCNONE; 34976e8c5aSLisandro Dalcin } else if (hasop) { /* likely is a parallel matrix run on one processor */ 35566e8bf2SBarry Smith *type = PCBJACOBI; 36566e8bf2SBarry Smith } else { 37566e8bf2SBarry Smith *type = PCNONE; 38566e8bf2SBarry Smith } 39566e8bf2SBarry Smith } else { 40976e8c5aSLisandro Dalcin if (hasop) { 41566e8bf2SBarry Smith *type = PCBJACOBI; 42566e8bf2SBarry Smith } else { 43566e8bf2SBarry Smith *type = PCNONE; 44566e8bf2SBarry Smith } 45566e8bf2SBarry Smith } 46566e8bf2SBarry Smith } else { 47566e8bf2SBarry Smith if (size == 1) { 48566e8bf2SBarry Smith *type = PCILU; 49566e8bf2SBarry Smith } else { 50566e8bf2SBarry Smith *type = PCBJACOBI; 51566e8bf2SBarry Smith } 52566e8bf2SBarry Smith } 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 54566e8bf2SBarry Smith } 554b9ad928SBarry Smith 5688584be7SBarry Smith /*@ 57f1580f4eSBarry Smith PCReset - Resets a PC context to the pcsetupcalled = 0 state and removes any allocated `Vec`s and `Mat`s 5888584be7SBarry Smith 59c3339decSBarry Smith Collective 6088584be7SBarry Smith 6188584be7SBarry Smith Input Parameter: 6288584be7SBarry Smith . pc - the preconditioner context 6388584be7SBarry Smith 6488584be7SBarry Smith Level: developer 6588584be7SBarry Smith 66f1580f4eSBarry Smith Note: 67f1580f4eSBarry Smith This allows a `PC` to be reused for a different sized linear system but using the same options that have been previously set in the PC 6888584be7SBarry Smith 69f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()` 7088584be7SBarry Smith @*/ 71d71ae5a4SJacob Faibussowitsch PetscErrorCode PCReset(PC pc) 72d71ae5a4SJacob Faibussowitsch { 7388584be7SBarry Smith PetscFunctionBegin; 7488584be7SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 75dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, reset); 769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleright)); 779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 802fa5cd67SKarl Rupp 810ce6c5a2SBarry Smith pc->setupcalled = 0; 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8388584be7SBarry Smith } 8488584be7SBarry Smith 851fb7b255SJunchao Zhang /*@C 86f1580f4eSBarry Smith PCDestroy - Destroys `PC` context that was created with `PCCreate()`. 874b9ad928SBarry Smith 88c3339decSBarry Smith Collective 894b9ad928SBarry Smith 904b9ad928SBarry Smith Input Parameter: 914b9ad928SBarry Smith . pc - the preconditioner context 924b9ad928SBarry Smith 934b9ad928SBarry Smith Level: developer 944b9ad928SBarry Smith 95f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()` 964b9ad928SBarry Smith @*/ 97d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDestroy(PC *pc) 98d71ae5a4SJacob Faibussowitsch { 994b9ad928SBarry Smith PetscFunctionBegin; 1003ba16761SJacob Faibussowitsch if (!*pc) PetscFunctionReturn(PETSC_SUCCESS); 1016bf464f9SBarry Smith PetscValidHeaderSpecific((*pc), PC_CLASSID, 1); 1029371c9d4SSatish Balay if (--((PetscObject)(*pc))->refct > 0) { 1039371c9d4SSatish Balay *pc = NULL; 1043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1059371c9d4SSatish Balay } 1064b9ad928SBarry Smith 1079566063dSJacob Faibussowitsch PetscCall(PCReset(*pc)); 108241cb3abSLisandro Dalcin 109e04113cfSBarry Smith /* if memory was published with SAWs then destroy it */ 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsViewOff((PetscObject)*pc)); 111dbbe0bcdSBarry Smith PetscTryTypeMethod((*pc), destroy); 1129566063dSJacob Faibussowitsch PetscCall(DMDestroy(&(*pc)->dm)); 1139566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(pc)); 1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1154b9ad928SBarry Smith } 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith /*@C 118c5eb9154SBarry Smith PCGetDiagonalScale - Indicates if the preconditioner applies an additional left and right 1194b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1204b9ad928SBarry Smith 121c3339decSBarry Smith Logically Collective 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith Input Parameter: 1244b9ad928SBarry Smith . pc - the preconditioner context 1254b9ad928SBarry Smith 1264b9ad928SBarry Smith Output Parameter: 127f1580f4eSBarry Smith . flag - `PETSC_TRUE` if it applies the scaling 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith Level: developer 1304b9ad928SBarry Smith 131f1580f4eSBarry Smith Note: 132f1580f4eSBarry Smith If this returns `PETSC_TRUE` then the system solved via the Krylov method is 133f1580f4eSBarry Smith .vb 134f1580f4eSBarry Smith D M A D^{-1} y = D M b for left preconditioning or 135f1580f4eSBarry Smith D A M D^{-1} z = D b for right preconditioning 136f1580f4eSBarry Smith .ve 1374b9ad928SBarry Smith 138f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCSetDiagonalScale()` 1394b9ad928SBarry Smith @*/ 140d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetDiagonalScale(PC pc, PetscBool *flag) 141d71ae5a4SJacob Faibussowitsch { 1424b9ad928SBarry Smith PetscFunctionBegin; 1430700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1444f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1454b9ad928SBarry Smith *flag = pc->diagonalscale; 1463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1474b9ad928SBarry Smith } 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith /*@ 150c5eb9154SBarry Smith PCSetDiagonalScale - Indicates the left scaling to use to apply an additional left and right 1514b9ad928SBarry Smith scaling as needed by certain time-stepping codes. 1524b9ad928SBarry Smith 153c3339decSBarry Smith Logically Collective 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith Input Parameters: 1564b9ad928SBarry Smith + pc - the preconditioner context 1574b9ad928SBarry Smith - s - scaling vector 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith Level: intermediate 1604b9ad928SBarry Smith 16195452b02SPatrick Sanan Notes: 16295452b02SPatrick Sanan The system solved via the Krylov method is 163f1580f4eSBarry Smith .vb 164f1580f4eSBarry Smith D M A D^{-1} y = D M b for left preconditioning or 165f1580f4eSBarry Smith D A M D^{-1} z = D b for right preconditioning 166f1580f4eSBarry Smith .ve 1674b9ad928SBarry Smith 168f1580f4eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}. 1694b9ad928SBarry Smith 170db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleRight()`, `PCGetDiagonalScale()` 1714b9ad928SBarry Smith @*/ 172d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetDiagonalScale(PC pc, Vec s) 173d71ae5a4SJacob Faibussowitsch { 1744b9ad928SBarry Smith PetscFunctionBegin; 1750700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1760700a824SBarry Smith PetscValidHeaderSpecific(s, VEC_CLASSID, 2); 1774b9ad928SBarry Smith pc->diagonalscale = PETSC_TRUE; 1782fa5cd67SKarl Rupp 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)s)); 1809566063dSJacob Faibussowitsch PetscCall(VecDestroy(&pc->diagonalscaleleft)); 1812fa5cd67SKarl Rupp 1824b9ad928SBarry Smith pc->diagonalscaleleft = s; 1832fa5cd67SKarl Rupp 1849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(s, &pc->diagonalscaleright)); 1859566063dSJacob Faibussowitsch PetscCall(VecCopy(s, pc->diagonalscaleright)); 1869566063dSJacob Faibussowitsch PetscCall(VecReciprocal(pc->diagonalscaleright)); 1873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1884b9ad928SBarry Smith } 1894b9ad928SBarry Smith 190bc08b0f1SBarry Smith /*@ 191c5eb9154SBarry Smith PCDiagonalScaleLeft - Scales a vector by the left scaling as needed by certain time-stepping codes. 1924b9ad928SBarry Smith 193c3339decSBarry Smith Logically Collective 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith Input Parameters: 1964b9ad928SBarry Smith + pc - the preconditioner context 1974b9ad928SBarry Smith . in - input vector 198a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 1994b9ad928SBarry Smith 2004b9ad928SBarry Smith Level: intermediate 2014b9ad928SBarry Smith 20295452b02SPatrick Sanan Notes: 20395452b02SPatrick Sanan The system solved via the Krylov method is 204f1580f4eSBarry Smith .vb 205f1580f4eSBarry Smith D M A D^{-1} y = D M b for left preconditioning or 206f1580f4eSBarry Smith D A M D^{-1} z = D b for right preconditioning 207f1580f4eSBarry Smith .ve 2084b9ad928SBarry Smith 209f1580f4eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}. 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith If diagonal scaling is turned off and in is not out then in is copied to out 2124b9ad928SBarry Smith 213db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleSet()`, `PCDiagonalScaleRight()`, `PCDiagonalScale()` 2144b9ad928SBarry Smith @*/ 215d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleLeft(PC pc, Vec in, Vec out) 216d71ae5a4SJacob Faibussowitsch { 2174b9ad928SBarry Smith PetscFunctionBegin; 2180700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2190700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2200700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2214b9ad928SBarry Smith if (pc->diagonalscale) { 2229566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleleft, in)); 2234b9ad928SBarry Smith } else if (in != out) { 2249566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 2254b9ad928SBarry Smith } 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2274b9ad928SBarry Smith } 2284b9ad928SBarry Smith 229bc08b0f1SBarry Smith /*@ 2304b9ad928SBarry Smith PCDiagonalScaleRight - Scales a vector by the right scaling as needed by certain time-stepping codes. 2314b9ad928SBarry Smith 232c3339decSBarry Smith Logically Collective 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith Input Parameters: 2354b9ad928SBarry Smith + pc - the preconditioner context 2364b9ad928SBarry Smith . in - input vector 237a2b725a8SWilliam Gropp - out - scaled vector (maybe the same as in) 2384b9ad928SBarry Smith 2394b9ad928SBarry Smith Level: intermediate 2404b9ad928SBarry Smith 24195452b02SPatrick Sanan Notes: 24295452b02SPatrick Sanan The system solved via the Krylov method is 243f1580f4eSBarry Smith .vb 244f1580f4eSBarry Smith D M A D^{-1} y = D M b for left preconditioning or 245f1580f4eSBarry Smith D A M D^{-1} z = D b for right preconditioning 246f1580f4eSBarry Smith .ve 2474b9ad928SBarry Smith 248f1580f4eSBarry Smith `PCDiagonalScaleLeft()` scales a vector by D. `PCDiagonalScaleRight()` scales a vector by D^{-1}. 2494b9ad928SBarry Smith 2504b9ad928SBarry Smith If diagonal scaling is turned off and in is not out then in is copied to out 2514b9ad928SBarry Smith 252db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetUp()`, `PCDiagonalScaleLeft()`, `PCDiagonalScaleSet()`, `PCDiagonalScale()` 2534b9ad928SBarry Smith @*/ 254d71ae5a4SJacob Faibussowitsch PetscErrorCode PCDiagonalScaleRight(PC pc, Vec in, Vec out) 255d71ae5a4SJacob Faibussowitsch { 2564b9ad928SBarry Smith PetscFunctionBegin; 2570700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 2580700a824SBarry Smith PetscValidHeaderSpecific(in, VEC_CLASSID, 2); 2590700a824SBarry Smith PetscValidHeaderSpecific(out, VEC_CLASSID, 3); 2604b9ad928SBarry Smith if (pc->diagonalscale) { 2619566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(out, pc->diagonalscaleright, in)); 2624b9ad928SBarry Smith } else if (in != out) { 2639566063dSJacob Faibussowitsch PetscCall(VecCopy(in, out)); 2644b9ad928SBarry Smith } 2653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2664b9ad928SBarry Smith } 2674b9ad928SBarry Smith 26849517cdeSBarry Smith /*@ 26949517cdeSBarry Smith PCSetUseAmat - Sets a flag to indicate that when the preconditioner needs to apply (part of) the 270f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 271f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 27249517cdeSBarry Smith 273c3339decSBarry Smith Logically Collective 27449517cdeSBarry Smith 27549517cdeSBarry Smith Input Parameters: 27649517cdeSBarry Smith + pc - the preconditioner context 277f1580f4eSBarry Smith - flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 27849517cdeSBarry Smith 27949517cdeSBarry Smith Options Database Key: 280147403d9SBarry Smith . -pc_use_amat <true,false> - use the amat to apply the operator 28149517cdeSBarry Smith 28220f4b53cSBarry Smith Level: intermediate 28320f4b53cSBarry Smith 284f1580f4eSBarry Smith Note: 28549517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 28649517cdeSBarry Smith preconditioner are identical, this routine is does nothing. 28749517cdeSBarry Smith 288f1580f4eSBarry Smith .seealso: `PC`, `PCGetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE` 28949517cdeSBarry Smith @*/ 290d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUseAmat(PC pc, PetscBool flg) 291d71ae5a4SJacob Faibussowitsch { 29249517cdeSBarry Smith PetscFunctionBegin; 29349517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 29449517cdeSBarry Smith pc->useAmat = flg; 2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 29649517cdeSBarry Smith } 29749517cdeSBarry Smith 298422a814eSBarry Smith /*@ 299f1580f4eSBarry Smith PCSetErrorIfFailure - Causes `PC` to generate an error if a FPE, for example a zero pivot, is detected. 300422a814eSBarry Smith 301c3339decSBarry Smith Logically Collective 302422a814eSBarry Smith 303422a814eSBarry Smith Input Parameters: 304422a814eSBarry Smith + pc - iterative context obtained from PCCreate() 305f1580f4eSBarry Smith - flg - `PETSC_TRUE` indicates you want the error generated 306422a814eSBarry Smith 307422a814eSBarry Smith Level: advanced 308422a814eSBarry Smith 309422a814eSBarry Smith Notes: 310f1580f4eSBarry Smith Normally PETSc continues if a linear solver fails due to a failed setup of a preconditioner, you can call `KSPGetConvergedReason()` after a `KSPSolve()` 311422a814eSBarry Smith to determine if it has converged or failed. Or use -ksp_error_if_not_converged to cause the program to terminate as soon as lack of convergence is 312422a814eSBarry Smith detected. 313422a814eSBarry Smith 314422a814eSBarry Smith This is propagated into KSPs used by this PC, which then propagate it into PCs used by those KSPs 315422a814eSBarry Smith 316f1580f4eSBarry Smith .seealso: `PC`, `KSPSetErrorIfNotConverged()`, `PCGetInitialGuessNonzero()`, `PCSetInitialGuessKnoll()`, `PCGetInitialGuessKnoll()` 317422a814eSBarry Smith @*/ 318d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetErrorIfFailure(PC pc, PetscBool flg) 319d71ae5a4SJacob Faibussowitsch { 320422a814eSBarry Smith PetscFunctionBegin; 321422a814eSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 322422a814eSBarry Smith PetscValidLogicalCollectiveBool(pc, flg, 2); 323422a814eSBarry Smith pc->erroriffailure = flg; 3243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 325422a814eSBarry Smith } 326422a814eSBarry Smith 32749517cdeSBarry Smith /*@ 32849517cdeSBarry Smith PCGetUseAmat - Gets a flag to indicate that when the preconditioner needs to apply (part of) the 329f1580f4eSBarry Smith operator during the preconditioning process it applies the Amat provided to `TSSetRHSJacobian()`, 330f1580f4eSBarry Smith `TSSetIJacobian()`, `SNESSetJacobian()`, `KSPSetOperators()` or `PCSetOperators()` not the Pmat. 33149517cdeSBarry Smith 332c3339decSBarry Smith Logically Collective 33349517cdeSBarry Smith 33449517cdeSBarry Smith Input Parameter: 33549517cdeSBarry Smith . pc - the preconditioner context 33649517cdeSBarry Smith 33749517cdeSBarry Smith Output Parameter: 338f1580f4eSBarry Smith . flg - `PETSC_TRUE` to use the Amat, `PETSC_FALSE` to use the Pmat (default is false) 33949517cdeSBarry Smith 34020f4b53cSBarry Smith Level: intermediate 34120f4b53cSBarry Smith 342f1580f4eSBarry Smith Note: 34349517cdeSBarry Smith For the common case in which the linear system matrix and the matrix used to construct the 34449517cdeSBarry Smith preconditioner are identical, this routine is does nothing. 34549517cdeSBarry Smith 346f1580f4eSBarry Smith .seealso: `PC`, `PCSetUseAmat()`, `PCBJACOBI`, `PGMG`, `PCFIELDSPLIT`, `PCCOMPOSITE` 34749517cdeSBarry Smith @*/ 348d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetUseAmat(PC pc, PetscBool *flg) 349d71ae5a4SJacob Faibussowitsch { 35049517cdeSBarry Smith PetscFunctionBegin; 35149517cdeSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 35249517cdeSBarry Smith *flg = pc->useAmat; 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35449517cdeSBarry Smith } 35549517cdeSBarry Smith 356f39d8e23SSatish Balay /*@ 3573821be0aSBarry Smith PCSetKSPNestLevel - sets the amount of nesting the `KSP` that contains this `PC` has 3583821be0aSBarry Smith 3593821be0aSBarry Smith Collective 3603821be0aSBarry Smith 3613821be0aSBarry Smith Input Parameters: 3623821be0aSBarry Smith + pc - the `PC` 3633821be0aSBarry Smith - level - the nest level 3643821be0aSBarry Smith 3653821be0aSBarry Smith Level: developer 3663821be0aSBarry Smith 3673821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPGetNestLevel()`, `PCGetKSPNestLevel()`, `KSPSetNestLevel()` 3683821be0aSBarry Smith @*/ 3693821be0aSBarry Smith PetscErrorCode PCSetKSPNestLevel(PC pc, PetscInt level) 3703821be0aSBarry Smith { 3713821be0aSBarry Smith PetscFunctionBegin; 3727a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3737a99bfcaSBarry Smith PetscValidLogicalCollectiveInt(pc, level, 2); 3743821be0aSBarry Smith pc->kspnestlevel = level; 3753821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 3763821be0aSBarry Smith } 3773821be0aSBarry Smith 3783821be0aSBarry Smith /*@ 3793821be0aSBarry Smith PCGetKSPNestLevel - gets the amount of nesting the `KSP` that contains this `PC` has 3803821be0aSBarry Smith 3817a99bfcaSBarry Smith Not Collective 3823821be0aSBarry Smith 3833821be0aSBarry Smith Input Parameter: 3843821be0aSBarry Smith . pc - the `PC` 3853821be0aSBarry Smith 3863821be0aSBarry Smith Output Parameter: 3873821be0aSBarry Smith . level - the nest level 3883821be0aSBarry Smith 3893821be0aSBarry Smith Level: developer 3903821be0aSBarry Smith 3913821be0aSBarry Smith .seealso: [](ch_ksp), `KSPSetUp()`, `KSPSolve()`, `KSPDestroy()`, `KSP`, `KSPGMRES`, `KSPType`, `KSPSetNestLevel()`, `PCSetKSPNestLevel()`, `KSPGetNestLevel()` 3923821be0aSBarry Smith @*/ 3933821be0aSBarry Smith PetscErrorCode PCGetKSPNestLevel(PC pc, PetscInt *level) 3943821be0aSBarry Smith { 3953821be0aSBarry Smith PetscFunctionBegin; 3967a99bfcaSBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 3977a99bfcaSBarry Smith PetscAssertPointer(level, 2); 3983821be0aSBarry Smith *level = pc->kspnestlevel; 3993821be0aSBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 4003821be0aSBarry Smith } 4013821be0aSBarry Smith 4023821be0aSBarry Smith /*@ 403f1580f4eSBarry Smith PCCreate - Creates a preconditioner context, `PC` 4044b9ad928SBarry Smith 405d083f849SBarry Smith Collective 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith Input Parameter: 4084b9ad928SBarry Smith . comm - MPI communicator 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith Output Parameter: 4117a99bfcaSBarry Smith . newpc - location to put the preconditioner context 4124b9ad928SBarry Smith 41320f4b53cSBarry Smith Level: developer 41420f4b53cSBarry Smith 415f1580f4eSBarry Smith Note: 416f1580f4eSBarry Smith The default preconditioner for sparse matrices is `PCILU` or `PCICC` with 0 fill on one process and block Jacobi (`PCBJACOBI`) with `PCILU` or `PCICC` 417f1580f4eSBarry Smith in parallel. For dense matrices it is always `PCNONE`. 4184b9ad928SBarry Smith 419f1580f4eSBarry Smith .seealso: `PC`, `PCSetUp()`, `PCApply()`, `PCDestroy()` 4204b9ad928SBarry Smith @*/ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCreate(MPI_Comm comm, PC *newpc) 422d71ae5a4SJacob Faibussowitsch { 4234b9ad928SBarry Smith PC pc; 4244b9ad928SBarry Smith 4254b9ad928SBarry Smith PetscFunctionBegin; 4264f572ea9SToby Isaac PetscAssertPointer(newpc, 2); 4270a545947SLisandro Dalcin *newpc = NULL; 4289566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 4294b9ad928SBarry Smith 4309566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(pc, PC_CLASSID, "PC", "Preconditioner", "PC", comm, PCDestroy, PCView)); 431a7d21a52SLisandro Dalcin 4320a545947SLisandro Dalcin pc->mat = NULL; 4330a545947SLisandro Dalcin pc->pmat = NULL; 4344b9ad928SBarry Smith pc->setupcalled = 0; 4350c24e6a1SHong Zhang pc->setfromoptionscalled = 0; 4360a545947SLisandro Dalcin pc->data = NULL; 4374b9ad928SBarry Smith pc->diagonalscale = PETSC_FALSE; 4380a545947SLisandro Dalcin pc->diagonalscaleleft = NULL; 4390a545947SLisandro Dalcin pc->diagonalscaleright = NULL; 4404b9ad928SBarry Smith 4410a545947SLisandro Dalcin pc->modifysubmatrices = NULL; 4420a545947SLisandro Dalcin pc->modifysubmatricesP = NULL; 4432fa5cd67SKarl Rupp 444a7d21a52SLisandro Dalcin *newpc = pc; 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4464b9ad928SBarry Smith } 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith /*@ 4494b9ad928SBarry Smith PCApply - Applies the preconditioner to a vector. 4504b9ad928SBarry Smith 451c3339decSBarry Smith Collective 4524b9ad928SBarry Smith 4534b9ad928SBarry Smith Input Parameters: 4544b9ad928SBarry Smith + pc - the preconditioner context 4554b9ad928SBarry Smith - x - input vector 4564b9ad928SBarry Smith 4574b9ad928SBarry Smith Output Parameter: 4584b9ad928SBarry Smith . y - output vector 4594b9ad928SBarry Smith 4604b9ad928SBarry Smith Level: developer 4614b9ad928SBarry Smith 462f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()`, `PCApplyBAorAB()` 4634b9ad928SBarry Smith @*/ 464d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply(PC pc, Vec x, Vec y) 465d71ae5a4SJacob Faibussowitsch { 466106e20bfSBarry Smith PetscInt m, n, mv, nv; 4674b9ad928SBarry Smith 4684b9ad928SBarry Smith PetscFunctionBegin; 4690700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 4700700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 4710700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 4727827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 473e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 474540bdfbdSVaclav Hapla /* use pmat to check vector sizes since for KSPLSQR the pmat may be of a different size than mat */ 4759566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m, &n)); 4769566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(x, &mv)); 4779566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(y, &nv)); 478540bdfbdSVaclav Hapla /* check pmat * y = x is feasible */ 47963a3b9bcSJacob Faibussowitsch PetscCheck(mv == m, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local rows %" PetscInt_FMT " does not equal input vector size %" PetscInt_FMT, m, mv); 48063a3b9bcSJacob Faibussowitsch PetscCheck(nv == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Preconditioner number of local columns %" PetscInt_FMT " does not equal output vector size %" PetscInt_FMT, n, nv); 4819566063dSJacob Faibussowitsch PetscCall(VecSetErrorIfLocked(y, 3)); 482106e20bfSBarry Smith 4839566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 4849566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 4859566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 486dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, apply, x, y); 4879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 488e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 4899566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4914b9ad928SBarry Smith } 4924b9ad928SBarry Smith 4934b9ad928SBarry Smith /*@ 494f1580f4eSBarry Smith PCMatApply - Applies the preconditioner to multiple vectors stored as a `MATDENSE`. Like `PCApply()`, Y and X must be different matrices. 495c41ea70eSPierre Jolivet 496c3339decSBarry Smith Collective 497c677e75fSPierre Jolivet 498c677e75fSPierre Jolivet Input Parameters: 499c677e75fSPierre Jolivet + pc - the preconditioner context 500c677e75fSPierre Jolivet - X - block of input vectors 501c677e75fSPierre Jolivet 502c677e75fSPierre Jolivet Output Parameter: 503c677e75fSPierre Jolivet . Y - block of output vectors 504c677e75fSPierre Jolivet 505c677e75fSPierre Jolivet Level: developer 506c677e75fSPierre Jolivet 507f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `KSPMatSolve()` 508c677e75fSPierre Jolivet @*/ 509d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMatApply(PC pc, Mat X, Mat Y) 510d71ae5a4SJacob Faibussowitsch { 511c677e75fSPierre Jolivet Mat A; 512c677e75fSPierre Jolivet Vec cy, cx; 513bd82155bSPierre Jolivet PetscInt m1, M1, m2, M2, n1, N1, n2, N2, m3, M3, n3, N3; 514c677e75fSPierre Jolivet PetscBool match; 515c677e75fSPierre Jolivet 516c677e75fSPierre Jolivet PetscFunctionBegin; 517c677e75fSPierre Jolivet PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 518e2b200f6SPierre Jolivet PetscValidHeaderSpecific(X, MAT_CLASSID, 2); 519e2b200f6SPierre Jolivet PetscValidHeaderSpecific(Y, MAT_CLASSID, 3); 520e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, X, 2); 521e2b200f6SPierre Jolivet PetscCheckSameComm(pc, 1, Y, 3); 5227827d75bSBarry Smith PetscCheck(Y != X, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "Y and X must be different matrices"); 5239566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, NULL, &A)); 5249566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m3, &n3)); 5259566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m2, &n2)); 5269566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Y, &m1, &n1)); 5279566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M3, &N3)); 5289566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, &M2, &N2)); 5299566063dSJacob Faibussowitsch PetscCall(MatGetSize(Y, &M1, &N1)); 53063a3b9bcSJacob Faibussowitsch PetscCheck(n1 == n2 && N1 == N2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible number of columns between block of input vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and block of output vectors (n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")", n2, N2, n1, N1); 53163a3b9bcSJacob Faibussowitsch PetscCheck(m2 == m3 && M2 == M3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of input vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m2, M2, m3, M3, n3, N3); 53263a3b9bcSJacob Faibussowitsch PetscCheck(m1 == n3 && M1 == N3, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible layout between block of output vectors (m,M) = (%" PetscInt_FMT ",%" PetscInt_FMT ") and Pmat (m,M)x(n,N) = (%" PetscInt_FMT ",%" PetscInt_FMT ")x(%" PetscInt_FMT ",%" PetscInt_FMT ")", m1, M1, m3, M3, n3, N3); 5339566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)Y, &match, MATSEQDENSE, MATMPIDENSE, "")); 53428b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of output vectors not stored in a dense Mat"); 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, "")); 53628b400f6SJacob Faibussowitsch PetscCheck(match, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of input vectors not stored in a dense Mat"); 5379566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 538c677e75fSPierre Jolivet if (pc->ops->matapply) { 5399566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_MatApply, pc, X, Y, 0)); 540dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, matapply, X, Y); 5419566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_MatApply, pc, X, Y, 0)); 542c677e75fSPierre Jolivet } else { 5439566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "PC type %s applying column by column\n", ((PetscObject)pc)->type_name)); 544bd82155bSPierre Jolivet for (n1 = 0; n1 < N1; ++n1) { 5459566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, n1, &cx)); 5469566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, n1, &cy)); 5479566063dSJacob Faibussowitsch PetscCall(PCApply(pc, cx, cy)); 5489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, n1, &cy)); 5499566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, n1, &cx)); 550c677e75fSPierre Jolivet } 551c677e75fSPierre Jolivet } 5523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 553c677e75fSPierre Jolivet } 554c677e75fSPierre Jolivet 555c677e75fSPierre Jolivet /*@ 5564b9ad928SBarry Smith PCApplySymmetricLeft - Applies the left part of a symmetric preconditioner to a vector. 5574b9ad928SBarry Smith 558c3339decSBarry Smith Collective 5594b9ad928SBarry Smith 5604b9ad928SBarry Smith Input Parameters: 5614b9ad928SBarry Smith + pc - the preconditioner context 5624b9ad928SBarry Smith - x - input vector 5634b9ad928SBarry Smith 5644b9ad928SBarry Smith Output Parameter: 5654b9ad928SBarry Smith . y - output vector 5664b9ad928SBarry Smith 56720f4b53cSBarry Smith Level: developer 56820f4b53cSBarry Smith 569f1580f4eSBarry Smith Note: 570f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 5714b9ad928SBarry Smith 572f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricRight()` 5734b9ad928SBarry Smith @*/ 574d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricLeft(PC pc, Vec x, Vec y) 575d71ae5a4SJacob Faibussowitsch { 5764b9ad928SBarry Smith PetscFunctionBegin; 5770700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 5780700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 5790700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 5807827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 581e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 5829566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 5839566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 5849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricLeft, pc, x, y, 0)); 585dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricleft, x, y); 5869566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricLeft, pc, x, y, 0)); 5879566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 588e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 5893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5904b9ad928SBarry Smith } 5914b9ad928SBarry Smith 5924b9ad928SBarry Smith /*@ 5934b9ad928SBarry Smith PCApplySymmetricRight - Applies the right part of a symmetric preconditioner to a vector. 5944b9ad928SBarry Smith 595c3339decSBarry Smith Collective 5964b9ad928SBarry Smith 5974b9ad928SBarry Smith Input Parameters: 5984b9ad928SBarry Smith + pc - the preconditioner context 5994b9ad928SBarry Smith - x - input vector 6004b9ad928SBarry Smith 6014b9ad928SBarry Smith Output Parameter: 6024b9ad928SBarry Smith . y - output vector 6034b9ad928SBarry Smith 6044b9ad928SBarry Smith Level: developer 6054b9ad928SBarry Smith 606f1580f4eSBarry Smith Note: 607f1580f4eSBarry Smith Currently, this routine is implemented only for `PCICC` and `PCJACOBI` preconditioners. 6084b9ad928SBarry Smith 609f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplySymmetricLeft()` 6104b9ad928SBarry Smith @*/ 611d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplySymmetricRight(PC pc, Vec x, Vec y) 612d71ae5a4SJacob Faibussowitsch { 6134b9ad928SBarry Smith PetscFunctionBegin; 6140700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6150700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6160700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6177827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 618e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6199566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6209566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6219566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplySymmetricRight, pc, x, y, 0)); 622dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applysymmetricright, x, y); 6239566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplySymmetricRight, pc, x, y, 0)); 6249566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 625e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6274b9ad928SBarry Smith } 6284b9ad928SBarry Smith 6294b9ad928SBarry Smith /*@ 6304b9ad928SBarry Smith PCApplyTranspose - Applies the transpose of preconditioner to a vector. 6314b9ad928SBarry Smith 632c3339decSBarry Smith Collective 6334b9ad928SBarry Smith 6344b9ad928SBarry Smith Input Parameters: 6354b9ad928SBarry Smith + pc - the preconditioner context 6364b9ad928SBarry Smith - x - input vector 6374b9ad928SBarry Smith 6384b9ad928SBarry Smith Output Parameter: 6394b9ad928SBarry Smith . y - output vector 6404b9ad928SBarry Smith 64120f4b53cSBarry Smith Level: developer 64220f4b53cSBarry Smith 643f1580f4eSBarry Smith Note: 64495452b02SPatrick Sanan For complex numbers this applies the non-Hermitian transpose. 6454c97465dSBarry Smith 646feefa0e1SJacob Faibussowitsch Developer Notes: 647f1580f4eSBarry Smith We need to implement a `PCApplyHermitianTranspose()` 6484c97465dSBarry Smith 649f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyBAorAB()`, `PCApplyBAorABTranspose()`, `PCApplyTransposeExists()` 6504b9ad928SBarry Smith @*/ 651d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTranspose(PC pc, Vec x, Vec y) 652d71ae5a4SJacob Faibussowitsch { 6534b9ad928SBarry Smith PetscFunctionBegin; 6540700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6550700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 2); 6560700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 6577827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 658e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 2, PETSC_TRUE)); 6599566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 6609566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(x)); 6619566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_Apply, pc, x, y, 0)); 662dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applytranspose, x, y); 6639566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_Apply, pc, x, y, 0)); 6649566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(x)); 665e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 3, PETSC_FALSE)); 6663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6674b9ad928SBarry Smith } 6684b9ad928SBarry Smith 6694b9ad928SBarry Smith /*@ 6709cc050e5SLisandro Dalcin PCApplyTransposeExists - Test whether the preconditioner has a transpose apply operation 6714b9ad928SBarry Smith 672c3339decSBarry Smith Collective 6734b9ad928SBarry Smith 674f1580f4eSBarry Smith Input Parameter: 6754b9ad928SBarry Smith . pc - the preconditioner context 6764b9ad928SBarry Smith 6774b9ad928SBarry Smith Output Parameter: 678f1580f4eSBarry Smith . flg - `PETSC_TRUE` if a transpose operation is defined 6794b9ad928SBarry Smith 6804b9ad928SBarry Smith Level: developer 6814b9ad928SBarry Smith 682f1580f4eSBarry Smith .seealso: `PC`, `PCApplyTranspose()` 6834b9ad928SBarry Smith @*/ 684d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyTransposeExists(PC pc, PetscBool *flg) 685d71ae5a4SJacob Faibussowitsch { 6864b9ad928SBarry Smith PetscFunctionBegin; 6870700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 6884f572ea9SToby Isaac PetscAssertPointer(flg, 2); 6896ce5e81cSLisandro Dalcin if (pc->ops->applytranspose) *flg = PETSC_TRUE; 6906ce5e81cSLisandro Dalcin else *flg = PETSC_FALSE; 6913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6924b9ad928SBarry Smith } 6934b9ad928SBarry Smith 6944b9ad928SBarry Smith /*@ 6951a2f9f99SBarry Smith PCApplyBAorAB - Applies the preconditioner and operator to a vector. y = B*A*x or y = A*B*x. 6964b9ad928SBarry Smith 697c3339decSBarry Smith Collective 6984b9ad928SBarry Smith 6994b9ad928SBarry Smith Input Parameters: 7004b9ad928SBarry Smith + pc - the preconditioner context 701f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 7024b9ad928SBarry Smith . x - input vector 7034b9ad928SBarry Smith - work - work vector 7044b9ad928SBarry Smith 7054b9ad928SBarry Smith Output Parameter: 7064b9ad928SBarry Smith . y - output vector 7074b9ad928SBarry Smith 7084b9ad928SBarry Smith Level: developer 7094b9ad928SBarry Smith 710f1580f4eSBarry Smith Note: 711f1580f4eSBarry Smith If the `PC` has had `PCSetDiagonalScale()` set then D M A D^{-1} for left preconditioning or D A M D^{-1} is actually applied. Note that the 712f1580f4eSBarry Smith specific `KSPSolve()` method must also be written to handle the post-solve "correction" for the diagonal scaling. 71313e0d083SBarry Smith 714f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorABTranspose()` 7154b9ad928SBarry Smith @*/ 716d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorAB(PC pc, PCSide side, Vec x, Vec y, Vec work) 717d71ae5a4SJacob Faibussowitsch { 7184b9ad928SBarry Smith PetscFunctionBegin; 7190700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 720a69c7061SStefano Zampini PetscValidLogicalCollectiveEnum(pc, side, 2); 7210700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 7220700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 7230700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 724a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, x, 3); 725a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, y, 4); 726a69c7061SStefano Zampini PetscCheckSameComm(pc, 1, work, 5); 7277827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 7287827d75bSBarry Smith PetscCheck(side == PC_LEFT || side == PC_SYMMETRIC || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right, left, or symmetric"); 7297827d75bSBarry Smith PetscCheck(!pc->diagonalscale || side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot include diagonal scaling with symmetric preconditioner application"); 730e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 7314b9ad928SBarry Smith 7329566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 7334b9ad928SBarry Smith if (pc->diagonalscale) { 7344b9ad928SBarry Smith if (pc->ops->applyBA) { 7354b9ad928SBarry Smith Vec work2; /* this is expensive, but to fix requires a second work vector argument to PCApplyBAorAB() */ 7369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x, &work2)); 7379566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, work2)); 738dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, work2, y, work); 7399566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&work2)); 7414b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7429566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7439566063dSJacob Faibussowitsch PetscCall(PCApply(pc, y, work)); 7449566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7459566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7464b9ad928SBarry Smith } else if (side == PC_LEFT) { 7479566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleRight(pc, x, y)); 7489566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, y, work)); 7499566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7509566063dSJacob Faibussowitsch PetscCall(PCDiagonalScaleLeft(pc, y, y)); 7517827d75bSBarry Smith } else PetscCheck(side != PC_SYMMETRIC, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot provide diagonal scaling with symmetric application of preconditioner"); 7524b9ad928SBarry Smith } else { 7534b9ad928SBarry Smith if (pc->ops->applyBA) { 754dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBA, side, x, y, work); 7554b9ad928SBarry Smith } else if (side == PC_RIGHT) { 7569566063dSJacob Faibussowitsch PetscCall(PCApply(pc, x, work)); 7579566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7584b9ad928SBarry Smith } else if (side == PC_LEFT) { 7599566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, x, work)); 7609566063dSJacob Faibussowitsch PetscCall(PCApply(pc, work, y)); 7614b9ad928SBarry Smith } else if (side == PC_SYMMETRIC) { 7624b9ad928SBarry Smith /* There's an extra copy here; maybe should provide 2 work vectors instead? */ 7639566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricRight(pc, x, work)); 7649566063dSJacob Faibussowitsch PetscCall(MatMult(pc->mat, work, y)); 7659566063dSJacob Faibussowitsch PetscCall(VecCopy(y, work)); 7669566063dSJacob Faibussowitsch PetscCall(PCApplySymmetricLeft(pc, work, y)); 7674b9ad928SBarry Smith } 7684b9ad928SBarry Smith } 769e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 7703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7714b9ad928SBarry Smith } 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith /*@ 7744b9ad928SBarry Smith PCApplyBAorABTranspose - Applies the transpose of the preconditioner 7754b9ad928SBarry Smith and operator to a vector. That is, applies tr(B) * tr(A) with left preconditioning, 7769b3038f0SBarry Smith NOT tr(B*A) = tr(A)*tr(B). 7774b9ad928SBarry Smith 778c3339decSBarry Smith Collective 7794b9ad928SBarry Smith 7804b9ad928SBarry Smith Input Parameters: 7814b9ad928SBarry Smith + pc - the preconditioner context 782f1580f4eSBarry Smith . side - indicates the preconditioner side, one of `PC_LEFT`, `PC_RIGHT`, or `PC_SYMMETRIC` 7834b9ad928SBarry Smith . x - input vector 7844b9ad928SBarry Smith - work - work vector 7854b9ad928SBarry Smith 7864b9ad928SBarry Smith Output Parameter: 7874b9ad928SBarry Smith . y - output vector 7884b9ad928SBarry Smith 78920f4b53cSBarry Smith Level: developer 79020f4b53cSBarry Smith 791f1580f4eSBarry Smith Note: 792f1580f4eSBarry Smith This routine is used internally so that the same Krylov code can be used to solve A x = b and A' x = b, with a preconditioner 7939b3038f0SBarry Smith defined by B'. This is why this has the funny form that it computes tr(B) * tr(A) 7949b3038f0SBarry Smith 795f1580f4eSBarry Smith .seealso: `PC`, `PCApply()`, `PCApplyTranspose()`, `PCApplyBAorAB()` 7964b9ad928SBarry Smith @*/ 797d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyBAorABTranspose(PC pc, PCSide side, Vec x, Vec y, Vec work) 798d71ae5a4SJacob Faibussowitsch { 7994b9ad928SBarry Smith PetscFunctionBegin; 8000700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8010700a824SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 8020700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 4); 8030700a824SBarry Smith PetscValidHeaderSpecific(work, VEC_CLASSID, 5); 8047827d75bSBarry Smith PetscCheck(x != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "x and y must be different vectors"); 805e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(x, 3, PETSC_TRUE)); 8064b9ad928SBarry Smith if (pc->ops->applyBAtranspose) { 807dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyBAtranspose, side, x, y, work); 808e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8104b9ad928SBarry Smith } 8117827d75bSBarry Smith PetscCheck(side == PC_LEFT || side == PC_RIGHT, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Side must be right or left"); 8124b9ad928SBarry Smith 8139566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 8144b9ad928SBarry Smith if (side == PC_RIGHT) { 8159566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, x, work)); 8169566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, work, y)); 8174b9ad928SBarry Smith } else if (side == PC_LEFT) { 8189566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->mat, x, work)); 8199566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(pc, work, y)); 8204b9ad928SBarry Smith } 8214b9ad928SBarry Smith /* add support for PC_SYMMETRIC */ 822e0f629ddSJacob Faibussowitsch if (pc->erroriffailure) PetscCall(VecValidValues_Internal(y, 4, PETSC_FALSE)); 8233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8244b9ad928SBarry Smith } 8254b9ad928SBarry Smith 8264b9ad928SBarry Smith /*@ 8274b9ad928SBarry Smith PCApplyRichardsonExists - Determines whether a particular preconditioner has a 8284b9ad928SBarry Smith built-in fast application of Richardson's method. 8294b9ad928SBarry Smith 8304b9ad928SBarry Smith Not Collective 8314b9ad928SBarry Smith 8324b9ad928SBarry Smith Input Parameter: 8334b9ad928SBarry Smith . pc - the preconditioner 8344b9ad928SBarry Smith 8354b9ad928SBarry Smith Output Parameter: 836f1580f4eSBarry Smith . exists - `PETSC_TRUE` or `PETSC_FALSE` 8374b9ad928SBarry Smith 8384b9ad928SBarry Smith Level: developer 8394b9ad928SBarry Smith 840f1580f4eSBarry Smith .seealso: `PC`, `PCRICHARDSON`, `PCApplyRichardson()` 8414b9ad928SBarry Smith @*/ 842d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardsonExists(PC pc, PetscBool *exists) 843d71ae5a4SJacob Faibussowitsch { 8444b9ad928SBarry Smith PetscFunctionBegin; 8450700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8464f572ea9SToby Isaac PetscAssertPointer(exists, 2); 8474b9ad928SBarry Smith if (pc->ops->applyrichardson) *exists = PETSC_TRUE; 8484b9ad928SBarry Smith else *exists = PETSC_FALSE; 8493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8504b9ad928SBarry Smith } 8514b9ad928SBarry Smith 8524b9ad928SBarry Smith /*@ 8534b9ad928SBarry Smith PCApplyRichardson - Applies several steps of Richardson iteration with 8544b9ad928SBarry Smith the particular preconditioner. This routine is usually used by the 8554b9ad928SBarry Smith Krylov solvers and not the application code directly. 8564b9ad928SBarry Smith 857c3339decSBarry Smith Collective 8584b9ad928SBarry Smith 8594b9ad928SBarry Smith Input Parameters: 8604b9ad928SBarry Smith + pc - the preconditioner context 8617319c654SBarry Smith . b - the right hand side 8624b9ad928SBarry Smith . w - one work vector 8634b9ad928SBarry Smith . rtol - relative decrease in residual norm convergence criteria 86470441072SBarry Smith . abstol - absolute residual norm convergence criteria 8654b9ad928SBarry Smith . dtol - divergence residual norm increase criteria 8667319c654SBarry Smith . its - the number of iterations to apply. 8677319c654SBarry Smith - guesszero - if the input x contains nonzero initial guess 8684b9ad928SBarry Smith 869d8d19677SJose E. Roman Output Parameters: 8704d0a8057SBarry Smith + outits - number of iterations actually used (for SOR this always equals its) 8714d0a8057SBarry Smith . reason - the reason the apply terminated 872f1580f4eSBarry Smith - y - the solution (also contains initial guess if guesszero is `PETSC_FALSE` 8734b9ad928SBarry Smith 87420f4b53cSBarry Smith Level: developer 87520f4b53cSBarry Smith 8764b9ad928SBarry Smith Notes: 8774b9ad928SBarry Smith Most preconditioners do not support this function. Use the command 878f1580f4eSBarry Smith `PCApplyRichardsonExists()` to determine if one does. 8794b9ad928SBarry Smith 880f1580f4eSBarry Smith Except for the `PCMG` this routine ignores the convergence tolerances 8814b9ad928SBarry Smith and always runs for the number of iterations 8824b9ad928SBarry Smith 883f1580f4eSBarry Smith .seealso: `PC`, `PCApplyRichardsonExists()` 8844b9ad928SBarry Smith @*/ 885d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApplyRichardson(PC pc, Vec b, Vec y, Vec w, PetscReal rtol, PetscReal abstol, PetscReal dtol, PetscInt its, PetscBool guesszero, PetscInt *outits, PCRichardsonConvergedReason *reason) 886d71ae5a4SJacob Faibussowitsch { 8874b9ad928SBarry Smith PetscFunctionBegin; 8880700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8890700a824SBarry Smith PetscValidHeaderSpecific(b, VEC_CLASSID, 2); 8900700a824SBarry Smith PetscValidHeaderSpecific(y, VEC_CLASSID, 3); 8910700a824SBarry Smith PetscValidHeaderSpecific(w, VEC_CLASSID, 4); 8927827d75bSBarry Smith PetscCheck(b != y, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_IDN, "b and y must be different vectors"); 8939566063dSJacob Faibussowitsch PetscCall(PCSetUp(pc)); 894dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, applyrichardson, b, y, w, rtol, abstol, dtol, its, guesszero, outits, reason); 8953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8964b9ad928SBarry Smith } 8974b9ad928SBarry Smith 898422a814eSBarry Smith /*@ 899f1580f4eSBarry Smith PCSetFailedReason - Sets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 9001b2b9847SBarry Smith 901c3339decSBarry Smith Logically Collective 9021b2b9847SBarry Smith 903d8d19677SJose E. Roman Input Parameters: 9041b2b9847SBarry Smith + pc - the preconditioner context 9051b2b9847SBarry Smith - reason - the reason it failedx 9061b2b9847SBarry Smith 9071b2b9847SBarry Smith Level: advanced 9081b2b9847SBarry Smith 909f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCFailedReason` 9101b2b9847SBarry Smith @*/ 911d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFailedReason(PC pc, PCFailedReason reason) 912d71ae5a4SJacob Faibussowitsch { 9131b2b9847SBarry Smith PetscFunctionBegin; 9146479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9151b2b9847SBarry Smith pc->failedreason = reason; 9163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9171b2b9847SBarry Smith } 9181b2b9847SBarry Smith 9191b2b9847SBarry Smith /*@ 920f1580f4eSBarry Smith PCGetFailedReason - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail 921422a814eSBarry Smith 922c3339decSBarry Smith Logically Collective 923422a814eSBarry Smith 924422a814eSBarry Smith Input Parameter: 925422a814eSBarry Smith . pc - the preconditioner context 926422a814eSBarry Smith 927422a814eSBarry Smith Output Parameter: 9281b2b9847SBarry Smith . reason - the reason it failed 929422a814eSBarry Smith 930422a814eSBarry Smith Level: advanced 931422a814eSBarry Smith 932f1580f4eSBarry Smith Note: 933f1580f4eSBarry Smith This is the maximum over reason over all ranks in the PC communicator. It is only valid after 9346479e835SStefano Zampini a call `KSPCheckDot()` or `KSPCheckNorm()` inside a `KSPSolve()` or `PCReduceFailedReason()`. 9356479e835SStefano Zampini It is not valid immediately after a `PCSetUp()` or `PCApply()`, then use `PCGetFailedReasonRank()` 9361b2b9847SBarry Smith 937feefa0e1SJacob Faibussowitsch .seealso: `PC`, ``PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReasonRank()`, `PCSetFailedReason()` 938422a814eSBarry Smith @*/ 939d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReason(PC pc, PCFailedReason *reason) 940d71ae5a4SJacob Faibussowitsch { 941422a814eSBarry Smith PetscFunctionBegin; 9426479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9439754764cSHong Zhang if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9440ae0b32dSHong Zhang else *reason = pc->failedreason; 9453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 946422a814eSBarry Smith } 947422a814eSBarry Smith 9481b2b9847SBarry Smith /*@ 949f1580f4eSBarry Smith PCGetFailedReasonRank - Gets the reason a `PCSetUp()` failed or `PC_NOERROR` if it did not fail on this MPI rank 9501b2b9847SBarry Smith 951f1580f4eSBarry Smith Not Collective 9521b2b9847SBarry Smith 9531b2b9847SBarry Smith Input Parameter: 9541b2b9847SBarry Smith . pc - the preconditioner context 9551b2b9847SBarry Smith 9561b2b9847SBarry Smith Output Parameter: 9571b2b9847SBarry Smith . reason - the reason it failed 9581b2b9847SBarry Smith 95920f4b53cSBarry Smith Level: advanced 96020f4b53cSBarry Smith 961f1580f4eSBarry Smith Note: 962f1580f4eSBarry Smith Different ranks may have different reasons or no reason, see `PCGetFailedReason()` 9631b2b9847SBarry Smith 9646479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()`, `PCReduceFailedReason()` 9651b2b9847SBarry Smith @*/ 966d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetFailedReasonRank(PC pc, PCFailedReason *reason) 967d71ae5a4SJacob Faibussowitsch { 9681b2b9847SBarry Smith PetscFunctionBegin; 9696479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9701b2b9847SBarry Smith if (pc->setupcalled < 0) *reason = (PCFailedReason)pc->setupcalled; 9711b2b9847SBarry Smith else *reason = pc->failedreason; 9723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9731b2b9847SBarry Smith } 974422a814eSBarry Smith 9756479e835SStefano Zampini /*@ 9766479e835SStefano Zampini PCReduceFailedReason - Reduce the failed reason among the MPI processes that share the `PC` 9776479e835SStefano Zampini 9786479e835SStefano Zampini Collective 9796479e835SStefano Zampini 9806479e835SStefano Zampini Input Parameter: 9816479e835SStefano Zampini . pc - the preconditioner context 9826479e835SStefano Zampini 9836479e835SStefano Zampini Level: advanced 9846479e835SStefano Zampini 9856479e835SStefano Zampini Note: 9866479e835SStefano Zampini Different MPI ranks may have different reasons or no reason, see `PCGetFailedReason()`. This routine 9876479e835SStefano Zampini makes them have a common value (failure if any MPI process had a failure). 9886479e835SStefano Zampini 9896479e835SStefano Zampini .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()`, `PCGetFailedReason()`, `PCSetFailedReason()` 9906479e835SStefano Zampini @*/ 9916479e835SStefano Zampini PetscErrorCode PCReduceFailedReason(PC pc) 9926479e835SStefano Zampini { 9936479e835SStefano Zampini PetscInt buf; 9946479e835SStefano Zampini 9956479e835SStefano Zampini PetscFunctionBegin; 9966479e835SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9976479e835SStefano Zampini buf = (PetscInt)pc->failedreason; 9986479e835SStefano Zampini PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &buf, 1, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)pc))); 9996479e835SStefano Zampini pc->failedreason = (PCFailedReason)buf; 10006479e835SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 10016479e835SStefano Zampini } 10026479e835SStefano Zampini 1003c00cb57fSBarry Smith /* Next line needed to deactivate KSP_Solve logging */ 1004c00cb57fSBarry Smith #include <petsc/private/kspimpl.h> 1005c00cb57fSBarry Smith 10064b9ad928SBarry Smith /* 10074b9ad928SBarry Smith a setupcall of 0 indicates never setup, 100823ee1639SBarry Smith 1 indicates has been previously setup 1009422a814eSBarry Smith -1 indicates a PCSetUp() was attempted and failed 10104b9ad928SBarry Smith */ 10114b9ad928SBarry Smith /*@ 10124b9ad928SBarry Smith PCSetUp - Prepares for the use of a preconditioner. 10134b9ad928SBarry Smith 1014c3339decSBarry Smith Collective 10154b9ad928SBarry Smith 10164b9ad928SBarry Smith Input Parameter: 10174b9ad928SBarry Smith . pc - the preconditioner context 10184b9ad928SBarry Smith 10194b9ad928SBarry Smith Level: developer 10204b9ad928SBarry Smith 1021f1580f4eSBarry Smith .seealso: `PC`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10224b9ad928SBarry Smith @*/ 1023d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUp(PC pc) 1024d71ae5a4SJacob Faibussowitsch { 1025566e8bf2SBarry Smith const char *def; 1026fc9b51d3SBarry Smith PetscObjectState matstate, matnonzerostate; 10274b9ad928SBarry Smith 10284b9ad928SBarry Smith PetscFunctionBegin; 10290700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 103028b400f6SJacob Faibussowitsch PetscCheck(pc->mat, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Matrix must be set first"); 10316ce5e81cSLisandro Dalcin 103223ee1639SBarry Smith if (pc->setupcalled && pc->reusepreconditioner) { 10339566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Leaving PC with identical preconditioner since reuse preconditioner is set\n")); 10343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10354b9ad928SBarry Smith } 10364b9ad928SBarry Smith 10379566063dSJacob Faibussowitsch PetscCall(PetscObjectStateGet((PetscObject)pc->pmat, &matstate)); 10389566063dSJacob Faibussowitsch PetscCall(MatGetNonzeroState(pc->pmat, &matnonzerostate)); 103923ee1639SBarry Smith if (!pc->setupcalled) { 1040*5e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC for first time\n")); 104123ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 10429df67409SStefano Zampini } else if (matstate == pc->matstate) PetscFunctionReturn(PETSC_SUCCESS); 10439df67409SStefano Zampini else { 10449df67409SStefano Zampini if (matnonzerostate != pc->matnonzerostate) { 10459566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc, "Setting up PC with different nonzero pattern\n")); 104623ee1639SBarry Smith pc->flag = DIFFERENT_NONZERO_PATTERN; 104723ee1639SBarry Smith } else { 1048*5e62d202SMark Adams //PetscCall(PetscInfo(pc, "Setting up PC with same nonzero pattern\n")); 104923ee1639SBarry Smith pc->flag = SAME_NONZERO_PATTERN; 105023ee1639SBarry Smith } 105123ee1639SBarry Smith } 105223ee1639SBarry Smith pc->matstate = matstate; 1053fc9b51d3SBarry Smith pc->matnonzerostate = matnonzerostate; 105423ee1639SBarry Smith 10557adad957SLisandro Dalcin if (!((PetscObject)pc)->type_name) { 10569566063dSJacob Faibussowitsch PetscCall(PCGetDefaultType_Private(pc, &def)); 10579566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, def)); 1058566e8bf2SBarry Smith } 10594b9ad928SBarry Smith 10609566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->pmat, pc->erroriffailure)); 10619566063dSJacob Faibussowitsch PetscCall(MatSetErrorIfFailure(pc->mat, pc->erroriffailure)); 10629566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUp, pc, 0, 0, 0)); 10634b9ad928SBarry Smith if (pc->ops->setup) { 1064c00cb57fSBarry Smith /* do not log solves and applications of preconditioners while constructing preconditioners; perhaps they should be logged separately from the regular solves */ 10659566063dSJacob Faibussowitsch PetscCall(KSPInitializePackage()); 10669566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePush(KSP_Solve)); 10679566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePush(PC_Apply)); 1068dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setup); 10699566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePop(KSP_Solve)); 10709566063dSJacob Faibussowitsch PetscCall(PetscLogEventDeactivatePop(PC_Apply)); 10714b9ad928SBarry Smith } 10729566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUp, pc, 0, 0, 0)); 1073422a814eSBarry Smith if (!pc->setupcalled) pc->setupcalled = 1; 10743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10754b9ad928SBarry Smith } 10764b9ad928SBarry Smith 10774b9ad928SBarry Smith /*@ 10784b9ad928SBarry Smith PCSetUpOnBlocks - Sets up the preconditioner for each block in 10794b9ad928SBarry Smith the block Jacobi, block Gauss-Seidel, and overlapping Schwarz 10804b9ad928SBarry Smith methods. 10814b9ad928SBarry Smith 1082c3339decSBarry Smith Collective 10834b9ad928SBarry Smith 1084f1580f4eSBarry Smith Input Parameter: 10854b9ad928SBarry Smith . pc - the preconditioner context 10864b9ad928SBarry Smith 10874b9ad928SBarry Smith Level: developer 10884b9ad928SBarry Smith 1089f1580f4eSBarry Smith Note: 1090f1580f4eSBarry Smith For nested preconditioners such as `PCBJACOBI` `PCSetUp()` is not called on each sub-`KSP` when `PCSetUp()` is 1091f1580f4eSBarry Smith called on the outer `PC`, this routine ensures it is called. 1092f1580f4eSBarry Smith 1093feefa0e1SJacob Faibussowitsch .seealso: `PC`, `PCSetUp()`, `PCCreate()`, `PCApply()`, `PCDestroy()` 10944b9ad928SBarry Smith @*/ 1095d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetUpOnBlocks(PC pc) 1096d71ae5a4SJacob Faibussowitsch { 10974b9ad928SBarry Smith PetscFunctionBegin; 10980700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 10993ba16761SJacob Faibussowitsch if (!pc->ops->setuponblocks) PetscFunctionReturn(PETSC_SUCCESS); 11009566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_SetUpOnBlocks, pc, 0, 0, 0)); 1101dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, setuponblocks); 11029566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_SetUpOnBlocks, pc, 0, 0, 0)); 11033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11044b9ad928SBarry Smith } 11054b9ad928SBarry Smith 11064b9ad928SBarry Smith /*@C 11074b9ad928SBarry Smith PCSetModifySubMatrices - Sets a user-defined routine for modifying the 110804c3f3b8SBarry Smith submatrices that arise within certain subdomain-based preconditioners such as `PCASM` 11094b9ad928SBarry Smith 1110c3339decSBarry Smith Logically Collective 11114b9ad928SBarry Smith 11124b9ad928SBarry Smith Input Parameters: 11134b9ad928SBarry Smith + pc - the preconditioner context 11144b9ad928SBarry Smith . func - routine for modifying the submatrices 111504c3f3b8SBarry Smith - ctx - optional user-defined context (may be `NULL`) 11164b9ad928SBarry Smith 111720f4b53cSBarry Smith Calling sequence of `func`: 111820f4b53cSBarry Smith + pc - the preconditioner context 111904c3f3b8SBarry Smith . nsub - number of index sets 112020f4b53cSBarry Smith . row - an array of index sets that contain the global row numbers 11214b9ad928SBarry Smith that comprise each local submatrix 11224b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11234b9ad928SBarry Smith that comprise each local submatrix 11244b9ad928SBarry Smith . submat - array of local submatrices 11254b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 112604c3f3b8SBarry Smith user-defined func routine (may be `NULL`) 11274b9ad928SBarry Smith 112820f4b53cSBarry Smith Level: advanced 112920f4b53cSBarry Smith 11304b9ad928SBarry Smith Notes: 113104c3f3b8SBarry Smith The basic submatrices are extracted from the preconditioner matrix as 113204c3f3b8SBarry Smith usual; the user can then alter these (for example, to set different boundary 113304c3f3b8SBarry Smith conditions for each submatrix) before they are used for the local solves. 113404c3f3b8SBarry Smith 1135f1580f4eSBarry Smith `PCSetModifySubMatrices()` MUST be called before `KSPSetUp()` and 1136f1580f4eSBarry Smith `KSPSolve()`. 11374b9ad928SBarry Smith 1138f1580f4eSBarry Smith A routine set by `PCSetModifySubMatrices()` is currently called within 1139f1580f4eSBarry Smith the block Jacobi (`PCBJACOBI`) and additive Schwarz (`PCASM`) 11404b9ad928SBarry Smith preconditioners. All other preconditioners ignore this routine. 11414b9ad928SBarry Smith 1142f1580f4eSBarry Smith .seealso: `PC`, `PCBJACOBI`, `PCASM`, `PCModifySubMatrices()` 11434b9ad928SBarry Smith @*/ 114404c3f3b8SBarry Smith PetscErrorCode PCSetModifySubMatrices(PC pc, PetscErrorCode (*func)(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx), void *ctx) 1145d71ae5a4SJacob Faibussowitsch { 11464b9ad928SBarry Smith PetscFunctionBegin; 11470700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11484b9ad928SBarry Smith pc->modifysubmatrices = func; 11494b9ad928SBarry Smith pc->modifysubmatricesP = ctx; 11503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11514b9ad928SBarry Smith } 11524b9ad928SBarry Smith 11534b9ad928SBarry Smith /*@C 11544b9ad928SBarry Smith PCModifySubMatrices - Calls an optional user-defined routine within 1155f1580f4eSBarry Smith certain preconditioners if one has been set with `PCSetModifySubMatrices()`. 11564b9ad928SBarry Smith 1157c3339decSBarry Smith Collective 11584b9ad928SBarry Smith 11594b9ad928SBarry Smith Input Parameters: 11604b9ad928SBarry Smith + pc - the preconditioner context 11614b9ad928SBarry Smith . nsub - the number of local submatrices 11624b9ad928SBarry Smith . row - an array of index sets that contain the global row numbers 11634b9ad928SBarry Smith that comprise each local submatrix 11644b9ad928SBarry Smith . col - an array of index sets that contain the global column numbers 11654b9ad928SBarry Smith that comprise each local submatrix 11664b9ad928SBarry Smith . submat - array of local submatrices 11674b9ad928SBarry Smith - ctx - optional user-defined context for private data for the 11684b9ad928SBarry Smith user-defined routine (may be null) 11694b9ad928SBarry Smith 11704b9ad928SBarry Smith Output Parameter: 11714b9ad928SBarry Smith . submat - array of local submatrices (the entries of which may 11724b9ad928SBarry Smith have been modified) 11734b9ad928SBarry Smith 117420f4b53cSBarry Smith Level: developer 117520f4b53cSBarry Smith 117604c3f3b8SBarry Smith Note: 11774b9ad928SBarry Smith The user should NOT generally call this routine, as it will 117804c3f3b8SBarry Smith automatically be called within certain preconditioners. 11794b9ad928SBarry Smith 1180f1580f4eSBarry Smith .seealso: `PC`, `PCSetModifySubMatrices()` 11814b9ad928SBarry Smith @*/ 1182d71ae5a4SJacob Faibussowitsch PetscErrorCode PCModifySubMatrices(PC pc, PetscInt nsub, const IS row[], const IS col[], Mat submat[], void *ctx) 1183d71ae5a4SJacob Faibussowitsch { 11844b9ad928SBarry Smith PetscFunctionBegin; 11850700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11863ba16761SJacob Faibussowitsch if (!pc->modifysubmatrices) PetscFunctionReturn(PETSC_SUCCESS); 11879566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ModifySubMatrices, pc, 0, 0, 0)); 11889566063dSJacob Faibussowitsch PetscCall((*pc->modifysubmatrices)(pc, nsub, row, col, submat, ctx)); 11899566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ModifySubMatrices, pc, 0, 0, 0)); 11903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11914b9ad928SBarry Smith } 11924b9ad928SBarry Smith 11934b9ad928SBarry Smith /*@ 11944b9ad928SBarry Smith PCSetOperators - Sets the matrix associated with the linear system and 11954b9ad928SBarry Smith a (possibly) different one associated with the preconditioner. 11964b9ad928SBarry Smith 1197c3339decSBarry Smith Logically Collective 11984b9ad928SBarry Smith 11994b9ad928SBarry Smith Input Parameters: 12004b9ad928SBarry Smith + pc - the preconditioner context 1201e5d3d808SBarry Smith . Amat - the matrix that defines the linear system 1202d1e9a80fSBarry Smith - Pmat - the matrix to be used in constructing the preconditioner, usually the same as Amat. 12034b9ad928SBarry Smith 120420f4b53cSBarry Smith Level: intermediate 1205189c0b8aSBarry Smith 120620f4b53cSBarry Smith Notes: 120720f4b53cSBarry Smith Passing a `NULL` for `Amat` or `Pmat` removes the matrix that is currently used. 120820f4b53cSBarry Smith 120920f4b53cSBarry Smith If you wish to replace either `Amat` or `Pmat` but leave the other one untouched then 1210f1580f4eSBarry Smith first call `KSPGetOperators()` to get the one you wish to keep, call `PetscObjectReference()` 1211f1580f4eSBarry Smith on it and then pass it back in in your call to `KSPSetOperators()`. 1212189c0b8aSBarry Smith 12134b9ad928SBarry Smith More Notes about Repeated Solution of Linear Systems: 121420f4b53cSBarry Smith PETSc does NOT reset the matrix entries of either `Amat` or `Pmat` 12154b9ad928SBarry Smith to zero after a linear solve; the user is completely responsible for 1216f1580f4eSBarry Smith matrix assembly. See the routine `MatZeroEntries()` if desiring to 12174b9ad928SBarry Smith zero all elements of a matrix. 12184b9ad928SBarry Smith 1219f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()` 12204b9ad928SBarry Smith @*/ 1221d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOperators(PC pc, Mat Amat, Mat Pmat) 1222d71ae5a4SJacob Faibussowitsch { 12233b3f6333SBarry Smith PetscInt m1, n1, m2, n2; 12244b9ad928SBarry Smith 12254b9ad928SBarry Smith PetscFunctionBegin; 12260700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 12270700a824SBarry Smith if (Amat) PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 12280700a824SBarry Smith if (Pmat) PetscValidHeaderSpecific(Pmat, MAT_CLASSID, 3); 12293fc8bf9cSBarry Smith if (Amat) PetscCheckSameComm(pc, 1, Amat, 2); 12303fc8bf9cSBarry Smith if (Pmat) PetscCheckSameComm(pc, 1, Pmat, 3); 123131641f1aSBarry Smith if (pc->setupcalled && pc->mat && pc->pmat && Amat && Pmat) { 12329566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Amat, &m1, &n1)); 12339566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat, &m2, &n2)); 123463a3b9bcSJacob Faibussowitsch PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Amat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1); 12359566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Pmat, &m1, &n1)); 12369566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &m2, &n2)); 123763a3b9bcSJacob Faibussowitsch PetscCheck(m1 == m2 && n1 == n2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot change local size of Pmat after use old sizes %" PetscInt_FMT " %" PetscInt_FMT " new sizes %" PetscInt_FMT " %" PetscInt_FMT, m2, n2, m1, n1); 12383b3f6333SBarry Smith } 12394b9ad928SBarry Smith 124023ee1639SBarry Smith if (Pmat != pc->pmat) { 124123ee1639SBarry Smith /* changing the operator that defines the preconditioner thus reneed to clear current states so new preconditioner is built */ 124223ee1639SBarry Smith pc->matnonzerostate = -1; 124323ee1639SBarry Smith pc->matstate = -1; 124423ee1639SBarry Smith } 124523ee1639SBarry Smith 1246906ed7ccSBarry Smith /* reference first in case the matrices are the same */ 12479566063dSJacob Faibussowitsch if (Amat) PetscCall(PetscObjectReference((PetscObject)Amat)); 12489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->mat)); 12499566063dSJacob Faibussowitsch if (Pmat) PetscCall(PetscObjectReference((PetscObject)Pmat)); 12509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 12514b9ad928SBarry Smith pc->mat = Amat; 12524b9ad928SBarry Smith pc->pmat = Pmat; 12533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1254e56f5c9eSBarry Smith } 1255e56f5c9eSBarry Smith 125623ee1639SBarry Smith /*@ 125723ee1639SBarry Smith PCSetReusePreconditioner - reuse the current preconditioner even if the operator in the preconditioner has changed. 125823ee1639SBarry Smith 1259c3339decSBarry Smith Logically Collective 126023ee1639SBarry Smith 126123ee1639SBarry Smith Input Parameters: 126223ee1639SBarry Smith + pc - the preconditioner context 1263f1580f4eSBarry Smith - flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 126423ee1639SBarry Smith 12652b26979fSBarry Smith Level: intermediate 12662b26979fSBarry Smith 1267f1580f4eSBarry Smith Note: 1268f1580f4eSBarry Smith Normally if a matrix inside a `PC` changes the `PC` automatically updates itself using information from the changed matrix. This option 1269f1580f4eSBarry Smith prevents this. 1270f1580f4eSBarry Smith 1271f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCGetReusePreconditioner()`, `KSPSetReusePreconditioner()` 127223ee1639SBarry Smith @*/ 1273d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetReusePreconditioner(PC pc, PetscBool flag) 1274d71ae5a4SJacob Faibussowitsch { 127523ee1639SBarry Smith PetscFunctionBegin; 127623ee1639SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1277f9177522SStefano Zampini PetscValidLogicalCollectiveBool(pc, flag, 2); 127823ee1639SBarry Smith pc->reusepreconditioner = flag; 12793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12804b9ad928SBarry Smith } 12814b9ad928SBarry Smith 1282c60c7ad4SBarry Smith /*@ 1283f1580f4eSBarry Smith PCGetReusePreconditioner - Determines if the `PC` reuses the current preconditioner even if the operator in the preconditioner has changed. 1284c60c7ad4SBarry Smith 1285bba28a21SBarry Smith Not Collective 1286c60c7ad4SBarry Smith 1287c60c7ad4SBarry Smith Input Parameter: 1288c60c7ad4SBarry Smith . pc - the preconditioner context 1289c60c7ad4SBarry Smith 1290c60c7ad4SBarry Smith Output Parameter: 1291f1580f4eSBarry Smith . flag - `PETSC_TRUE` do not compute a new preconditioner, `PETSC_FALSE` do compute a new preconditioner 1292c60c7ad4SBarry Smith 1293d0418729SSatish Balay Level: intermediate 1294d0418729SSatish Balay 1295f1580f4eSBarry Smith .seealso: `PC`, `PCGetOperators()`, `MatZeroEntries()`, `PCSetReusePreconditioner()` 1296c60c7ad4SBarry Smith @*/ 1297d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetReusePreconditioner(PC pc, PetscBool *flag) 1298d71ae5a4SJacob Faibussowitsch { 1299c60c7ad4SBarry Smith PetscFunctionBegin; 1300c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 13014f572ea9SToby Isaac PetscAssertPointer(flag, 2); 1302c60c7ad4SBarry Smith *flag = pc->reusepreconditioner; 13033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1304c60c7ad4SBarry Smith } 1305c60c7ad4SBarry Smith 1306487a658cSBarry Smith /*@ 13074b9ad928SBarry Smith PCGetOperators - Gets the matrix associated with the linear system and 13084b9ad928SBarry Smith possibly a different one associated with the preconditioner. 13094b9ad928SBarry Smith 131020f4b53cSBarry Smith Not Collective, though parallel mats are returned if pc is parallel 13114b9ad928SBarry Smith 13124b9ad928SBarry Smith Input Parameter: 13134b9ad928SBarry Smith . pc - the preconditioner context 13144b9ad928SBarry Smith 13154b9ad928SBarry Smith Output Parameters: 1316e5d3d808SBarry Smith + Amat - the matrix defining the linear system 131723ee1639SBarry Smith - Pmat - the matrix from which the preconditioner is constructed, usually the same as Amat. 13184b9ad928SBarry Smith 13194b9ad928SBarry Smith Level: intermediate 13204b9ad928SBarry Smith 1321f1580f4eSBarry Smith Note: 132295452b02SPatrick Sanan Does not increase the reference count of the matrices, so you should not destroy them 1323298cc208SBarry Smith 1324f1580f4eSBarry Smith Alternative usage: If the operators have NOT been set with `KSPSetOperators()`/`PCSetOperators()` then the operators 1325f1580f4eSBarry Smith are created in `PC` and returned to the user. In this case, if both operators 132673950996SBarry Smith mat and pmat are requested, two DIFFERENT operators will be returned. If 132773950996SBarry Smith only one is requested both operators in the PC will be the same (i.e. as 1328f1580f4eSBarry Smith if one had called `KSPSetOperators()`/`PCSetOperators()` with the same argument for both Mats). 132973950996SBarry Smith The user must set the sizes of the returned matrices and their type etc just 1330f1580f4eSBarry Smith as if the user created them with `MatCreate()`. For example, 133173950996SBarry Smith 1332f1580f4eSBarry Smith .vb 1333f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,NULL); is equivalent to 1334f1580f4eSBarry Smith set size, type, etc of Amat 133573950996SBarry Smith 1336f1580f4eSBarry Smith MatCreate(comm,&mat); 1337f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Amat); 1338f1580f4eSBarry Smith PetscObjectDereference((PetscObject)mat); 1339f1580f4eSBarry Smith set size, type, etc of Amat 1340f1580f4eSBarry Smith .ve 134173950996SBarry Smith 134273950996SBarry Smith and 134373950996SBarry Smith 1344f1580f4eSBarry Smith .vb 1345f1580f4eSBarry Smith KSP/PCGetOperators(ksp/pc,&Amat,&Pmat); is equivalent to 1346f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 134773950996SBarry Smith 1348f1580f4eSBarry Smith MatCreate(comm,&Amat); 1349f1580f4eSBarry Smith MatCreate(comm,&Pmat); 1350f1580f4eSBarry Smith KSP/PCSetOperators(ksp/pc,Amat,Pmat); 1351f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Amat); 1352f1580f4eSBarry Smith PetscObjectDereference((PetscObject)Pmat); 1353f1580f4eSBarry Smith set size, type, etc of Amat and Pmat 1354f1580f4eSBarry Smith .ve 135573950996SBarry Smith 1356f1580f4eSBarry Smith The rationale for this support is so that when creating a `TS`, `SNES`, or `KSP` the hierarchy 1357b8d709abSRichard Tran Mills of underlying objects (i.e. `SNES`, `KSP`, `PC`, `Mat`) and their lifespans can be completely 1358f1580f4eSBarry Smith managed by the top most level object (i.e. the `TS`, `SNES`, or `KSP`). Another way to look 1359f1580f4eSBarry Smith at this is when you create a `SNES` you do not NEED to create a `KSP` and attach it to 1360f1580f4eSBarry Smith the `SNES` object (the `SNES` object manages it for you). Similarly when you create a KSP 1361f1580f4eSBarry Smith you do not need to attach a `PC` to it (the `KSP` object manages the `PC` object for you). 1362f1580f4eSBarry Smith Thus, why should YOU have to create the `Mat` and attach it to the `SNES`/`KSP`/`PC`, when 136373950996SBarry Smith it can be created for you? 136473950996SBarry Smith 1365f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperatorsSet()` 13664b9ad928SBarry Smith @*/ 1367d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperators(PC pc, Mat *Amat, Mat *Pmat) 1368d71ae5a4SJacob Faibussowitsch { 13694b9ad928SBarry Smith PetscFunctionBegin; 13700700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1371e5d3d808SBarry Smith if (Amat) { 137273950996SBarry Smith if (!pc->mat) { 13739fca8c71SStefano Zampini if (pc->pmat && !Pmat) { /* Pmat has been set, but user did not request it, so use for Amat */ 13749a4708feSJed Brown pc->mat = pc->pmat; 13759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 1376e5d3d808SBarry Smith } else { /* both Amat and Pmat are empty */ 13779566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->mat)); 1378e5d3d808SBarry Smith if (!Pmat) { /* user did NOT request Pmat, so make same as Amat */ 137973950996SBarry Smith pc->pmat = pc->mat; 13809566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 138173950996SBarry Smith } 138273950996SBarry Smith } 13839a4708feSJed Brown } 1384e5d3d808SBarry Smith *Amat = pc->mat; 138573950996SBarry Smith } 1386e5d3d808SBarry Smith if (Pmat) { 138773950996SBarry Smith if (!pc->pmat) { 1388e5d3d808SBarry Smith if (pc->mat && !Amat) { /* Amat has been set but was not requested, so use for pmat */ 13899a4708feSJed Brown pc->pmat = pc->mat; 13909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->pmat)); 13919a4708feSJed Brown } else { 13929566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc), &pc->pmat)); 1393e5d3d808SBarry Smith if (!Amat) { /* user did NOT request Amat, so make same as Pmat */ 139473950996SBarry Smith pc->mat = pc->pmat; 13959566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)pc->mat)); 139673950996SBarry Smith } 139773950996SBarry Smith } 13989a4708feSJed Brown } 1399e5d3d808SBarry Smith *Pmat = pc->pmat; 140073950996SBarry Smith } 14013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14024b9ad928SBarry Smith } 14034b9ad928SBarry Smith 1404906ed7ccSBarry Smith /*@C 1405906ed7ccSBarry Smith PCGetOperatorsSet - Determines if the matrix associated with the linear system and 1406f1580f4eSBarry Smith possibly a different one associated with the preconditioner have been set in the `PC`. 1407906ed7ccSBarry Smith 140820f4b53cSBarry Smith Not Collective, though the results on all processes should be the same 1409906ed7ccSBarry Smith 1410906ed7ccSBarry Smith Input Parameter: 1411906ed7ccSBarry Smith . pc - the preconditioner context 1412906ed7ccSBarry Smith 1413906ed7ccSBarry Smith Output Parameters: 1414906ed7ccSBarry Smith + mat - the matrix associated with the linear system was set 1415906ed7ccSBarry Smith - pmat - matrix associated with the preconditioner was set, usually the same 1416906ed7ccSBarry Smith 1417906ed7ccSBarry Smith Level: intermediate 1418906ed7ccSBarry Smith 1419f1580f4eSBarry Smith .seealso: `PC`, `PCSetOperators()`, `KSPGetOperators()`, `KSPSetOperators()`, `PCGetOperators()` 1420906ed7ccSBarry Smith @*/ 1421d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOperatorsSet(PC pc, PetscBool *mat, PetscBool *pmat) 1422d71ae5a4SJacob Faibussowitsch { 1423906ed7ccSBarry Smith PetscFunctionBegin; 14240700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1425906ed7ccSBarry Smith if (mat) *mat = (pc->mat) ? PETSC_TRUE : PETSC_FALSE; 1426906ed7ccSBarry Smith if (pmat) *pmat = (pc->pmat) ? PETSC_TRUE : PETSC_FALSE; 14273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1428906ed7ccSBarry Smith } 1429906ed7ccSBarry Smith 1430f39d8e23SSatish Balay /*@ 1431a4fd02acSBarry Smith PCFactorGetMatrix - Gets the factored matrix from the 1432f1580f4eSBarry Smith preconditioner context. This routine is valid only for the `PCLU`, 1433f1580f4eSBarry Smith `PCILU`, `PCCHOLESKY`, and `PCICC` methods. 14344b9ad928SBarry Smith 1435c3339decSBarry Smith Not Collective though mat is parallel if pc is parallel 14364b9ad928SBarry Smith 1437f1580f4eSBarry Smith Input Parameter: 14384b9ad928SBarry Smith . pc - the preconditioner context 14394b9ad928SBarry Smith 1440feefa0e1SJacob Faibussowitsch Output Parameters: 14414b9ad928SBarry Smith . mat - the factored matrix 14424b9ad928SBarry Smith 14434b9ad928SBarry Smith Level: advanced 14444b9ad928SBarry Smith 1445f1580f4eSBarry Smith Note: 144695452b02SPatrick Sanan Does not increase the reference count for the matrix so DO NOT destroy it 14479405f653SBarry Smith 1448f1580f4eSBarry Smith .seealso: `PC`, `PCLU`, `PCILU`, `PCCHOLESKY`, `PCICC` 14494b9ad928SBarry Smith @*/ 1450d71ae5a4SJacob Faibussowitsch PetscErrorCode PCFactorGetMatrix(PC pc, Mat *mat) 1451d71ae5a4SJacob Faibussowitsch { 14524b9ad928SBarry Smith PetscFunctionBegin; 14530700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14544f572ea9SToby Isaac PetscAssertPointer(mat, 2); 14557def7855SStefano Zampini PetscCall(PCFactorSetUpMatSolverType(pc)); 1456dbbe0bcdSBarry Smith PetscUseTypeMethod(pc, getfactoredmatrix, mat); 14573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14584b9ad928SBarry Smith } 14594b9ad928SBarry Smith 14604b9ad928SBarry Smith /*@C 14614b9ad928SBarry Smith PCSetOptionsPrefix - Sets the prefix used for searching for all 1462f1580f4eSBarry Smith `PC` options in the database. 14634b9ad928SBarry Smith 1464c3339decSBarry Smith Logically Collective 14654b9ad928SBarry Smith 14664b9ad928SBarry Smith Input Parameters: 14674b9ad928SBarry Smith + pc - the preconditioner context 1468f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 14694b9ad928SBarry Smith 1470f1580f4eSBarry Smith Note: 14714b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 14724b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 14734b9ad928SBarry Smith hyphen. 14744b9ad928SBarry Smith 14754b9ad928SBarry Smith Level: advanced 14764b9ad928SBarry Smith 1477f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCAppendOptionsPrefix()`, `PCGetOptionsPrefix()` 14784b9ad928SBarry Smith @*/ 1479d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetOptionsPrefix(PC pc, const char prefix[]) 1480d71ae5a4SJacob Faibussowitsch { 14814b9ad928SBarry Smith PetscFunctionBegin; 14820700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 14839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pc, prefix)); 14843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14854b9ad928SBarry Smith } 14864b9ad928SBarry Smith 14874b9ad928SBarry Smith /*@C 14884b9ad928SBarry Smith PCAppendOptionsPrefix - Appends to the prefix used for searching for all 1489f1580f4eSBarry Smith `PC` options in the database. 14904b9ad928SBarry Smith 1491c3339decSBarry Smith Logically Collective 14924b9ad928SBarry Smith 14934b9ad928SBarry Smith Input Parameters: 14944b9ad928SBarry Smith + pc - the preconditioner context 1495f1580f4eSBarry Smith - prefix - the prefix string to prepend to all `PC` option requests 14964b9ad928SBarry Smith 1497f1580f4eSBarry Smith Note: 14984b9ad928SBarry Smith A hyphen (-) must NOT be given at the beginning of the prefix name. 14994b9ad928SBarry Smith The first character of all runtime options is AUTOMATICALLY the 15004b9ad928SBarry Smith hyphen. 15014b9ad928SBarry Smith 15024b9ad928SBarry Smith Level: advanced 15034b9ad928SBarry Smith 1504f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCGetOptionsPrefix()` 15054b9ad928SBarry Smith @*/ 1506d71ae5a4SJacob Faibussowitsch PetscErrorCode PCAppendOptionsPrefix(PC pc, const char prefix[]) 1507d71ae5a4SJacob Faibussowitsch { 15084b9ad928SBarry Smith PetscFunctionBegin; 15090700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15109566063dSJacob Faibussowitsch PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pc, prefix)); 15113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15124b9ad928SBarry Smith } 15134b9ad928SBarry Smith 15144b9ad928SBarry Smith /*@C 15154b9ad928SBarry Smith PCGetOptionsPrefix - Gets the prefix used for searching for all 15164b9ad928SBarry Smith PC options in the database. 15174b9ad928SBarry Smith 15184b9ad928SBarry Smith Not Collective 15194b9ad928SBarry Smith 1520f1580f4eSBarry Smith Input Parameter: 15214b9ad928SBarry Smith . pc - the preconditioner context 15224b9ad928SBarry Smith 1523f1580f4eSBarry Smith Output Parameter: 15244b9ad928SBarry Smith . prefix - pointer to the prefix string used, is returned 15254b9ad928SBarry Smith 1526feefa0e1SJacob Faibussowitsch Fortran Notes: 1527f1580f4eSBarry Smith The user should pass in a string 'prefix' of 15284b9ad928SBarry Smith sufficient length to hold the prefix. 15294b9ad928SBarry Smith 15304b9ad928SBarry Smith Level: advanced 15314b9ad928SBarry Smith 1532f1580f4eSBarry Smith .seealso: `PC`, `PCSetFromOptions`, `PCSetOptionsPrefix()`, `PCAppendOptionsPrefix()` 15334b9ad928SBarry Smith @*/ 1534d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetOptionsPrefix(PC pc, const char *prefix[]) 1535d71ae5a4SJacob Faibussowitsch { 15364b9ad928SBarry Smith PetscFunctionBegin; 15370700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15384f572ea9SToby Isaac PetscAssertPointer(prefix, 2); 15399566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc, prefix)); 15403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15414b9ad928SBarry Smith } 15424b9ad928SBarry Smith 15438066bbecSBarry Smith /* 15448066bbecSBarry Smith Indicates the right hand side will be changed by KSPSolve(), this occurs for a few 15458066bbecSBarry Smith preconditioners including BDDC and Eisentat that transform the equations before applying 15468066bbecSBarry Smith the Krylov methods 15478066bbecSBarry Smith */ 1548d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC pc, PetscBool *change) 1549d71ae5a4SJacob Faibussowitsch { 1550a06fd7f2SStefano Zampini PetscFunctionBegin; 1551a06fd7f2SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15524f572ea9SToby Isaac PetscAssertPointer(change, 2); 1553a06fd7f2SStefano Zampini *change = PETSC_FALSE; 1554cac4c232SBarry Smith PetscTryMethod(pc, "PCPreSolveChangeRHS_C", (PC, PetscBool *), (pc, change)); 15553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1556a06fd7f2SStefano Zampini } 1557a06fd7f2SStefano Zampini 15584b9ad928SBarry Smith /*@ 15594b9ad928SBarry Smith PCPreSolve - Optional pre-solve phase, intended for any 15604b9ad928SBarry Smith preconditioner-specific actions that must be performed before 15614b9ad928SBarry Smith the iterative solve itself. 15624b9ad928SBarry Smith 1563c3339decSBarry Smith Collective 15644b9ad928SBarry Smith 15654b9ad928SBarry Smith Input Parameters: 15664b9ad928SBarry Smith + pc - the preconditioner context 15674b9ad928SBarry Smith - ksp - the Krylov subspace context 15684b9ad928SBarry Smith 15694b9ad928SBarry Smith Level: developer 15704b9ad928SBarry Smith 1571feefa0e1SJacob Faibussowitsch Example Usage: 15724b9ad928SBarry Smith .vb 15734b9ad928SBarry Smith PCPreSolve(pc,ksp); 157423ce1328SBarry Smith KSPSolve(ksp,b,x); 15754b9ad928SBarry Smith PCPostSolve(pc,ksp); 15764b9ad928SBarry Smith .ve 15774b9ad928SBarry Smith 15784b9ad928SBarry Smith Notes: 1579f1580f4eSBarry Smith The pre-solve phase is distinct from the `PCSetUp()` phase. 15804b9ad928SBarry Smith 1581f1580f4eSBarry Smith `KSPSolve()` calls this directly, so is rarely called by the user. 15824b9ad928SBarry Smith 1583f1580f4eSBarry Smith .seealso: `PC`, `PCPostSolve()` 15844b9ad928SBarry Smith @*/ 1585d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPreSolve(PC pc, KSP ksp) 1586d71ae5a4SJacob Faibussowitsch { 15874b9ad928SBarry Smith Vec x, rhs; 15884b9ad928SBarry Smith 15894b9ad928SBarry Smith PetscFunctionBegin; 15900700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 15910700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1592d9a03883SBarry Smith pc->presolvedone++; 15937827d75bSBarry Smith PetscCheck(pc->presolvedone <= 2, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot embed PCPreSolve() more than twice"); 15949566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 15959566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 15964b9ad928SBarry Smith 1597dbbe0bcdSBarry Smith if (pc->ops->presolve) PetscUseTypeMethod(pc, presolve, ksp, rhs, x); 15981baa6e33SBarry Smith else if (pc->presolve) PetscCall((pc->presolve)(pc, ksp)); 15993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16004b9ad928SBarry Smith } 16014b9ad928SBarry Smith 1602f560b561SHong Zhang /*@C 1603f1580f4eSBarry Smith PCSetPreSolve - Sets function used by `PCPreSolve()` which is intended for any 1604f560b561SHong Zhang preconditioner-specific actions that must be performed before 1605f560b561SHong Zhang the iterative solve itself. 1606f560b561SHong Zhang 1607c3339decSBarry Smith Logically Collective 1608f560b561SHong Zhang 1609f560b561SHong Zhang Input Parameters: 1610f560b561SHong Zhang + pc - the preconditioner object 1611f560b561SHong Zhang - presolve - the function to call before the solve 1612f560b561SHong Zhang 161320f4b53cSBarry Smith Calling sequence of `presolve`: 161420f4b53cSBarry Smith + pc - the `PC` context 161520f4b53cSBarry Smith - ksp - the `KSP` context 1616f560b561SHong Zhang 1617f560b561SHong Zhang Level: developer 1618f560b561SHong Zhang 1619db781477SPatrick Sanan .seealso: `PC`, `PCSetUp()`, `PCPreSolve()` 1620f560b561SHong Zhang @*/ 162104c3f3b8SBarry Smith PetscErrorCode PCSetPreSolve(PC pc, PetscErrorCode (*presolve)(PC pc, KSP ksp)) 1622d71ae5a4SJacob Faibussowitsch { 1623f560b561SHong Zhang PetscFunctionBegin; 1624f560b561SHong Zhang PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1625f560b561SHong Zhang pc->presolve = presolve; 16263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1627f560b561SHong Zhang } 1628f560b561SHong Zhang 16294b9ad928SBarry Smith /*@ 16304b9ad928SBarry Smith PCPostSolve - Optional post-solve phase, intended for any 16314b9ad928SBarry Smith preconditioner-specific actions that must be performed after 16324b9ad928SBarry Smith the iterative solve itself. 16334b9ad928SBarry Smith 1634c3339decSBarry Smith Collective 16354b9ad928SBarry Smith 16364b9ad928SBarry Smith Input Parameters: 16374b9ad928SBarry Smith + pc - the preconditioner context 16384b9ad928SBarry Smith - ksp - the Krylov subspace context 16394b9ad928SBarry Smith 1640feefa0e1SJacob Faibussowitsch Example Usage: 16414b9ad928SBarry Smith .vb 16424b9ad928SBarry Smith PCPreSolve(pc,ksp); 164323ce1328SBarry Smith KSPSolve(ksp,b,x); 16444b9ad928SBarry Smith PCPostSolve(pc,ksp); 16454b9ad928SBarry Smith .ve 16464b9ad928SBarry Smith 16474b9ad928SBarry Smith Note: 1648f1580f4eSBarry Smith `KSPSolve()` calls this routine directly, so it is rarely called by the user. 16494b9ad928SBarry Smith 16504b9ad928SBarry Smith Level: developer 16514b9ad928SBarry Smith 1652f1580f4eSBarry Smith .seealso: `PC`, `PCSetPostSolve()`, `PCSetPresolve()`, `PCPreSolve()`, `KSPSolve()` 16534b9ad928SBarry Smith @*/ 1654d71ae5a4SJacob Faibussowitsch PetscErrorCode PCPostSolve(PC pc, KSP ksp) 1655d71ae5a4SJacob Faibussowitsch { 16564b9ad928SBarry Smith Vec x, rhs; 16574b9ad928SBarry Smith 16584b9ad928SBarry Smith PetscFunctionBegin; 16590700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16600700a824SBarry Smith PetscValidHeaderSpecific(ksp, KSP_CLASSID, 2); 1661d9a03883SBarry Smith pc->presolvedone--; 16629566063dSJacob Faibussowitsch PetscCall(KSPGetSolution(ksp, &x)); 16639566063dSJacob Faibussowitsch PetscCall(KSPGetRhs(ksp, &rhs)); 1664dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, postsolve, ksp, rhs, x); 16653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16664b9ad928SBarry Smith } 16674b9ad928SBarry Smith 166855849f57SBarry Smith /*@C 1669f1580f4eSBarry Smith PCLoad - Loads a `PC` that has been stored in binary with `PCView()`. 167055849f57SBarry Smith 1671c3339decSBarry Smith Collective 167255849f57SBarry Smith 167355849f57SBarry Smith Input Parameters: 1674f1580f4eSBarry Smith + newdm - the newly loaded `PC`, this needs to have been created with `PCCreate()` or 1675f1580f4eSBarry Smith some related function before a call to `PCLoad()`. 1676f1580f4eSBarry Smith - viewer - binary file viewer, obtained from `PetscViewerBinaryOpen()` 167755849f57SBarry Smith 167855849f57SBarry Smith Level: intermediate 167955849f57SBarry Smith 1680f1580f4eSBarry Smith Note: 168155849f57SBarry Smith The type is determined by the data in the file, any type set into the PC before this call is ignored. 168255849f57SBarry Smith 1683f1580f4eSBarry Smith .seealso: `PC`, `PetscViewerBinaryOpen()`, `PCView()`, `MatLoad()`, `VecLoad()` 168455849f57SBarry Smith @*/ 1685d71ae5a4SJacob Faibussowitsch PetscErrorCode PCLoad(PC newdm, PetscViewer viewer) 1686d71ae5a4SJacob Faibussowitsch { 168755849f57SBarry Smith PetscBool isbinary; 1688060da220SMatthew G. Knepley PetscInt classid; 168955849f57SBarry Smith char type[256]; 169055849f57SBarry Smith 169155849f57SBarry Smith PetscFunctionBegin; 169255849f57SBarry Smith PetscValidHeaderSpecific(newdm, PC_CLASSID, 1); 169355849f57SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 16949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 169528b400f6SJacob Faibussowitsch PetscCheck(isbinary, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 169655849f57SBarry Smith 16979566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, &classid, 1, NULL, PETSC_INT)); 16987827d75bSBarry Smith PetscCheck(classid == PC_FILE_CLASSID, PetscObjectComm((PetscObject)newdm), PETSC_ERR_ARG_WRONG, "Not PC next in file"); 16999566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(viewer, type, 256, NULL, PETSC_CHAR)); 17009566063dSJacob Faibussowitsch PetscCall(PCSetType(newdm, type)); 1701dbbe0bcdSBarry Smith PetscTryTypeMethod(newdm, load, viewer); 17023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 170355849f57SBarry Smith } 170455849f57SBarry Smith 17059804daf3SBarry Smith #include <petscdraw.h> 1706e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1707e04113cfSBarry Smith #include <petscviewersaws.h> 17080acecf5bSBarry Smith #endif 1709fe2efc57SMark 1710fe2efc57SMark /*@C 1711f1580f4eSBarry Smith PCViewFromOptions - View from the `PC` based on options in the database 1712fe2efc57SMark 1713c3339decSBarry Smith Collective 1714fe2efc57SMark 1715fe2efc57SMark Input Parameters: 171620f4b53cSBarry Smith + A - the `PC` context 1717f1580f4eSBarry Smith . obj - Optional object that provides the options prefix 1718736c3998SJose E. Roman - name - command line option 1719fe2efc57SMark 1720fe2efc57SMark Level: intermediate 1721f1580f4eSBarry Smith 1722db781477SPatrick Sanan .seealso: `PC`, `PCView`, `PetscObjectViewFromOptions()`, `PCCreate()` 1723fe2efc57SMark @*/ 1724d71ae5a4SJacob Faibussowitsch PetscErrorCode PCViewFromOptions(PC A, PetscObject obj, const char name[]) 1725d71ae5a4SJacob Faibussowitsch { 1726fe2efc57SMark PetscFunctionBegin; 1727fe2efc57SMark PetscValidHeaderSpecific(A, PC_CLASSID, 1); 17289566063dSJacob Faibussowitsch PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name)); 17293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1730fe2efc57SMark } 1731fe2efc57SMark 17324b9ad928SBarry Smith /*@C 1733f1580f4eSBarry Smith PCView - Prints information about the `PC` 17344b9ad928SBarry Smith 1735c3339decSBarry Smith Collective 17364b9ad928SBarry Smith 17374b9ad928SBarry Smith Input Parameters: 1738feefa0e1SJacob Faibussowitsch + pc - the `PC` context 17394b9ad928SBarry Smith - viewer - optional visualization context 17404b9ad928SBarry Smith 174120f4b53cSBarry Smith Level: developer 174220f4b53cSBarry Smith 1743f1580f4eSBarry Smith Notes: 17444b9ad928SBarry Smith The available visualization contexts include 1745f1580f4eSBarry Smith + `PETSC_VIEWER_STDOUT_SELF` - standard output (default) 1746f1580f4eSBarry Smith - `PETSC_VIEWER_STDOUT_WORLD` - synchronized standard 17474b9ad928SBarry Smith output where only the first processor opens 17484b9ad928SBarry Smith the file. All other processors send their 17494b9ad928SBarry Smith data to the first processor to print. 17504b9ad928SBarry Smith 17514b9ad928SBarry Smith The user can open an alternative visualization contexts with 1752f1580f4eSBarry Smith `PetscViewerASCIIOpen()` (output to a specified file). 17534b9ad928SBarry Smith 1754f1580f4eSBarry Smith .seealso: `PC`, `PetscViewer`, `KSPView()`, `PetscViewerASCIIOpen()` 17554b9ad928SBarry Smith @*/ 1756d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView(PC pc, PetscViewer viewer) 1757d71ae5a4SJacob Faibussowitsch { 175819fd82e9SBarry Smith PCType cstr; 17596cd81132SPierre Jolivet PetscViewerFormat format; 17606cd81132SPierre Jolivet PetscBool iascii, isstring, isbinary, isdraw, pop = PETSC_FALSE; 1761e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1762536b137fSBarry Smith PetscBool issaws; 17630acecf5bSBarry Smith #endif 17644b9ad928SBarry Smith 17654b9ad928SBarry Smith PetscFunctionBegin; 17660700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 176748a46eb9SPierre Jolivet if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)pc), &viewer)); 17680700a824SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1769c9780b6fSBarry Smith PetscCheckSameComm(pc, 1, viewer, 2); 17704b9ad928SBarry Smith 17719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 17729566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 17739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary)); 17749566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 1775e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 17769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSAWS, &issaws)); 17770acecf5bSBarry Smith #endif 1778219b1045SBarry Smith 177932077d6dSBarry Smith if (iascii) { 17809566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)pc, viewer)); 178148a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(PetscViewerASCIIPrintf(viewer, " PC has not been set up so information may be incomplete\n")); 17829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 1783dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 17849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 1785834dbeb0SBarry Smith if (pc->mat) { 17866cd81132SPierre Jolivet PetscCall(PetscViewerGetFormat(viewer, &format)); 17876cd81132SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 17889566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_INFO)); 17896cd81132SPierre Jolivet pop = PETSC_TRUE; 17906cd81132SPierre Jolivet } 17914b9ad928SBarry Smith if (pc->pmat == pc->mat) { 17929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix = precond matrix:\n")); 17939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 17949566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 17959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 17964b9ad928SBarry Smith } else { 1797834dbeb0SBarry Smith if (pc->pmat) { 17989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix followed by preconditioner matrix:\n")); 17994b9ad928SBarry Smith } else { 18009566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " linear system matrix:\n")); 18014b9ad928SBarry Smith } 18029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 18039566063dSJacob Faibussowitsch PetscCall(MatView(pc->mat, viewer)); 18049566063dSJacob Faibussowitsch if (pc->pmat) PetscCall(MatView(pc->pmat, viewer)); 18059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 18064b9ad928SBarry Smith } 18076cd81132SPierre Jolivet if (pop) PetscCall(PetscViewerPopFormat(viewer)); 18084b9ad928SBarry Smith } 18094b9ad928SBarry Smith } else if (isstring) { 18109566063dSJacob Faibussowitsch PetscCall(PCGetType(pc, &cstr)); 18119566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " PCType: %-7.7s", cstr)); 1812dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18139566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18149566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18155b0b0462SBarry Smith } else if (isbinary) { 181655849f57SBarry Smith PetscInt classid = PC_FILE_CLASSID; 181755849f57SBarry Smith MPI_Comm comm; 181855849f57SBarry Smith PetscMPIInt rank; 181955849f57SBarry Smith char type[256]; 182055849f57SBarry Smith 18219566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 18229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 1823dd400576SPatrick Sanan if (rank == 0) { 18249566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, &classid, 1, PETSC_INT)); 18259566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(type, ((PetscObject)pc)->type_name, 256)); 18269566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(viewer, type, 256, PETSC_CHAR)); 182755849f57SBarry Smith } 1828dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 1829219b1045SBarry Smith } else if (isdraw) { 1830219b1045SBarry Smith PetscDraw draw; 1831d9884438SBarry Smith char str[25]; 183289fd9fafSBarry Smith PetscReal x, y, bottom, h; 1833d9884438SBarry Smith PetscInt n; 1834219b1045SBarry Smith 18359566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 18369566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 18371d840656SPeter Brune if (pc->mat) { 18389566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat, &n, NULL)); 183963a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s (%" PetscInt_FMT ")", ((PetscObject)pc)->type_name, n)); 18401d840656SPeter Brune } else { 18419566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(str, 25, "PC: %s", ((PetscObject)pc)->type_name)); 18421d840656SPeter Brune } 18439566063dSJacob Faibussowitsch PetscCall(PetscDrawStringBoxed(draw, x, y, PETSC_DRAW_RED, PETSC_DRAW_BLACK, str, NULL, &h)); 184489fd9fafSBarry Smith bottom = y - h; 18459566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 1846dbbe0bcdSBarry Smith PetscTryTypeMethod(pc, view, viewer); 18479566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 1848e04113cfSBarry Smith #if defined(PETSC_HAVE_SAWS) 1849536b137fSBarry Smith } else if (issaws) { 1850d45a07a7SBarry Smith PetscMPIInt rank; 1851d45a07a7SBarry Smith 18529566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)pc)); 18539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 185448a46eb9SPierre Jolivet if (!((PetscObject)pc)->amsmem && rank == 0) PetscCall(PetscObjectViewSAWs((PetscObject)pc, viewer)); 18559566063dSJacob Faibussowitsch if (pc->mat) PetscCall(MatView(pc->mat, viewer)); 18569566063dSJacob Faibussowitsch if (pc->pmat && pc->pmat != pc->mat) PetscCall(MatView(pc->pmat, viewer)); 18570acecf5bSBarry Smith #endif 18584b9ad928SBarry Smith } 18593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18604b9ad928SBarry Smith } 18614b9ad928SBarry Smith 18624b9ad928SBarry Smith /*@C 18631c84c290SBarry Smith PCRegister - Adds a method to the preconditioner package. 18641c84c290SBarry Smith 18651c84c290SBarry Smith Not collective 18661c84c290SBarry Smith 18671c84c290SBarry Smith Input Parameters: 186820f4b53cSBarry Smith + sname - name of a new user-defined solver 186920f4b53cSBarry Smith - function - routine to create method context 18701c84c290SBarry Smith 1871feefa0e1SJacob Faibussowitsch Example Usage: 18721c84c290SBarry Smith .vb 1873bdf89e91SBarry Smith PCRegister("my_solver", MySolverCreate); 18741c84c290SBarry Smith .ve 18751c84c290SBarry Smith 18761c84c290SBarry Smith Then, your solver can be chosen with the procedural interface via 18771c84c290SBarry Smith $ PCSetType(pc, "my_solver") 18781c84c290SBarry Smith or at runtime via the option 18791c84c290SBarry Smith $ -pc_type my_solver 18804b9ad928SBarry Smith 18814b9ad928SBarry Smith Level: advanced 18821c84c290SBarry Smith 188320f4b53cSBarry Smith Note: 188420f4b53cSBarry Smith `PCRegister()` may be called multiple times to add several user-defined preconditioners. 188520f4b53cSBarry Smith 1886f1580f4eSBarry Smith .seealso: `PC`, `PCType`, `PCRegisterAll()` 18874b9ad928SBarry Smith @*/ 1888d71ae5a4SJacob Faibussowitsch PetscErrorCode PCRegister(const char sname[], PetscErrorCode (*function)(PC)) 1889d71ae5a4SJacob Faibussowitsch { 18904b9ad928SBarry Smith PetscFunctionBegin; 18919566063dSJacob Faibussowitsch PetscCall(PCInitializePackage()); 18929566063dSJacob Faibussowitsch PetscCall(PetscFunctionListAdd(&PCList, sname, function)); 18933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18944b9ad928SBarry Smith } 18954b9ad928SBarry Smith 1896d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMult_PC(Mat A, Vec X, Vec Y) 1897d71ae5a4SJacob Faibussowitsch { 1898186a3e20SStefano Zampini PC pc; 1899186a3e20SStefano Zampini 1900186a3e20SStefano Zampini PetscFunctionBegin; 19019566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(A, &pc)); 19029566063dSJacob Faibussowitsch PetscCall(PCApply(pc, X, Y)); 19033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1904186a3e20SStefano Zampini } 1905186a3e20SStefano Zampini 1906feefa0e1SJacob Faibussowitsch /*@C 19070bacdadaSStefano Zampini PCComputeOperator - Computes the explicit preconditioned operator. 19084b9ad928SBarry Smith 1909c3339decSBarry Smith Collective 19104b9ad928SBarry Smith 1911d8d19677SJose E. Roman Input Parameters: 1912186a3e20SStefano Zampini + pc - the preconditioner object 1913186a3e20SStefano Zampini - mattype - the matrix type to be used for the operator 19144b9ad928SBarry Smith 19154b9ad928SBarry Smith Output Parameter: 1916a5b23f4aSJose E. Roman . mat - the explicit preconditioned operator 19174b9ad928SBarry Smith 191820f4b53cSBarry Smith Level: advanced 191920f4b53cSBarry Smith 1920f1580f4eSBarry Smith Note: 1921186a3e20SStefano Zampini This computation is done by applying the operators to columns of the identity matrix. 1922186a3e20SStefano Zampini This routine is costly in general, and is recommended for use only with relatively small systems. 1923186a3e20SStefano Zampini Currently, this routine uses a dense matrix format when mattype == NULL 19244b9ad928SBarry Smith 1925f1580f4eSBarry Smith .seealso: `PC`, `KSPComputeOperator()`, `MatType` 19264b9ad928SBarry Smith @*/ 1927d71ae5a4SJacob Faibussowitsch PetscErrorCode PCComputeOperator(PC pc, MatType mattype, Mat *mat) 1928d71ae5a4SJacob Faibussowitsch { 1929186a3e20SStefano Zampini PetscInt N, M, m, n; 1930186a3e20SStefano Zampini Mat A, Apc; 19314b9ad928SBarry Smith 19324b9ad928SBarry Smith PetscFunctionBegin; 19330700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 19344f572ea9SToby Isaac PetscAssertPointer(mat, 3); 19359566063dSJacob Faibussowitsch PetscCall(PCGetOperators(pc, &A, NULL)); 19369566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &m, &n)); 19379566063dSJacob Faibussowitsch PetscCall(MatGetSize(A, &M, &N)); 19389566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PetscObjectComm((PetscObject)pc), m, n, M, N, pc, &Apc)); 19399566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(Apc, MATOP_MULT, (void (*)(void))MatMult_PC)); 19409566063dSJacob Faibussowitsch PetscCall(MatComputeOperator(Apc, mattype, mat)); 19419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Apc)); 19423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19434b9ad928SBarry Smith } 19444b9ad928SBarry Smith 19456c237d78SBarry Smith /*@ 19466c237d78SBarry Smith PCSetCoordinates - sets the coordinates of all the nodes on the local process 19476c237d78SBarry Smith 1948c3339decSBarry Smith Collective 19496c237d78SBarry Smith 19506c237d78SBarry Smith Input Parameters: 19516c237d78SBarry Smith + pc - the solver context 19526c237d78SBarry Smith . dim - the dimension of the coordinates 1, 2, or 3 195314893cbeSStefano Zampini . nloc - the blocked size of the coordinates array 195414893cbeSStefano Zampini - coords - the coordinates array 19556c237d78SBarry Smith 19566c237d78SBarry Smith Level: intermediate 19576c237d78SBarry Smith 1958f1580f4eSBarry Smith Note: 195920f4b53cSBarry Smith `coords` is an array of the dim coordinates for the nodes on 196020f4b53cSBarry Smith the local processor, of size `dim`*`nloc`. 196114893cbeSStefano Zampini If there are 108 equation on a processor 19626c237d78SBarry Smith for a displacement finite element discretization of elasticity (so 196314893cbeSStefano Zampini that there are nloc = 36 = 108/3 nodes) then the array must have 108 19646c237d78SBarry Smith double precision values (ie, 3 * 36). These x y z coordinates 19656c237d78SBarry Smith should be ordered for nodes 0 to N-1 like so: [ 0.x, 0.y, 0.z, 1.x, 19666c237d78SBarry Smith ... , N-1.z ]. 19676c237d78SBarry Smith 1968f1580f4eSBarry Smith .seealso: `PC`, `MatSetNearNullSpace()` 19696c237d78SBarry Smith @*/ 1970d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetCoordinates(PC pc, PetscInt dim, PetscInt nloc, PetscReal coords[]) 1971d71ae5a4SJacob Faibussowitsch { 19726c237d78SBarry Smith PetscFunctionBegin; 197332954da3SStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 197432954da3SStefano Zampini PetscValidLogicalCollectiveInt(pc, dim, 2); 197522794d57SStefano Zampini PetscTryMethod(pc, "PCSetCoordinates_C", (PC, PetscInt, PetscInt, PetscReal[]), (pc, dim, nloc, coords)); 19763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19776c237d78SBarry Smith } 1978fd2dd295SFande Kong 1979fd2dd295SFande Kong /*@ 1980fd2dd295SFande Kong PCGetInterpolations - Gets interpolation matrices for all levels (except level 0) 1981fd2dd295SFande Kong 1982c3339decSBarry Smith Logically Collective 1983fd2dd295SFande Kong 1984d8d19677SJose E. Roman Input Parameter: 1985d8d19677SJose E. Roman . pc - the precondition context 1986fd2dd295SFande Kong 1987d8d19677SJose E. Roman Output Parameters: 1988d8d19677SJose E. Roman + num_levels - the number of levels 1989d8d19677SJose E. Roman - interpolations - the interpolation matrices (size of num_levels-1) 1990fd2dd295SFande Kong 1991fd2dd295SFande Kong Level: advanced 1992fd2dd295SFande Kong 1993feefa0e1SJacob Faibussowitsch Developer Notes: 1994f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 1995fd2dd295SFande Kong 1996f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetInterpolation()`, `PCGetCoarseOperators()` 1997fd2dd295SFande Kong @*/ 1998d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetInterpolations(PC pc, PetscInt *num_levels, Mat *interpolations[]) 1999d71ae5a4SJacob Faibussowitsch { 2000fd2dd295SFande Kong PetscFunctionBegin; 2001fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20024f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20034f572ea9SToby Isaac PetscAssertPointer(interpolations, 3); 2004cac4c232SBarry Smith PetscUseMethod(pc, "PCGetInterpolations_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, interpolations)); 20053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2006fd2dd295SFande Kong } 2007fd2dd295SFande Kong 2008fd2dd295SFande Kong /*@ 2009fd2dd295SFande Kong PCGetCoarseOperators - Gets coarse operator matrices for all levels (except the finest level) 2010fd2dd295SFande Kong 2011c3339decSBarry Smith Logically Collective 2012fd2dd295SFande Kong 2013d8d19677SJose E. Roman Input Parameter: 2014d8d19677SJose E. Roman . pc - the precondition context 2015fd2dd295SFande Kong 2016d8d19677SJose E. Roman Output Parameters: 2017d8d19677SJose E. Roman + num_levels - the number of levels 2018d8d19677SJose E. Roman - coarseOperators - the coarse operator matrices (size of num_levels-1) 2019fd2dd295SFande Kong 2020fd2dd295SFande Kong Level: advanced 2021fd2dd295SFande Kong 2022feefa0e1SJacob Faibussowitsch Developer Notes: 2023f1580f4eSBarry Smith Why is this here instead of in `PCMG` etc? 2024fd2dd295SFande Kong 2025f1580f4eSBarry Smith .seealso: `PC`, `PCMG`, `PCMGGetRestriction()`, `PCMGSetInterpolation()`, `PCMGGetRScale()`, `PCMGGetInterpolation()`, `PCGetInterpolations()` 2026fd2dd295SFande Kong @*/ 2027d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGetCoarseOperators(PC pc, PetscInt *num_levels, Mat *coarseOperators[]) 2028d71ae5a4SJacob Faibussowitsch { 2029fd2dd295SFande Kong PetscFunctionBegin; 2030fd2dd295SFande Kong PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 20314f572ea9SToby Isaac PetscAssertPointer(num_levels, 2); 20324f572ea9SToby Isaac PetscAssertPointer(coarseOperators, 3); 2033cac4c232SBarry Smith PetscUseMethod(pc, "PCGetCoarseOperators_C", (PC, PetscInt *, Mat *[]), (pc, num_levels, coarseOperators)); 20343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2035fd2dd295SFande Kong } 2036