xref: /petsc/src/ksp/pc/impls/mat/pcmat.c (revision 562efe2ef922487c6beae96ba39e19afd4eefbe6)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
24b9ad928SBarry Smith 
3345a4b08SToby Isaac typedef enum {
4345a4b08SToby Isaac   PCMATOP_MULT,
5345a4b08SToby Isaac   PCMATOP_MULT_TRANSPOSE,
6345a4b08SToby Isaac   PCMATOP_MULT_HERMITIAN_TRANSPOSE,
7345a4b08SToby Isaac   PCMATOP_SOLVE,
8345a4b08SToby Isaac   PCMATOP_SOLVE_TRANSPOSE,
9345a4b08SToby Isaac } PCMatOperation;
10345a4b08SToby Isaac 
11345a4b08SToby Isaac const char *const PCMatOpTypes[] = {"Mult", "MultTranspose", "MultHermitianTranspose", "Solve", "SolveTranspose", NULL};
12345a4b08SToby Isaac 
13345a4b08SToby Isaac typedef struct _PCMAT {
14345a4b08SToby Isaac   PCMatOperation apply;
15345a4b08SToby Isaac } PC_Mat;
16345a4b08SToby Isaac 
17d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Mat(PC pc, Vec x, Vec y)
18d71ae5a4SJacob Faibussowitsch {
19345a4b08SToby Isaac   PC_Mat *pcmat = (PC_Mat *)pc->data;
20345a4b08SToby Isaac 
214b9ad928SBarry Smith   PetscFunctionBegin;
22345a4b08SToby Isaac   switch (pcmat->apply) {
23345a4b08SToby Isaac   case PCMATOP_MULT:
249566063dSJacob Faibussowitsch     PetscCall(MatMult(pc->pmat, x, y));
25345a4b08SToby Isaac     break;
26345a4b08SToby Isaac   case PCMATOP_MULT_TRANSPOSE:
27345a4b08SToby Isaac     PetscCall(MatMultTranspose(pc->pmat, x, y));
28345a4b08SToby Isaac     break;
29345a4b08SToby Isaac   case PCMATOP_SOLVE:
30345a4b08SToby Isaac     PetscCall(MatSolve(pc->pmat, x, y));
31345a4b08SToby Isaac     break;
32345a4b08SToby Isaac   case PCMATOP_SOLVE_TRANSPOSE:
33345a4b08SToby Isaac     PetscCall(MatSolveTranspose(pc->pmat, x, y));
34345a4b08SToby Isaac     break;
35345a4b08SToby Isaac   case PCMATOP_MULT_HERMITIAN_TRANSPOSE:
36345a4b08SToby Isaac     PetscCall(MatMultHermitianTranspose(pc->pmat, x, y));
37345a4b08SToby Isaac     break;
38345a4b08SToby Isaac   }
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
404b9ad928SBarry Smith }
414b9ad928SBarry Smith 
42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_Mat(PC pc, Mat X, Mat Y)
43d71ae5a4SJacob Faibussowitsch {
44345a4b08SToby Isaac   PC_Mat *pcmat = (PC_Mat *)pc->data;
45345a4b08SToby Isaac 
467b6e2003SPierre Jolivet   PetscFunctionBegin;
47345a4b08SToby Isaac   switch (pcmat->apply) {
48345a4b08SToby Isaac   case PCMATOP_MULT:
499566063dSJacob Faibussowitsch     PetscCall(MatMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
50345a4b08SToby Isaac     break;
51345a4b08SToby Isaac   case PCMATOP_MULT_TRANSPOSE:
52345a4b08SToby Isaac     PetscCall(MatTransposeMatMult(pc->pmat, X, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
53345a4b08SToby Isaac     break;
54345a4b08SToby Isaac   case PCMATOP_SOLVE:
55345a4b08SToby Isaac     PetscCall(MatMatSolve(pc->pmat, X, Y));
56345a4b08SToby Isaac     break;
57345a4b08SToby Isaac   case PCMATOP_SOLVE_TRANSPOSE:
58345a4b08SToby Isaac     PetscCall(MatMatSolveTranspose(pc->pmat, X, Y));
59345a4b08SToby Isaac     break;
60345a4b08SToby Isaac   case PCMATOP_MULT_HERMITIAN_TRANSPOSE: {
61345a4b08SToby Isaac     Mat W;
62345a4b08SToby Isaac 
63345a4b08SToby Isaac     PetscCall(MatDuplicate(X, MAT_COPY_VALUES, &W));
64345a4b08SToby Isaac     PetscCall(MatConjugate(W));
65345a4b08SToby Isaac     PetscCall(MatTransposeMatMult(pc->pmat, W, MAT_REUSE_MATRIX, PETSC_DEFAULT, &Y));
66345a4b08SToby Isaac     PetscCall(MatConjugate(Y));
67345a4b08SToby Isaac     PetscCall(MatDestroy(&W));
68345a4b08SToby Isaac     break;
69345a4b08SToby Isaac   }
70345a4b08SToby Isaac   }
713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
727b6e2003SPierre Jolivet }
737b6e2003SPierre Jolivet 
74d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Mat(PC pc, Vec x, Vec y)
75d71ae5a4SJacob Faibussowitsch {
76345a4b08SToby Isaac   PC_Mat *pcmat = (PC_Mat *)pc->data;
77345a4b08SToby Isaac 
784b9ad928SBarry Smith   PetscFunctionBegin;
79345a4b08SToby Isaac   switch (pcmat->apply) {
80345a4b08SToby Isaac   case PCMATOP_MULT:
819566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(pc->pmat, x, y));
82345a4b08SToby Isaac     break;
83345a4b08SToby Isaac   case PCMATOP_MULT_TRANSPOSE:
84345a4b08SToby Isaac     PetscCall(MatMult(pc->pmat, x, y));
85345a4b08SToby Isaac     break;
86345a4b08SToby Isaac   case PCMATOP_SOLVE:
87345a4b08SToby Isaac     PetscCall(MatSolveTranspose(pc->pmat, x, y));
88345a4b08SToby Isaac     break;
89345a4b08SToby Isaac   case PCMATOP_SOLVE_TRANSPOSE:
90345a4b08SToby Isaac     PetscCall(MatSolve(pc->pmat, x, y));
91345a4b08SToby Isaac     break;
92345a4b08SToby Isaac   case PCMATOP_MULT_HERMITIAN_TRANSPOSE: {
93345a4b08SToby Isaac     Vec w;
94345a4b08SToby Isaac 
95345a4b08SToby Isaac     PetscCall(VecDuplicate(x, &w));
96345a4b08SToby Isaac     PetscCall(VecCopy(x, w));
97345a4b08SToby Isaac     PetscCall(VecConjugate(w));
98345a4b08SToby Isaac     PetscCall(MatMult(pc->pmat, w, y));
99345a4b08SToby Isaac     PetscCall(VecConjugate(y));
100345a4b08SToby Isaac     PetscCall(VecDestroy(&w));
101345a4b08SToby Isaac     break;
102345a4b08SToby Isaac   }
103345a4b08SToby Isaac   }
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1054b9ad928SBarry Smith }
1064b9ad928SBarry Smith 
107d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Mat(PC pc)
108d71ae5a4SJacob Faibussowitsch {
1094b9ad928SBarry Smith   PetscFunctionBegin;
110345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", NULL));
111345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", NULL));
112345a4b08SToby Isaac   PetscCall(PetscFree(pc->data));
113345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
114345a4b08SToby Isaac }
115345a4b08SToby Isaac 
116345a4b08SToby Isaac /*@
117345a4b08SToby Isaac   PCMatSetApplyOperation - Set which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`.
118345a4b08SToby Isaac 
119345a4b08SToby Isaac   Logically collective
120345a4b08SToby Isaac 
121345a4b08SToby Isaac   Input Parameters:
122345a4b08SToby Isaac + pc    - An instance of `PCMAT`
123345a4b08SToby Isaac - matop - The selected `MatOperation`
124345a4b08SToby Isaac 
125345a4b08SToby Isaac   Level: intermediate
126345a4b08SToby Isaac 
127345a4b08SToby Isaac   Note:
128345a4b08SToby Isaac   If you have a matrix type that implements an exact inverse that isn't a factorization,
129345a4b08SToby Isaac   you can use `PCMatSetApplyOperation(pc, MATOP_SOLVE)`.
130345a4b08SToby Isaac 
131*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMAT`, `PCMatGetApplyOperation()`, `PCApply()`, `MatOperation`
132345a4b08SToby Isaac @*/
133345a4b08SToby Isaac PetscErrorCode PCMatSetApplyOperation(PC pc, MatOperation matop)
134345a4b08SToby Isaac {
135345a4b08SToby Isaac   PetscFunctionBegin;
136345a4b08SToby Isaac   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
137345a4b08SToby Isaac   PetscTryMethod((PetscObject)pc, "PCMatSetApplyOperation_C", (PC, MatOperation), (pc, matop));
138345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
139345a4b08SToby Isaac }
140345a4b08SToby Isaac 
141345a4b08SToby Isaac /*@
142345a4b08SToby Isaac   PCMatGetApplyOperation - Get which matrix operation of the preconditioning matrix implements `PCApply()` for `PCMAT`.
143345a4b08SToby Isaac 
144345a4b08SToby Isaac   Logically collective
145345a4b08SToby Isaac 
146345a4b08SToby Isaac   Input Parameter:
147345a4b08SToby Isaac . pc - An instance of `PCMAT`
148345a4b08SToby Isaac 
149345a4b08SToby Isaac   Output Parameter:
150345a4b08SToby Isaac . matop - The `MatOperation`
151345a4b08SToby Isaac 
152345a4b08SToby Isaac   Level: intermediate
153345a4b08SToby Isaac 
154*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCMAT`, `PCMatSetApplyOperation()`, `PCApply()`, `MatOperation`
155345a4b08SToby Isaac @*/
156345a4b08SToby Isaac PetscErrorCode PCMatGetApplyOperation(PC pc, MatOperation *matop)
157345a4b08SToby Isaac {
158345a4b08SToby Isaac   PetscFunctionBegin;
159345a4b08SToby Isaac   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1604f572ea9SToby Isaac   PetscAssertPointer(matop, 2);
161345a4b08SToby Isaac   PetscUseMethod((PetscObject)pc, "PCMatGetApplyOperation_C", (PC, MatOperation *), (pc, matop));
162345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
163345a4b08SToby Isaac }
164345a4b08SToby Isaac 
165345a4b08SToby Isaac static PetscErrorCode PCMatSetApplyOperation_Mat(PC pc, MatOperation matop)
166345a4b08SToby Isaac {
167345a4b08SToby Isaac   PC_Mat        *pcmat = (PC_Mat *)pc->data;
168345a4b08SToby Isaac   PCMatOperation pcmatop;
169345a4b08SToby Isaac 
170345a4b08SToby Isaac   PetscFunctionBegin;
171345a4b08SToby Isaac 
172345a4b08SToby Isaac   // clang-format off
173345a4b08SToby Isaac #define MATOP_TO_PCMATOP_CASE(var, OP) case MATOP_##OP: (var) = PCMATOP_##OP; break
174345a4b08SToby Isaac   switch (matop) {
175345a4b08SToby Isaac   MATOP_TO_PCMATOP_CASE(pcmatop, MULT);
176345a4b08SToby Isaac   MATOP_TO_PCMATOP_CASE(pcmatop, MULT_TRANSPOSE);
177345a4b08SToby Isaac   MATOP_TO_PCMATOP_CASE(pcmatop, MULT_HERMITIAN_TRANSPOSE);
178345a4b08SToby Isaac   MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE);
179345a4b08SToby Isaac   MATOP_TO_PCMATOP_CASE(pcmatop, SOLVE_TRANSPOSE);
180345a4b08SToby Isaac   default:
181345a4b08SToby Isaac     SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Unsupported MatOperation %d for PCMatSetApplyOperation()", (int)matop);
182345a4b08SToby Isaac   }
183345a4b08SToby Isaac #undef MATOP_TO_PCMATOP_CASE
184345a4b08SToby Isaac   // clang-format on
185345a4b08SToby Isaac 
186345a4b08SToby Isaac   pcmat->apply = pcmatop;
187345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
188345a4b08SToby Isaac }
189345a4b08SToby Isaac 
190345a4b08SToby Isaac static PetscErrorCode PCMatGetApplyOperation_Mat(PC pc, MatOperation *matop_p)
191345a4b08SToby Isaac {
192345a4b08SToby Isaac   PC_Mat      *pcmat = (PC_Mat *)pc->data;
193345a4b08SToby Isaac   MatOperation matop = MATOP_MULT;
194345a4b08SToby Isaac 
195345a4b08SToby Isaac   PetscFunctionBegin;
196345a4b08SToby Isaac 
197345a4b08SToby Isaac   // clang-format off
198345a4b08SToby Isaac #define PCMATOP_TO_MATOP_CASE(var, OP) case PCMATOP_##OP: (var) = MATOP_##OP; break
199345a4b08SToby Isaac   switch (pcmat->apply) {
200345a4b08SToby Isaac   PCMATOP_TO_MATOP_CASE(matop, MULT);
201345a4b08SToby Isaac   PCMATOP_TO_MATOP_CASE(matop, MULT_TRANSPOSE);
202345a4b08SToby Isaac   PCMATOP_TO_MATOP_CASE(matop, MULT_HERMITIAN_TRANSPOSE);
203345a4b08SToby Isaac   PCMATOP_TO_MATOP_CASE(matop, SOLVE);
204345a4b08SToby Isaac   PCMATOP_TO_MATOP_CASE(matop, SOLVE_TRANSPOSE);
205345a4b08SToby Isaac   }
206345a4b08SToby Isaac #undef PCMATOP_TO_MATOP_CASE
207345a4b08SToby Isaac   // clang-format on
208345a4b08SToby Isaac 
209345a4b08SToby Isaac   *matop_p = matop;
210345a4b08SToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
211345a4b08SToby Isaac }
212345a4b08SToby Isaac 
213345a4b08SToby Isaac static PetscErrorCode PCView_Mat(PC pc, PetscViewer viewer)
214345a4b08SToby Isaac {
215345a4b08SToby Isaac   PC_Mat   *pcmat = (PC_Mat *)pc->data;
216345a4b08SToby Isaac   PetscBool iascii;
217345a4b08SToby Isaac 
218345a4b08SToby Isaac   PetscFunctionBegin;
219345a4b08SToby Isaac   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
220345a4b08SToby Isaac   if (iascii) { PetscCall(PetscViewerASCIIPrintf(viewer, "PCApply() == Mat%s()\n", PCMatOpTypes[pcmat->apply])); }
2213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2224b9ad928SBarry Smith }
2234b9ad928SBarry Smith 
2244b9ad928SBarry Smith /*MC
225345a4b08SToby Isaac      PCMAT - A preconditioner obtained by applying an operation of the `pmat` provided in
226345a4b08SToby Isaac              in `PCSetOperators(pc,amat,pmat)` or `KSPSetOperators(ksp,amat,pmat)`.  By default, the operation is `MATOP_MULT`,
227345a4b08SToby Isaac              meaning that the `pmat` provides an approximate inverse of `amat`.
228345a4b08SToby Isaac              If some other operation of `pmat` implements the approximate inverse,
229345a4b08SToby Isaac              use `PCMatSetApplyOperation()` to select that operation.
2304b9ad928SBarry Smith 
2314b9ad928SBarry Smith    Level: intermediate
2324b9ad928SBarry Smith 
233*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCSHELL`, `MatOperation`, `PCMatSetApplyOperation()`, `PCMatGetApplyOperation()`
2344b9ad928SBarry Smith M*/
2354b9ad928SBarry Smith 
236d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Mat(PC pc)
237d71ae5a4SJacob Faibussowitsch {
238345a4b08SToby Isaac   PC_Mat *data;
239345a4b08SToby Isaac 
2404b9ad928SBarry Smith   PetscFunctionBegin;
241345a4b08SToby Isaac   PetscCall(PetscNew(&data));
242345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatSetApplyOperation_C", PCMatSetApplyOperation_Mat));
243345a4b08SToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCMatGetApplyOperation_C", PCMatGetApplyOperation_Mat));
244345a4b08SToby Isaac   data->apply                  = PCMATOP_MULT;
245345a4b08SToby Isaac   pc->data                     = data;
2464b9ad928SBarry Smith   pc->ops->apply               = PCApply_Mat;
2477b6e2003SPierre Jolivet   pc->ops->matapply            = PCMatApply_Mat;
2484b9ad928SBarry Smith   pc->ops->applytranspose      = PCApplyTranspose_Mat;
2490a545947SLisandro Dalcin   pc->ops->setup               = NULL;
2504b9ad928SBarry Smith   pc->ops->destroy             = PCDestroy_Mat;
2510a545947SLisandro Dalcin   pc->ops->setfromoptions      = NULL;
252345a4b08SToby Isaac   pc->ops->view                = PCView_Mat;
2530a545947SLisandro Dalcin   pc->ops->applyrichardson     = NULL;
2540a545947SLisandro Dalcin   pc->ops->applysymmetricleft  = NULL;
2550a545947SLisandro Dalcin   pc->ops->applysymmetricright = NULL;
2563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2574b9ad928SBarry Smith }
258