xref: /petsc/src/ksp/pc/impls/galerkin/galerkin.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
12a6744ebSBarry Smith /*
22a6744ebSBarry Smith       Defines a preconditioner defined by R^T S R
32a6744ebSBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcimpl.h>
5c6db04a5SJed Brown #include <petscksp.h> /*I "petscksp.h" I*/
62a6744ebSBarry Smith 
72a6744ebSBarry Smith typedef struct {
82a6744ebSBarry Smith   KSP ksp;
92a6744ebSBarry Smith   Mat R, P;
102a6744ebSBarry Smith   Vec b, x;
11b3402f20SBarry Smith   PetscErrorCode (*computeasub)(PC, Mat, Mat, Mat *, void *);
12b3402f20SBarry Smith   void *computeasub_ctx;
132a6744ebSBarry Smith } PC_Galerkin;
142a6744ebSBarry Smith 
PCApply_Galerkin(PC pc,Vec x,Vec y)15d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Galerkin(PC pc, Vec x, Vec y)
16d71ae5a4SJacob Faibussowitsch {
172a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
182a6744ebSBarry Smith 
192a6744ebSBarry Smith   PetscFunctionBegin;
202fa5cd67SKarl Rupp   if (jac->R) {
219566063dSJacob Faibussowitsch     PetscCall(MatRestrict(jac->R, x, jac->b));
222fa5cd67SKarl Rupp   } else {
239566063dSJacob Faibussowitsch     PetscCall(MatRestrict(jac->P, x, jac->b));
242fa5cd67SKarl Rupp   }
259566063dSJacob Faibussowitsch   PetscCall(KSPSolve(jac->ksp, jac->b, jac->x));
269566063dSJacob Faibussowitsch   PetscCall(KSPCheckSolve(jac->ksp, pc, jac->x));
272fa5cd67SKarl Rupp   if (jac->P) {
289566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(jac->P, jac->x, y));
292fa5cd67SKarl Rupp   } else {
309566063dSJacob Faibussowitsch     PetscCall(MatInterpolate(jac->R, jac->x, y));
312fa5cd67SKarl Rupp   }
323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
332a6744ebSBarry Smith }
342a6744ebSBarry Smith 
PCSetUp_Galerkin(PC pc)35d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Galerkin(PC pc)
36d71ae5a4SJacob Faibussowitsch {
372a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
38ace3abfcSBarry Smith   PetscBool    a;
39906ed7ccSBarry Smith   Vec         *xx, *yy;
402a6744ebSBarry Smith 
412a6744ebSBarry Smith   PetscFunctionBegin;
42b3402f20SBarry Smith   if (jac->computeasub) {
43b3402f20SBarry Smith     Mat Ap;
44b3402f20SBarry Smith     if (!pc->setupcalled) {
459566063dSJacob Faibussowitsch       PetscCall((*jac->computeasub)(pc, pc->pmat, NULL, &Ap, jac->computeasub_ctx));
469566063dSJacob Faibussowitsch       PetscCall(KSPSetOperators(jac->ksp, Ap, Ap));
479566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&Ap));
48b3402f20SBarry Smith     } else {
499566063dSJacob Faibussowitsch       PetscCall(KSPGetOperators(jac->ksp, NULL, &Ap));
509566063dSJacob Faibussowitsch       PetscCall((*jac->computeasub)(pc, pc->pmat, Ap, NULL, jac->computeasub_ctx));
51b3402f20SBarry Smith     }
52b3402f20SBarry Smith   }
53b3402f20SBarry Smith 
542a6744ebSBarry Smith   if (!jac->x) {
559566063dSJacob Faibussowitsch     PetscCall(KSPGetOperatorsSet(jac->ksp, &a, NULL));
5628b400f6SJacob Faibussowitsch     PetscCheck(a, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()");
579566063dSJacob Faibussowitsch     PetscCall(KSPCreateVecs(jac->ksp, 1, &xx, 1, &yy));
58906ed7ccSBarry Smith     jac->x = *xx;
59906ed7ccSBarry Smith     jac->b = *yy;
609566063dSJacob Faibussowitsch     PetscCall(PetscFree(xx));
619566063dSJacob Faibussowitsch     PetscCall(PetscFree(yy));
622a6744ebSBarry Smith   }
637827d75bSBarry Smith   PetscCheck(jac->R || jac->P, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction()/Interpolation()");
642a6744ebSBarry Smith   /* should check here that sizes of R/P match size of a */
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
662a6744ebSBarry Smith }
672a6744ebSBarry Smith 
PCReset_Galerkin(PC pc)68d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Galerkin(PC pc)
69d71ae5a4SJacob Faibussowitsch {
702a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
712a6744ebSBarry Smith 
722a6744ebSBarry Smith   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->R));
749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->P));
759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->x));
769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->b));
779566063dSJacob Faibussowitsch   PetscCall(KSPReset(jac->ksp));
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
79a06653b4SBarry Smith }
80a06653b4SBarry Smith 
PCDestroy_Galerkin(PC pc)81d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Galerkin(PC pc)
82d71ae5a4SJacob Faibussowitsch {
83a06653b4SBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
84a06653b4SBarry Smith 
85a06653b4SBarry Smith   PetscFunctionBegin;
869566063dSJacob Faibussowitsch   PetscCall(PCReset_Galerkin(pc));
879566063dSJacob Faibussowitsch   PetscCall(KSPDestroy(&jac->ksp));
889566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
902a6744ebSBarry Smith }
912a6744ebSBarry Smith 
PCView_Galerkin(PC pc,PetscViewer viewer)92d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Galerkin(PC pc, PetscViewer viewer)
93d71ae5a4SJacob Faibussowitsch {
942a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
959f196a02SMartin Diehl   PetscBool    isascii;
962a6744ebSBarry Smith 
972a6744ebSBarry Smith   PetscFunctionBegin;
989f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
999f196a02SMartin Diehl   if (isascii) {
1009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  KSP on Galerkin follow\n"));
1019566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  ---------------------------------\n"));
1022a6744ebSBarry Smith   }
1039566063dSJacob Faibussowitsch   PetscCall(KSPView(jac->ksp, viewer));
1043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1052a6744ebSBarry Smith }
1062a6744ebSBarry Smith 
PCGalerkinGetKSP_Galerkin(PC pc,KSP * ksp)107d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGalerkinGetKSP_Galerkin(PC pc, KSP *ksp)
108d71ae5a4SJacob Faibussowitsch {
1092a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
1102a6744ebSBarry Smith 
1112a6744ebSBarry Smith   PetscFunctionBegin;
1122a6744ebSBarry Smith   *ksp = jac->ksp;
1133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1142a6744ebSBarry Smith }
1152a6744ebSBarry Smith 
PCGalerkinSetRestriction_Galerkin(PC pc,Mat R)116d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGalerkinSetRestriction_Galerkin(PC pc, Mat R)
117d71ae5a4SJacob Faibussowitsch {
1182a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
1192a6744ebSBarry Smith 
1202a6744ebSBarry Smith   PetscFunctionBegin;
1219566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)R));
1229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->R));
1232a6744ebSBarry Smith   jac->R = R;
1243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1252a6744ebSBarry Smith }
1262a6744ebSBarry Smith 
PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P)127d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCGalerkinSetInterpolation_Galerkin(PC pc, Mat P)
128d71ae5a4SJacob Faibussowitsch {
1292a6744ebSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
1302a6744ebSBarry Smith 
1312a6744ebSBarry Smith   PetscFunctionBegin;
1329566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)P));
1339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&jac->P));
1342a6744ebSBarry Smith   jac->P = P;
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1362a6744ebSBarry Smith }
1372a6744ebSBarry Smith 
PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (* computeAsub)(PC,Mat,Mat,Mat *,void *),PetscCtx ctx)138*2a8381b2SBarry Smith static PetscErrorCode PCGalerkinSetComputeSubmatrix_Galerkin(PC pc, PetscErrorCode (*computeAsub)(PC, Mat, Mat, Mat *, void *), PetscCtx ctx)
139d71ae5a4SJacob Faibussowitsch {
140b3402f20SBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
141b3402f20SBarry Smith 
142b3402f20SBarry Smith   PetscFunctionBegin;
143b3402f20SBarry Smith   jac->computeasub     = computeAsub;
144b3402f20SBarry Smith   jac->computeasub_ctx = ctx;
1453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146b3402f20SBarry Smith }
147b3402f20SBarry Smith 
1482a6744ebSBarry Smith /*@
149f1580f4eSBarry Smith   PCGalerkinSetRestriction - Sets the restriction operator for the `PCGALERKIN` preconditioner
1502a6744ebSBarry Smith 
151c3339decSBarry Smith   Logically Collective
1522a6744ebSBarry Smith 
153d8d19677SJose E. Roman   Input Parameters:
1542a6744ebSBarry Smith + pc - the preconditioner context
1552a6744ebSBarry Smith - R  - the restriction operator
1562a6744ebSBarry Smith 
157feefa0e1SJacob Faibussowitsch   Level: intermediate
15820f4b53cSBarry Smith 
159f1580f4eSBarry Smith   Note:
160f1580f4eSBarry Smith   Either this or `PCGalerkinSetInterpolation()` or both must be called
1612a6744ebSBarry Smith 
162562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetType()`, `PCType`, `PCGALERKIN`,
163db781477SPatrick Sanan           `PCGalerkinSetInterpolation()`, `PCGalerkinGetKSP()`
1642a6744ebSBarry Smith @*/
PCGalerkinSetRestriction(PC pc,Mat R)165d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGalerkinSetRestriction(PC pc, Mat R)
166d71ae5a4SJacob Faibussowitsch {
1672a6744ebSBarry Smith   PetscFunctionBegin;
1680700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
169cac4c232SBarry Smith   PetscTryMethod(pc, "PCGalerkinSetRestriction_C", (PC, Mat), (pc, R));
1703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1712a6744ebSBarry Smith }
1722a6744ebSBarry Smith 
1732a6744ebSBarry Smith /*@
174f1580f4eSBarry Smith   PCGalerkinSetInterpolation - Sets the interpolation operator for the `PCGALERKIN` preconditioner
1752a6744ebSBarry Smith 
176c3339decSBarry Smith   Logically Collective
1772a6744ebSBarry Smith 
178d8d19677SJose E. Roman   Input Parameters:
1792a6744ebSBarry Smith + pc - the preconditioner context
180feefa0e1SJacob Faibussowitsch - P  - the interpolation operator
1812a6744ebSBarry Smith 
182feefa0e1SJacob Faibussowitsch   Level: intermediate
18320f4b53cSBarry Smith 
184f1580f4eSBarry Smith   Note:
185f1580f4eSBarry Smith   Either this or `PCGalerkinSetRestriction()` or both must be called
1862a6744ebSBarry Smith 
187562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetType()`, `PCType`, `PCGALERKIN`,
188db781477SPatrick Sanan           `PCGalerkinSetRestriction()`, `PCGalerkinGetKSP()`
1892a6744ebSBarry Smith @*/
PCGalerkinSetInterpolation(PC pc,Mat P)190d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGalerkinSetInterpolation(PC pc, Mat P)
191d71ae5a4SJacob Faibussowitsch {
1922a6744ebSBarry Smith   PetscFunctionBegin;
1930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
194cac4c232SBarry Smith   PetscTryMethod(pc, "PCGalerkinSetInterpolation_C", (PC, Mat), (pc, P));
1953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1962a6744ebSBarry Smith }
1972a6744ebSBarry Smith 
198feefa0e1SJacob Faibussowitsch /*@C
199b3402f20SBarry Smith   PCGalerkinSetComputeSubmatrix - Provide a routine that will be called to compute the Galerkin submatrix
200b3402f20SBarry Smith 
201b3402f20SBarry Smith   Logically Collective
202b3402f20SBarry Smith 
203d8d19677SJose E. Roman   Input Parameters:
204b3402f20SBarry Smith + pc          - the preconditioner context
205b3402f20SBarry Smith . computeAsub - routine that computes the submatrix from the global matrix
20620f4b53cSBarry Smith - ctx         - context used by the routine, or `NULL`
207b3402f20SBarry Smith 
20820f4b53cSBarry Smith   Calling sequence of `computeAsub`:
20904c3f3b8SBarry Smith + pc  - the `PCGALERKIN` preconditioner
210f1580f4eSBarry Smith . A   - the matrix in the `PCGALERKIN`
21120f4b53cSBarry Smith . Ap  - the computed submatrix from any previous computation, if `NULL` it has not previously been computed
212b3402f20SBarry Smith . cAp - the submatrix computed by this routine
213b3402f20SBarry Smith - ctx - optional user-defined function context
214b3402f20SBarry Smith 
215feefa0e1SJacob Faibussowitsch   Level: intermediate
216b3402f20SBarry Smith 
21795452b02SPatrick Sanan   Notes:
218f1580f4eSBarry Smith   Instead of providing this routine you can call `PCGalerkinGetKSP()` and then `KSPSetOperators()` to provide the submatrix,
219f1580f4eSBarry Smith   but that will not work for multiple `KSPSolve()`s with different matrices unless you call it for each solve.
220b3402f20SBarry Smith 
22120f4b53cSBarry Smith   This routine is called each time the outer matrix is changed. In the first call the Ap argument is `NULL` and the routine should create the
222b3402f20SBarry Smith   matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix.
223b3402f20SBarry Smith 
224feefa0e1SJacob Faibussowitsch   Developer Notes:
225f1580f4eSBarry Smith   If the user does not call this routine nor call `PCGalerkinGetKSP()` and `KSPSetOperators()` then `PCGALERKIN`
226f1580f4eSBarry Smith   could automatically compute the submatrix via calls to `MatGalerkin()` or `MatRARt()`
227b3402f20SBarry Smith 
228562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetType()`, `PCType`, `PCGALERKIN`,
229db781477SPatrick Sanan           `PCGalerkinSetRestriction()`, `PCGalerkinSetInterpolation()`, `PCGalerkinGetKSP()`
230b3402f20SBarry Smith @*/
PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (* computeAsub)(PC pc,Mat A,Mat Ap,Mat * cAp,PetscCtx ctx),PetscCtx ctx)231*2a8381b2SBarry Smith PetscErrorCode PCGalerkinSetComputeSubmatrix(PC pc, PetscErrorCode (*computeAsub)(PC pc, Mat A, Mat Ap, Mat *cAp, PetscCtx ctx), PetscCtx ctx)
232d71ae5a4SJacob Faibussowitsch {
233b3402f20SBarry Smith   PetscFunctionBegin;
234b3402f20SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
235cac4c232SBarry Smith   PetscTryMethod(pc, "PCGalerkinSetComputeSubmatrix_C", (PC, PetscErrorCode (*)(PC, Mat, Mat, Mat *, void *), void *), (pc, computeAsub, ctx));
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
237b3402f20SBarry Smith }
238b3402f20SBarry Smith 
239b3402f20SBarry Smith /*@
240f1580f4eSBarry Smith   PCGalerkinGetKSP - Gets the `KSP` object in the `PCGALERKIN`
2412a6744ebSBarry Smith 
2422a6744ebSBarry Smith   Not Collective
2432a6744ebSBarry Smith 
2442a6744ebSBarry Smith   Input Parameter:
2452a6744ebSBarry Smith . pc - the preconditioner context
2462a6744ebSBarry Smith 
247f1580f4eSBarry Smith   Output Parameter:
248f1580f4eSBarry Smith . ksp - the `KSP` object
2492a6744ebSBarry Smith 
250feefa0e1SJacob Faibussowitsch   Level: intermediate
2512a6744ebSBarry Smith 
252f1580f4eSBarry Smith   Note:
25304c3f3b8SBarry Smith   Once you have called this routine you can call `KSPSetOperators()` on the resulting `KSP` to provide the operator for the Galerkin problem,
254f1580f4eSBarry Smith   an alternative is to use `PCGalerkinSetComputeSubmatrix()` to provide a routine that computes the submatrix as needed.
255b3402f20SBarry Smith 
256562efe2eSBarry Smith .seealso: [](ch_ksp), `PC`, `PCCreate()`, `PCSetType()`, `PCType`, `PCGALERKIN`,
257db781477SPatrick Sanan           `PCGalerkinSetRestriction()`, `PCGalerkinSetInterpolation()`, `PCGalerkinSetComputeSubmatrix()`
2582a6744ebSBarry Smith @*/
PCGalerkinGetKSP(PC pc,KSP * ksp)259d71ae5a4SJacob Faibussowitsch PetscErrorCode PCGalerkinGetKSP(PC pc, KSP *ksp)
260d71ae5a4SJacob Faibussowitsch {
2612a6744ebSBarry Smith   PetscFunctionBegin;
2620700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
2634f572ea9SToby Isaac   PetscAssertPointer(ksp, 2);
264cac4c232SBarry Smith   PetscUseMethod(pc, "PCGalerkinGetKSP_C", (PC, KSP *), (pc, ksp));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2662a6744ebSBarry Smith }
2672a6744ebSBarry Smith 
PCSetFromOptions_Galerkin(PC pc,PetscOptionItems PetscOptionsObject)268ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_Galerkin(PC pc, PetscOptionItems PetscOptionsObject)
269d71ae5a4SJacob Faibussowitsch {
2704ac220ccSBarry Smith   PC_Galerkin *jac = (PC_Galerkin *)pc->data;
2714ac220ccSBarry Smith   const char  *prefix;
2724ac220ccSBarry Smith   PetscBool    flg;
2734ac220ccSBarry Smith 
2744ac220ccSBarry Smith   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(jac->ksp, &prefix));
2769566063dSJacob Faibussowitsch   PetscCall(PetscStrendswith(prefix, "galerkin_", &flg));
2774ac220ccSBarry Smith   if (!flg) {
2789566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
2799566063dSJacob Faibussowitsch     PetscCall(KSPSetOptionsPrefix(jac->ksp, prefix));
2809566063dSJacob Faibussowitsch     PetscCall(KSPAppendOptionsPrefix(jac->ksp, "galerkin_"));
2814ac220ccSBarry Smith   }
2824ac220ccSBarry Smith 
283d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Galerkin options");
2841baa6e33SBarry Smith   if (jac->ksp) PetscCall(KSPSetFromOptions(jac->ksp));
285d0609cedSBarry Smith   PetscOptionsHeadEnd();
2863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2874ac220ccSBarry Smith }
2882a6744ebSBarry Smith 
2892a6744ebSBarry Smith /*MC
2902a6744ebSBarry Smith      PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T)
2912a6744ebSBarry Smith 
2922a6744ebSBarry Smith     Level: intermediate
2932a6744ebSBarry Smith 
29420f4b53cSBarry Smith     Note:
29520f4b53cSBarry Smith     Use
29620f4b53cSBarry Smith .vb
29720f4b53cSBarry Smith      `PCGalerkinSetRestriction`(pc,R) and/or `PCGalerkinSetInterpolation`(pc,P)
29820f4b53cSBarry Smith      `PCGalerkinGetKSP`(pc,&ksp);
29920f4b53cSBarry Smith      `KSPSetOperators`(ksp,A,....)
30020f4b53cSBarry Smith      ...
30120f4b53cSBarry Smith .ve
30220f4b53cSBarry Smith 
303f1580f4eSBarry Smith     Developer Notes:
304f1580f4eSBarry Smith     If `KSPSetOperators()` has not been called on the inner `KSP` then `PCGALERKIN` could use `MatRARt()` or `MatPtAP()` to compute
3057817a140SBarry Smith     the operators automatically.
306f1580f4eSBarry Smith 
307f1580f4eSBarry Smith     Should there be a prefix for the inner `KSP`?
308f1580f4eSBarry Smith 
309f1580f4eSBarry Smith     There is no `KSPSetFromOptions_Galerkin()` that calls `KSPSetFromOptions()` on the inner `KSP`
3107817a140SBarry Smith 
311562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
312db781477SPatrick Sanan           `PCSHELL`, `PCKSP`, `PCGalerkinSetRestriction()`, `PCGalerkinSetInterpolation()`, `PCGalerkinGetKSP()`
3132a6744ebSBarry Smith M*/
3142a6744ebSBarry Smith 
PCCreate_Galerkin(PC pc)315d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Galerkin(PC pc)
316d71ae5a4SJacob Faibussowitsch {
3172a6744ebSBarry Smith   PC_Galerkin *jac;
3182a6744ebSBarry Smith 
3192a6744ebSBarry Smith   PetscFunctionBegin;
3204dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
3212fa5cd67SKarl Rupp 
3222a6744ebSBarry Smith   pc->ops->apply           = PCApply_Galerkin;
3232a6744ebSBarry Smith   pc->ops->setup           = PCSetUp_Galerkin;
324a06653b4SBarry Smith   pc->ops->reset           = PCReset_Galerkin;
3252a6744ebSBarry Smith   pc->ops->destroy         = PCDestroy_Galerkin;
3262a6744ebSBarry Smith   pc->ops->view            = PCView_Galerkin;
3274ac220ccSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Galerkin;
3284ac220ccSBarry Smith   pc->ops->applyrichardson = NULL;
3292a6744ebSBarry Smith 
3309566063dSJacob Faibussowitsch   PetscCall(KSPCreate(PetscObjectComm((PetscObject)pc), &jac->ksp));
3313821be0aSBarry Smith   PetscCall(KSPSetNestLevel(jac->ksp, pc->kspnestlevel));
3329566063dSJacob Faibussowitsch   PetscCall(KSPSetErrorIfNotConverged(jac->ksp, pc->erroriffailure));
3339566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)jac->ksp, (PetscObject)pc, 1));
3342a6744ebSBarry Smith 
3352a6744ebSBarry Smith   pc->data = (void *)jac;
3362a6744ebSBarry Smith 
3379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGalerkinSetRestriction_C", PCGalerkinSetRestriction_Galerkin));
3389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGalerkinSetInterpolation_C", PCGalerkinSetInterpolation_Galerkin));
3399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGalerkinGetKSP_C", PCGalerkinGetKSP_Galerkin));
3409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCGalerkinSetComputeSubmatrix_C", PCGalerkinSetComputeSubmatrix_Galerkin));
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3422a6744ebSBarry Smith }
343