14b9ad928SBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 34b9ad928SBarry Smith 4*345a4b08SToby Isaac typedef enum { 5*345a4b08SToby Isaac PCMATOP_MULT, 6*345a4b08SToby Isaac PCMATOP_MULT_TRANSPOSE, 7*345a4b08SToby Isaac PCMATOP_MULT_HERMITIAN_TRANSPOSE, 8*345a4b08SToby Isaac PCMATOP_SOLVE, 9*345a4b08SToby Isaac PCMATOP_SOLVE_TRANSPOSE, 10*345a4b08SToby Isaac } PCMatOperation; 11*345a4b08SToby Isaac 12*345a4b08SToby Isaac const char *const PCMatOpTypes[] = {"Mult", "MultTranspose", "MultHermitianTranspose", "Solve", "SolveTranspose", NULL}; 13*345a4b08SToby Isaac 14*345a4b08SToby Isaac typedef struct _PCMAT { 15*345a4b08SToby Isaac PCMatOperation apply; 16*345a4b08SToby Isaac } PC_Mat; 17*345a4b08SToby Isaac 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Mat(PC pc, Vec x, Vec y) 19d71ae5a4SJacob Faibussowitsch { 20*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 21*345a4b08SToby Isaac 224b9ad928SBarry Smith PetscFunctionBegin; 23*345a4b08SToby Isaac switch (pcmat->apply) { 24*345a4b08SToby Isaac case PCMATOP_MULT: 259566063dSJacob Faibussowitsch PetscCall(MatMult(pc->pmat, x, y)); 26*345a4b08SToby Isaac break; 27*345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 28*345a4b08SToby Isaac PetscCall(MatMultTranspose(pc->pmat, x, y)); 29*345a4b08SToby Isaac break; 30*345a4b08SToby Isaac case PCMATOP_SOLVE: 31*345a4b08SToby Isaac PetscCall(MatSolve(pc->pmat, x, y)); 32*345a4b08SToby Isaac break; 33*345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 34*345a4b08SToby Isaac PetscCall(MatSolveTranspose(pc->pmat, x, y)); 35*345a4b08SToby Isaac break; 36*345a4b08SToby Isaac case PCMATOP_MULT_HERMITIAN_TRANSPOSE: 37*345a4b08SToby Isaac PetscCall(MatMultHermitianTranspose(pc->pmat, x, y)); 38*345a4b08SToby Isaac break; 39*345a4b08SToby Isaac } 403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 414b9ad928SBarry Smith } 424b9ad928SBarry Smith 43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Mat(PC pc, Mat X, Mat Y) 44d71ae5a4SJacob Faibussowitsch { 45*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 46*345a4b08SToby Isaac 477b6e2003SPierre Jolivet PetscFunctionBegin; 48*345a4b08SToby Isaac switch (pcmat->apply) { 49*345a4b08SToby Isaac case PCMATOP_MULT: 509566063dSJacob Faibussowitsch PetscCall(MatMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 51*345a4b08SToby Isaac break; 52*345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 53*345a4b08SToby Isaac PetscCall(MatTransposeMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 54*345a4b08SToby Isaac break; 55*345a4b08SToby Isaac case PCMATOP_SOLVE: 56*345a4b08SToby Isaac PetscCall(MatMatSolve(pc->pmat, X, Y)); 57*345a4b08SToby Isaac break; 58*345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 59*345a4b08SToby Isaac PetscCall(MatMatSolveTranspose(pc->pmat, X, Y)); 60*345a4b08SToby Isaac break; 61*345a4b08SToby Isaac case PCMATOP_MULT_HERMITIAN_TRANSPOSE: { 62*345a4b08SToby Isaac Mat W; 63*345a4b08SToby Isaac 64*345a4b08SToby Isaac PetscCall(MatDuplicate(X, MAT_COPY_VALUES, &W)); 65*345a4b08SToby Isaac PetscCall(MatConjugate(W)); 66*345a4b08SToby Isaac PetscCall(MatTransposeMatMult(pc->pmat, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y)); 67*345a4b08SToby Isaac PetscCall(MatConjugate(Y)); 68*345a4b08SToby Isaac PetscCall(MatDestroy(&W)); 69*345a4b08SToby Isaac break; 70*345a4b08SToby Isaac } 71*345a4b08SToby Isaac } 723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 737b6e2003SPierre Jolivet } 747b6e2003SPierre Jolivet 75d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Mat(PC pc, Vec x, Vec y) 76d71ae5a4SJacob Faibussowitsch { 77*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 78*345a4b08SToby Isaac 794b9ad928SBarry Smith PetscFunctionBegin; 80*345a4b08SToby Isaac switch (pcmat->apply) { 81*345a4b08SToby Isaac case PCMATOP_MULT: 829566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(pc->pmat, x, y)); 83*345a4b08SToby Isaac break; 84*345a4b08SToby Isaac case PCMATOP_MULT_TRANSPOSE: 85*345a4b08SToby Isaac PetscCall(MatMult(pc->pmat, x, y)); 86*345a4b08SToby Isaac break; 87*345a4b08SToby Isaac case PCMATOP_SOLVE: 88*345a4b08SToby Isaac PetscCall(MatSolveTranspose(pc->pmat, x, y)); 89*345a4b08SToby Isaac break; 90*345a4b08SToby Isaac case PCMATOP_SOLVE_TRANSPOSE: 91*345a4b08SToby Isaac PetscCall(MatSolve(pc->pmat, x, y)); 92*345a4b08SToby Isaac break; 93*345a4b08SToby Isaac case PCMATOP_MULT_HERMITIAN_TRANSPOSE: { 94*345a4b08SToby Isaac Vec w; 95*345a4b08SToby Isaac 96*345a4b08SToby Isaac PetscCall(VecDuplicate(x, &w)); 97*345a4b08SToby Isaac PetscCall(VecCopy(x, w)); 98*345a4b08SToby Isaac PetscCall(VecConjugate(w)); 99*345a4b08SToby Isaac PetscCall(MatMult(pc->pmat, w, y)); 100*345a4b08SToby Isaac PetscCall(VecConjugate(y)); 101*345a4b08SToby Isaac PetscCall(VecDestroy(&w)); 102*345a4b08SToby Isaac break; 103*345a4b08SToby Isaac } 104*345a4b08SToby Isaac } 1053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1064b9ad928SBarry Smith } 1074b9ad928SBarry Smith 108d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Mat(PC pc) 109d71ae5a4SJacob Faibussowitsch { 1104b9ad928SBarry Smith PetscFunctionBegin; 111*345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", NULL)); 112*345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", NULL)); 113*345a4b08SToby Isaac PetscCall(PetscFree(pc->data)); 114*345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 115*345a4b08SToby Isaac } 116*345a4b08SToby Isaac 117*345a4b08SToby Isaac /*@ 118*345a4b08SToby Isaac PCMatSetApplyOperation - Set which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`. 119*345a4b08SToby Isaac 120*345a4b08SToby Isaac Logically collective 121*345a4b08SToby Isaac 122*345a4b08SToby Isaac Input Parameters: 123*345a4b08SToby Isaac + pc - An instance of `PCMAT` 124*345a4b08SToby Isaac - matop - The selected `MatOperation` 125*345a4b08SToby Isaac 126*345a4b08SToby Isaac Level: intermediate 127*345a4b08SToby Isaac 128*345a4b08SToby Isaac Note: 129*345a4b08SToby Isaac If you have a matrix type that implements an exact inverse that isn't a factorization, 130*345a4b08SToby Isaac you can use `PCMatSetApplyOperation(pc, MATOP_SOLVE)`. 131*345a4b08SToby Isaac 132*345a4b08SToby Isaac .seealso: `PCMAT`, `PCMatGetApplyOperation()`, `PCApply()`, `MatOperation` 133*345a4b08SToby Isaac @*/ 134*345a4b08SToby Isaac PetscErrorCode PCMatSetApplyOperation(PC pc, MatOperation matop) 135*345a4b08SToby Isaac { 136*345a4b08SToby Isaac PetscFunctionBegin; 137*345a4b08SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 138*345a4b08SToby Isaac PetscTryMethod((PetscObject)pc, "PCMatSetApplyOperation_C", (PC, MatOperation), (pc, matop)); 139*345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 140*345a4b08SToby Isaac } 141*345a4b08SToby Isaac 142*345a4b08SToby Isaac /*@ 143*345a4b08SToby Isaac PCMatGetApplyOperation - Get which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`. 144*345a4b08SToby Isaac 145*345a4b08SToby Isaac Logically collective 146*345a4b08SToby Isaac 147*345a4b08SToby Isaac Input Parameter: 148*345a4b08SToby Isaac . pc - An instance of `PCMAT` 149*345a4b08SToby Isaac 150*345a4b08SToby Isaac Output Parameter: 151*345a4b08SToby Isaac . matop - The `MatOperation` 152*345a4b08SToby Isaac 153*345a4b08SToby Isaac Level: intermediate 154*345a4b08SToby Isaac 155*345a4b08SToby Isaac .seealso: `PCMAT`, `PCMatSetApplyOperation()`, `PCApply()`, `MatOperation` 156*345a4b08SToby Isaac @*/ 157*345a4b08SToby Isaac PetscErrorCode PCMatGetApplyOperation(PC pc, MatOperation *matop) 158*345a4b08SToby Isaac { 159*345a4b08SToby Isaac PetscFunctionBegin; 160*345a4b08SToby Isaac PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 161*345a4b08SToby Isaac PetscValidPointer(matop, 2); 162*345a4b08SToby Isaac PetscUseMethod((PetscObject)pc, "PCMatGetApplyOperation_C", (PC, MatOperation *), (pc, matop)); 163*345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 164*345a4b08SToby Isaac } 165*345a4b08SToby Isaac 166*345a4b08SToby Isaac static PetscErrorCode PCMatSetApplyOperation_Mat(PC pc, MatOperation matop) 167*345a4b08SToby Isaac { 168*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 169*345a4b08SToby Isaac PCMatOperation pcmatop; 170*345a4b08SToby Isaac 171*345a4b08SToby Isaac PetscFunctionBegin; 172*345a4b08SToby Isaac 173*345a4b08SToby Isaac // clang-format off 174*345a4b08SToby Isaac #define MATOP_TO_PCMATOP_CASE(var, OP) case MATOP_##OP: (var) = PCMATOP_##OP; break 175*345a4b08SToby Isaac switch (matop) { 176*345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT); 177*345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT_TRANSPOSE); 178*345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, MULT_HERMITIAN_TRANSPOSE); 179*345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE); 180*345a4b08SToby Isaac MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE_TRANSPOSE); 181*345a4b08SToby Isaac default: 182*345a4b08SToby Isaac SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported MatOperation %d for PCMatSetApplyOperation()", (int)matop); 183*345a4b08SToby Isaac } 184*345a4b08SToby Isaac #undef MATOP_TO_PCMATOP_CASE 185*345a4b08SToby Isaac // clang-format on 186*345a4b08SToby Isaac 187*345a4b08SToby Isaac pcmat->apply = pcmatop; 188*345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 189*345a4b08SToby Isaac } 190*345a4b08SToby Isaac 191*345a4b08SToby Isaac static PetscErrorCode PCMatGetApplyOperation_Mat(PC pc, MatOperation *matop_p) 192*345a4b08SToby Isaac { 193*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 194*345a4b08SToby Isaac MatOperation matop = MATOP_MULT; 195*345a4b08SToby Isaac 196*345a4b08SToby Isaac PetscFunctionBegin; 197*345a4b08SToby Isaac 198*345a4b08SToby Isaac // clang-format off 199*345a4b08SToby Isaac #define PCMATOP_TO_MATOP_CASE(var, OP) case PCMATOP_##OP: (var) = MATOP_##OP; break 200*345a4b08SToby Isaac switch (pcmat->apply) { 201*345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT); 202*345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT_TRANSPOSE); 203*345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, MULT_HERMITIAN_TRANSPOSE); 204*345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, SOLVE); 205*345a4b08SToby Isaac PCMATOP_TO_MATOP_CASE(matop, SOLVE_TRANSPOSE); 206*345a4b08SToby Isaac } 207*345a4b08SToby Isaac #undef PCMATOP_TO_MATOP_CASE 208*345a4b08SToby Isaac // clang-format on 209*345a4b08SToby Isaac 210*345a4b08SToby Isaac *matop_p = matop; 211*345a4b08SToby Isaac PetscFunctionReturn(PETSC_SUCCESS); 212*345a4b08SToby Isaac } 213*345a4b08SToby Isaac 214*345a4b08SToby Isaac static PetscErrorCode PCView_Mat(PC pc, PetscViewer viewer) 215*345a4b08SToby Isaac { 216*345a4b08SToby Isaac PC_Mat *pcmat = (PC_Mat *)pc->data; 217*345a4b08SToby Isaac PetscBool iascii; 218*345a4b08SToby Isaac 219*345a4b08SToby Isaac PetscFunctionBegin; 220*345a4b08SToby Isaac PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 221*345a4b08SToby Isaac if (iascii) { PetscCall(PetscViewerASCIIPrintf(viewer, "PCApply() == Mat%s()\n", PCMatOpTypes[pcmat->apply])); } 2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2234b9ad928SBarry Smith } 2244b9ad928SBarry Smith 2254b9ad928SBarry Smith /*MC 226*345a4b08SToby Isaac PCMAT - A preconditioner obtained by applying an operation of the `pmat` provided in 227*345a4b08SToby Isaac in `PCSetOperators(pc,amat,pmat)` or `KSPSetOperators(ksp,amat,pmat)`. By default, the operation is `MATOP_MULT`, 228*345a4b08SToby Isaac meaning that the `pmat` provides an approximate inverse of `amat`. 229*345a4b08SToby Isaac If some other operation of `pmat` implements the approximate inverse, 230*345a4b08SToby Isaac use `PCMatSetApplyOperation()` to select that operation. 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith Level: intermediate 2334b9ad928SBarry Smith 234*345a4b08SToby Isaac .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSHELL`, `MatOperation`, `PCMatSetApplyOperation()`, `PCMatGetApplyOperation()` 2354b9ad928SBarry Smith M*/ 2364b9ad928SBarry Smith 237d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc) 238d71ae5a4SJacob Faibussowitsch { 239*345a4b08SToby Isaac PC_Mat *data; 240*345a4b08SToby Isaac 2414b9ad928SBarry Smith PetscFunctionBegin; 242*345a4b08SToby Isaac PetscCall(PetscNew(&data)); 243*345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", PCMatSetApplyOperation_Mat)); 244*345a4b08SToby Isaac PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", PCMatGetApplyOperation_Mat)); 245*345a4b08SToby Isaac data->apply = PCMATOP_MULT; 246*345a4b08SToby Isaac pc->data = data; 2474b9ad928SBarry Smith pc->ops->apply = PCApply_Mat; 2487b6e2003SPierre Jolivet pc->ops->matapply = PCMatApply_Mat; 2494b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Mat; 2500a545947SLisandro Dalcin pc->ops->setup = NULL; 2514b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Mat; 2520a545947SLisandro Dalcin pc->ops->setfromoptions = NULL; 253*345a4b08SToby Isaac pc->ops->view = PCView_Mat; 2540a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 2550a545947SLisandro Dalcin pc->ops->applysymmetricleft = NULL; 2560a545947SLisandro Dalcin pc->ops->applysymmetricright = NULL; 2573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2584b9ad928SBarry Smith } 259