1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 24b9ad928SBarry Smith 3345a4b08SToby Isaac typedef enum { 4*b4f8a55aSStefano Zampini PCMATOP_UNSPECIFIED, 5345a4b08SToby Isaac PCMATOP_MULT, 6345a4b08SToby Isaac PCMATOP_MULT_TRANSPOSE, 7345a4b08SToby Isaac PCMATOP_MULT_HERMITIAN_TRANSPOSE, 8345a4b08SToby Isaac PCMATOP_SOLVE, 9345a4b08SToby Isaac PCMATOP_SOLVE_TRANSPOSE, 10345a4b08SToby Isaac } PCMatOperation; 11345a4b08SToby Isaac 12*b4f8a55aSStefano Zampini const char *const PCMatOpTypes[] = {"Unspecified", "Mult", "MultTranspose", "MultHermitianTranspose", "Solve", "SolveTranspose", NULL}; 13345a4b08SToby Isaac 14345a4b08SToby Isaac typedef struct _PCMAT { 15345a4b08SToby Isaac PCMatOperation apply; 16345a4b08SToby Isaac } PC_Mat; 17345a4b08SToby Isaac 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Mat(PC pc, Vec x, Vec y) 19d71ae5a4SJacob Faibussowitsch { 20345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 21345a4b08SToby Isaac 224b9ad928SBarry Smith PetscFunctionBegin; 23345a4b08SToby Isaac switch (pcmat->apply) { 24345a4b08SToby Isaac case PCMATOP_MULT: 259566063dSJacob Faibussowitsch PetscCall(MatMult(pc->pmat, x, y)); 26345a4b08SToby Isaac break; 27345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 28345a4b08SToby Isaac PetscCall(MatMultTranspose(pc->pmat, x, y)); 29345a4b08SToby Isaac break; 30345a4b08SToby Isaac case PCMATOP_SOLVE: 31345a4b08SToby Isaac PetscCall(MatSolve(pc->pmat, x, y)); 32345a4b08SToby Isaac break; 33345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 34345a4b08SToby Isaac PetscCall(MatSolveTranspose(pc->pmat, x, y)); 35345a4b08SToby Isaac break; 36345a4b08SToby Isaac case PCMATOP_MULT_HERMITIAN_TRANSPOSE: 37345a4b08SToby Isaac PetscCall(MatMultHermitianTranspose(pc->pmat, x, y)); 38345a4b08SToby Isaac break; 39*b4f8a55aSStefano Zampini default: 40*b4f8a55aSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported %s case", PCMatOpTypes[pcmat->apply]); 41*b4f8a55aSStefano Zampini } 42*b4f8a55aSStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 43*b4f8a55aSStefano Zampini } 44*b4f8a55aSStefano Zampini 45*b4f8a55aSStefano Zampini static PetscErrorCode PCSetUp_Mat(PC pc) 46*b4f8a55aSStefano Zampini { 47*b4f8a55aSStefano Zampini PC_Mat *pcmat = (PC_Mat *)pc->data; 48*b4f8a55aSStefano Zampini 49*b4f8a55aSStefano Zampini PetscFunctionBegin; 50*b4f8a55aSStefano Zampini if (pcmat->apply == PCMATOP_UNSPECIFIED) { 51*b4f8a55aSStefano Zampini PetscBool hassolve; 52*b4f8a55aSStefano Zampini PetscCall(MatHasOperation(pc->pmat, MATOP_SOLVE, &hassolve)); 53*b4f8a55aSStefano Zampini if (hassolve) pcmat->apply = PCMATOP_SOLVE; 54*b4f8a55aSStefano Zampini else pcmat->apply = PCMATOP_MULT; 55345a4b08SToby Isaac } 563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 574b9ad928SBarry Smith } 584b9ad928SBarry Smith 59d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Mat(PC pc, Mat X, Mat Y) 60d71ae5a4SJacob Faibussowitsch { 61345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 62*b4f8a55aSStefano Zampini Mat W; 63345a4b08SToby Isaac 647b6e2003SPierre Jolivet PetscFunctionBegin; 65345a4b08SToby Isaac switch (pcmat->apply) { 66345a4b08SToby Isaac case PCMATOP_MULT: 679566063dSJacob Faibussowitsch PetscCall(MatMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 68345a4b08SToby Isaac break; 69345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 70345a4b08SToby Isaac PetscCall(MatTransposeMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 71345a4b08SToby Isaac break; 72345a4b08SToby Isaac case PCMATOP_SOLVE: 73345a4b08SToby Isaac PetscCall(MatMatSolve(pc->pmat, X, Y)); 74345a4b08SToby Isaac break; 75345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 76345a4b08SToby Isaac PetscCall(MatMatSolveTranspose(pc->pmat, X, Y)); 77345a4b08SToby Isaac break; 78*b4f8a55aSStefano Zampini case PCMATOP_MULT_HERMITIAN_TRANSPOSE: 79345a4b08SToby Isaac PetscCall(MatDuplicate(X, MAT_COPY_VALUES, &W)); 80345a4b08SToby Isaac PetscCall(MatConjugate(W)); 81345a4b08SToby Isaac PetscCall(MatTransposeMatMult(pc->pmat, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 82345a4b08SToby Isaac PetscCall(MatConjugate(Y)); 83345a4b08SToby Isaac PetscCall(MatDestroy(&W)); 84345a4b08SToby Isaac break; 85*b4f8a55aSStefano Zampini default: 86*b4f8a55aSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported %s case", PCMatOpTypes[pcmat->apply]); 87345a4b08SToby Isaac } 883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 897b6e2003SPierre Jolivet } 907b6e2003SPierre Jolivet 91d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Mat(PC pc, Vec x, Vec y) 92d71ae5a4SJacob Faibussowitsch { 93345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 94*b4f8a55aSStefano Zampini Vec w; 95345a4b08SToby Isaac 964b9ad928SBarry Smith PetscFunctionBegin; 97345a4b08SToby Isaac switch (pcmat->apply) { 98345a4b08SToby Isaac case PCMATOP_MULT: 999566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->pmat, x, y)); 100345a4b08SToby Isaac break; 101345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 102345a4b08SToby Isaac PetscCall(MatMult(pc->pmat, x, y)); 103345a4b08SToby Isaac break; 104345a4b08SToby Isaac case PCMATOP_SOLVE: 105345a4b08SToby Isaac PetscCall(MatSolveTranspose(pc->pmat, x, y)); 106345a4b08SToby Isaac break; 107345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 108345a4b08SToby Isaac PetscCall(MatSolve(pc->pmat, x, y)); 109345a4b08SToby Isaac break; 110*b4f8a55aSStefano Zampini case PCMATOP_MULT_HERMITIAN_TRANSPOSE: 111345a4b08SToby Isaac PetscCall(VecDuplicate(x, &w)); 112345a4b08SToby Isaac PetscCall(VecCopy(x, w)); 113345a4b08SToby Isaac PetscCall(VecConjugate(w)); 114345a4b08SToby Isaac PetscCall(MatMult(pc->pmat, w, y)); 115345a4b08SToby Isaac PetscCall(VecConjugate(y)); 116345a4b08SToby Isaac PetscCall(VecDestroy(&w)); 117345a4b08SToby Isaac break; 118*b4f8a55aSStefano Zampini default: 119*b4f8a55aSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported %s case", PCMatOpTypes[pcmat->apply]); 120345a4b08SToby Isaac } 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1224b9ad928SBarry Smith } 1234b9ad928SBarry Smith 124d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Mat(PC pc) 125d71ae5a4SJacob Faibussowitsch { 1264b9ad928SBarry Smith PetscFunctionBegin; 127345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", NULL)); 128345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", NULL)); 129345a4b08SToby Isaac PetscCall(PetscFree(pc->data)); 130345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 131345a4b08SToby Isaac } 132345a4b08SToby Isaac 133345a4b08SToby Isaac /*@ 134345a4b08SToby Isaac PCMatSetApplyOperation - Set which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`. 135345a4b08SToby Isaac 136345a4b08SToby Isaac Logically collective 137345a4b08SToby Isaac 138345a4b08SToby Isaac Input Parameters: 139345a4b08SToby Isaac + pc - An instance of `PCMAT` 140345a4b08SToby Isaac - matop - The selected `MatOperation` 141345a4b08SToby Isaac 142345a4b08SToby Isaac Level: intermediate 143345a4b08SToby Isaac 144345a4b08SToby Isaac Note: 145345a4b08SToby Isaac If you have a matrix type that implements an exact inverse that isn't a factorization, 146345a4b08SToby Isaac you can use `PCMatSetApplyOperation(pc, MATOP_SOLVE)`. 147345a4b08SToby Isaac 148345a4b08SToby Isaac .seealso: `PCMAT`, `PCMatGetApplyOperation()`, `PCApply()`, `MatOperation` 149345a4b08SToby Isaac @*/ 150345a4b08SToby Isaac PetscErrorCode PCMatSetApplyOperation(PC pc, MatOperation matop) 151345a4b08SToby Isaac { 152345a4b08SToby Isaac PetscFunctionBegin; 153345a4b08SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 154345a4b08SToby Isaac PetscTryMethod((PetscObject)pc, "PCMatSetApplyOperation_C", (PC, MatOperation), (pc, matop)); 155345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 156345a4b08SToby Isaac } 157345a4b08SToby Isaac 158345a4b08SToby Isaac /*@ 159345a4b08SToby Isaac PCMatGetApplyOperation - Get which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`. 160345a4b08SToby Isaac 161345a4b08SToby Isaac Logically collective 162345a4b08SToby Isaac 163345a4b08SToby Isaac Input Parameter: 164345a4b08SToby Isaac . pc - An instance of `PCMAT` 165345a4b08SToby Isaac 166345a4b08SToby Isaac Output Parameter: 167345a4b08SToby Isaac . matop - The `MatOperation` 168345a4b08SToby Isaac 169345a4b08SToby Isaac Level: intermediate 170345a4b08SToby Isaac 171345a4b08SToby Isaac .seealso: `PCMAT`, `PCMatSetApplyOperation()`, `PCApply()`, `MatOperation` 172345a4b08SToby Isaac @*/ 173345a4b08SToby Isaac PetscErrorCode PCMatGetApplyOperation(PC pc, MatOperation *matop) 174345a4b08SToby Isaac { 175345a4b08SToby Isaac PetscFunctionBegin; 176345a4b08SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1774f572ea9SToby Isaac PetscAssertPointer(matop, 2); 178345a4b08SToby Isaac PetscUseMethod((PetscObject)pc, "PCMatGetApplyOperation_C", (PC, MatOperation *), (pc, matop)); 179345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 180345a4b08SToby Isaac } 181345a4b08SToby Isaac 182345a4b08SToby Isaac static PetscErrorCode PCMatSetApplyOperation_Mat(PC pc, MatOperation matop) 183345a4b08SToby Isaac { 184345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 185345a4b08SToby Isaac PCMatOperation pcmatop; 186345a4b08SToby Isaac 187345a4b08SToby Isaac PetscFunctionBegin; 188345a4b08SToby Isaac 189345a4b08SToby Isaac // clang-format off 190345a4b08SToby Isaac #define MATOP_TO_PCMATOP_CASE(var, OP) case MATOP_##OP: (var) = PCMATOP_##OP; break 191345a4b08SToby Isaac switch (matop) { 192345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT); 193345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT_TRANSPOSE); 194345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT_HERMITIAN_TRANSPOSE); 195345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE); 196345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE_TRANSPOSE); 197345a4b08SToby Isaac default: 198345a4b08SToby Isaac SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported MatOperation %d for PCMatSetApplyOperation()", (int)matop); 199345a4b08SToby Isaac } 200345a4b08SToby Isaac #undef MATOP_TO_PCMATOP_CASE 201345a4b08SToby Isaac // clang-format on 202345a4b08SToby Isaac 203345a4b08SToby Isaac pcmat->apply = pcmatop; 204345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 205345a4b08SToby Isaac } 206345a4b08SToby Isaac 207345a4b08SToby Isaac static PetscErrorCode PCMatGetApplyOperation_Mat(PC pc, MatOperation *matop_p) 208345a4b08SToby Isaac { 209345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 210345a4b08SToby Isaac MatOperation matop = MATOP_MULT; 211345a4b08SToby Isaac 212345a4b08SToby Isaac PetscFunctionBegin; 213*b4f8a55aSStefano Zampini if (!pc->setupcalled) PetscCall(PCSetUp(pc)); 214345a4b08SToby Isaac 215345a4b08SToby Isaac // clang-format off 216345a4b08SToby Isaac #define PCMATOP_TO_MATOP_CASE(var, OP) case PCMATOP_##OP: (var) = MATOP_##OP; break 217345a4b08SToby Isaac switch (pcmat->apply) { 218345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT); 219345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT_TRANSPOSE); 220345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT_HERMITIAN_TRANSPOSE); 221345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, SOLVE); 222345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, SOLVE_TRANSPOSE); 223*b4f8a55aSStefano Zampini default: 224*b4f8a55aSStefano Zampini SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported %s case", PCMatOpTypes[pcmat->apply]); 225345a4b08SToby Isaac } 226345a4b08SToby Isaac #undef PCMATOP_TO_MATOP_CASE 227345a4b08SToby Isaac // clang-format on 228345a4b08SToby Isaac 229345a4b08SToby Isaac *matop_p = matop; 230345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 231345a4b08SToby Isaac } 232345a4b08SToby Isaac 233345a4b08SToby Isaac static PetscErrorCode PCView_Mat(PC pc, PetscViewer viewer) 234345a4b08SToby Isaac { 235345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 236345a4b08SToby Isaac PetscBool iascii; 237345a4b08SToby Isaac 238345a4b08SToby Isaac PetscFunctionBegin; 239345a4b08SToby Isaac PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 240345a4b08SToby Isaac if (iascii) { PetscCall(PetscViewerASCIIPrintf(viewer, "PCApply() == Mat%s()\n", PCMatOpTypes[pcmat->apply])); } 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2424b9ad928SBarry Smith } 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith /*MC 245345a4b08SToby Isaac PCMAT - A preconditioner obtained by applying an operation of the `pmat` provided in 246345a4b08SToby Isaac in `PCSetOperators(pc,amat,pmat)` or `KSPSetOperators(ksp,amat,pmat)`. By default, the operation is `MATOP_MULT`, 247345a4b08SToby Isaac meaning that the `pmat` provides an approximate inverse of `amat`. 248345a4b08SToby Isaac If some other operation of `pmat` implements the approximate inverse, 249345a4b08SToby Isaac use `PCMatSetApplyOperation()` to select that operation. 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith Level: intermediate 2524b9ad928SBarry Smith 253345a4b08SToby Isaac .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSHELL`, `MatOperation`, `PCMatSetApplyOperation()`, `PCMatGetApplyOperation()` 2544b9ad928SBarry Smith M*/ 2554b9ad928SBarry Smith 256d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc) 257d71ae5a4SJacob Faibussowitsch { 258345a4b08SToby Isaac PC_Mat *data; 259345a4b08SToby Isaac 2604b9ad928SBarry Smith PetscFunctionBegin; 261345a4b08SToby Isaac PetscCall(PetscNew(&data)); 262345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", PCMatSetApplyOperation_Mat)); 263345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", PCMatGetApplyOperation_Mat)); 264345a4b08SToby Isaac pc->data = data; 2654b9ad928SBarry Smith pc->ops->apply = PCApply_Mat; 2667b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Mat; 2674b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Mat; 268*b4f8a55aSStefano Zampini pc->ops->setup = PCSetUp_Mat; 2694b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Mat; 2700a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 271345a4b08SToby Isaac pc->ops->view = PCView_Mat; 2720a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2730a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 2740a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 2753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2764b9ad928SBarry Smith } 277