xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 78a3110e6e5b4f375af22805359ffc6186cef0f0)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith       Defines a preconditioner that can consist of a collection of PCs
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcimpl.h>
5c6db04a5SJed Brown #include <petscksp.h> /*I "petscksp.h" I*/
64b9ad928SBarry Smith 
74b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink;
84b9ad928SBarry Smith struct _PC_CompositeLink {
94b9ad928SBarry Smith   PC               pc;
104b9ad928SBarry Smith   PC_CompositeLink next;
11421e10b8SBarry Smith   PC_CompositeLink previous;
124b9ad928SBarry Smith };
134b9ad928SBarry Smith 
144b9ad928SBarry Smith typedef struct {
154b9ad928SBarry Smith   PC_CompositeLink head;
164b9ad928SBarry Smith   PCCompositeType  type;
174b9ad928SBarry Smith   Vec              work1;
184b9ad928SBarry Smith   Vec              work2;
194b9ad928SBarry Smith   PetscScalar      alpha;
20*78a3110eSAlexander   Mat              alpha_mat;
214b9ad928SBarry Smith } PC_Composite;
224b9ad928SBarry Smith 
23d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Composite_Multiplicative(PC pc, Vec x, Vec y)
24d71ae5a4SJacob Faibussowitsch {
254b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
264b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
274b9ad928SBarry Smith   Mat              mat  = pc->pmat;
284b9ad928SBarry Smith 
294b9ad928SBarry Smith   PetscFunctionBegin;
3028b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
31450d59ebSPatrick Farrell 
32450d59ebSPatrick Farrell   /* Set the reuse flag on children PCs */
33450d59ebSPatrick Farrell   while (next) {
349566063dSJacob Faibussowitsch     PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
35450d59ebSPatrick Farrell     next = next->next;
36450d59ebSPatrick Farrell   }
37450d59ebSPatrick Farrell   next = jac->head;
38450d59ebSPatrick Farrell 
394b9ad928SBarry Smith   if (next->next && !jac->work2) { /* allocate second work vector */
409566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(jac->work1, &jac->work2));
414b9ad928SBarry Smith   }
4249517cdeSBarry Smith   if (pc->useAmat) mat = pc->mat;
439566063dSJacob Faibussowitsch   PetscCall(PCApply(next->pc, x, y)); /* y <- B x */
444b9ad928SBarry Smith   while (next->next) {
454b9ad928SBarry Smith     next = next->next;
469566063dSJacob Faibussowitsch     PetscCall(MatMult(mat, y, jac->work1));               /* work1 <- A y */
479566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x)); /* work2 <- x - work1 */
489566063dSJacob Faibussowitsch     PetscCall(PCApply(next->pc, jac->work2, jac->work1)); /* work1 <- C work2 */
499566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, 1.0, jac->work1));               /* y <- y + work1 = B x + C (x - A B x) = (B + C (1 - A B)) x */
504b9ad928SBarry Smith   }
51421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
52421e10b8SBarry Smith     while (next->previous) {
53421e10b8SBarry Smith       next = next->previous;
549566063dSJacob Faibussowitsch       PetscCall(MatMult(mat, y, jac->work1));
559566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
569566063dSJacob Faibussowitsch       PetscCall(PCApply(next->pc, jac->work2, jac->work1));
579566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, 1.0, jac->work1));
58421e10b8SBarry Smith     }
59421e10b8SBarry Smith   }
603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
614b9ad928SBarry Smith }
624b9ad928SBarry Smith 
63d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc, Vec x, Vec y)
64d71ae5a4SJacob Faibussowitsch {
652533e041SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
662533e041SBarry Smith   PC_CompositeLink next = jac->head;
672533e041SBarry Smith   Mat              mat  = pc->pmat;
682533e041SBarry Smith 
692533e041SBarry Smith   PetscFunctionBegin;
7028b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
712533e041SBarry Smith   if (next->next && !jac->work2) { /* allocate second work vector */
729566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(jac->work1, &jac->work2));
732533e041SBarry Smith   }
7449517cdeSBarry Smith   if (pc->useAmat) mat = pc->mat;
752533e041SBarry Smith   /* locate last PC */
76ad540459SPierre Jolivet   while (next->next) next = next->next;
779566063dSJacob Faibussowitsch   PetscCall(PCApplyTranspose(next->pc, x, y));
782533e041SBarry Smith   while (next->previous) {
792533e041SBarry Smith     next = next->previous;
809566063dSJacob Faibussowitsch     PetscCall(MatMultTranspose(mat, y, jac->work1));
819566063dSJacob Faibussowitsch     PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
829566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1));
839566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, 1.0, jac->work1));
842533e041SBarry Smith   }
852533e041SBarry Smith   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
862533e041SBarry Smith     next = jac->head;
872533e041SBarry Smith     while (next->next) {
882533e041SBarry Smith       next = next->next;
899566063dSJacob Faibussowitsch       PetscCall(MatMultTranspose(mat, y, jac->work1));
909566063dSJacob Faibussowitsch       PetscCall(VecWAXPY(jac->work2, -1.0, jac->work1, x));
919566063dSJacob Faibussowitsch       PetscCall(PCApplyTranspose(next->pc, jac->work2, jac->work1));
929566063dSJacob Faibussowitsch       PetscCall(VecAXPY(y, 1.0, jac->work1));
932533e041SBarry Smith     }
942533e041SBarry Smith   }
953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
962533e041SBarry Smith }
972533e041SBarry Smith 
984b9ad928SBarry Smith /*
994b9ad928SBarry Smith     This is very special for a matrix of the form alpha I + R + S
1004b9ad928SBarry Smith     where first preconditioner is built from alpha I + S and second from
1014b9ad928SBarry Smith     alpha I + R
1024b9ad928SBarry Smith */
103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Composite_Special(PC pc, Vec x, Vec y)
104d71ae5a4SJacob Faibussowitsch {
1054b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
1064b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith   PetscFunctionBegin;
10928b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
1107827d75bSBarry Smith   PetscCheck(next->next && !next->next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Special composite preconditioners requires exactly two PCs");
1114b9ad928SBarry Smith 
112450d59ebSPatrick Farrell   /* Set the reuse flag on children PCs */
1139566063dSJacob Faibussowitsch   PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
1149566063dSJacob Faibussowitsch   PetscCall(PCSetReusePreconditioner(next->next->pc, pc->reusepreconditioner));
115450d59ebSPatrick Farrell 
1169566063dSJacob Faibussowitsch   PetscCall(PCApply(next->pc, x, jac->work1));
117*78a3110eSAlexander   if (jac->alpha_mat) {
118*78a3110eSAlexander     if (!jac->work2) PetscCall(VecDuplicate(jac->work1, &jac->work2));
119*78a3110eSAlexander     PetscCall(MatMult(jac->alpha_mat, jac->work1, jac->work2));
120*78a3110eSAlexander     PetscCall(PCApply(next->next->pc, jac->work2, y));
121*78a3110eSAlexander   } else PetscCall(PCApply(next->next->pc, jac->work1, y));
1223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1234b9ad928SBarry Smith }
1244b9ad928SBarry Smith 
125d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_Composite_Additive(PC pc, Vec x, Vec y)
126d71ae5a4SJacob Faibussowitsch {
1274b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
1284b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1294b9ad928SBarry Smith 
1304b9ad928SBarry Smith   PetscFunctionBegin;
13128b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
132450d59ebSPatrick Farrell 
133450d59ebSPatrick Farrell   /* Set the reuse flag on children PCs */
134450d59ebSPatrick Farrell   while (next) {
1359566063dSJacob Faibussowitsch     PetscCall(PCSetReusePreconditioner(next->pc, pc->reusepreconditioner));
136450d59ebSPatrick Farrell     next = next->next;
137450d59ebSPatrick Farrell   }
138450d59ebSPatrick Farrell   next = jac->head;
139450d59ebSPatrick Farrell 
1409566063dSJacob Faibussowitsch   PetscCall(PCApply(next->pc, x, y));
1414b9ad928SBarry Smith   while (next->next) {
1424b9ad928SBarry Smith     next = next->next;
1439566063dSJacob Faibussowitsch     PetscCall(PCApply(next->pc, x, jac->work1));
1449566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, 1.0, jac->work1));
1454b9ad928SBarry Smith   }
1463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1474b9ad928SBarry Smith }
1484b9ad928SBarry Smith 
149d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc, Vec x, Vec y)
150d71ae5a4SJacob Faibussowitsch {
1512533e041SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
1522533e041SBarry Smith   PC_CompositeLink next = jac->head;
1532533e041SBarry Smith 
1542533e041SBarry Smith   PetscFunctionBegin;
15528b400f6SJacob Faibussowitsch   PetscCheck(next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs");
1569566063dSJacob Faibussowitsch   PetscCall(PCApplyTranspose(next->pc, x, y));
1572533e041SBarry Smith   while (next->next) {
1582533e041SBarry Smith     next = next->next;
1599566063dSJacob Faibussowitsch     PetscCall(PCApplyTranspose(next->pc, x, jac->work1));
1609566063dSJacob Faibussowitsch     PetscCall(VecAXPY(y, 1.0, jac->work1));
1612533e041SBarry Smith   }
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1632533e041SBarry Smith }
1642533e041SBarry Smith 
165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_Composite(PC pc)
166d71ae5a4SJacob Faibussowitsch {
1674b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
1684b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1695a78d018SMatthew G. Knepley   DM               dm;
1704b9ad928SBarry Smith 
1714b9ad928SBarry Smith   PetscFunctionBegin;
17248a46eb9SPierre Jolivet   if (!jac->work1) PetscCall(MatCreateVecs(pc->pmat, &jac->work1, NULL));
1739566063dSJacob Faibussowitsch   PetscCall(PCGetDM(pc, &dm));
1744b9ad928SBarry Smith   while (next) {
17548a46eb9SPierre Jolivet     if (!next->pc->dm) PetscCall(PCSetDM(next->pc, dm));
17648a46eb9SPierre Jolivet     if (!next->pc->mat) PetscCall(PCSetOperators(next->pc, pc->mat, pc->pmat));
1774b9ad928SBarry Smith     next = next->next;
1784b9ad928SBarry Smith   }
1793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1804b9ad928SBarry Smith }
1814b9ad928SBarry Smith 
1824f853519SStefano Zampini static PetscErrorCode PCSetUpOnBlocks_Composite(PC pc)
1834f853519SStefano Zampini {
1844f853519SStefano Zampini   PC_Composite    *jac  = (PC_Composite *)pc->data;
1854f853519SStefano Zampini   PC_CompositeLink next = jac->head;
1864f853519SStefano Zampini   PCFailedReason   reason;
1874f853519SStefano Zampini 
1884f853519SStefano Zampini   PetscFunctionBegin;
1894f853519SStefano Zampini   while (next) {
1904f853519SStefano Zampini     PetscCall(PCSetUp(next->pc));
1914f853519SStefano Zampini     PetscCall(PCGetFailedReasonRank(next->pc, &reason));
1924f853519SStefano Zampini     if (reason) pc->failedreason = reason;
1934f853519SStefano Zampini     next = next->next;
1944f853519SStefano Zampini   }
1954f853519SStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
1964f853519SStefano Zampini }
1974f853519SStefano Zampini 
198d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_Composite(PC pc)
199d71ae5a4SJacob Faibussowitsch {
20069d2c0f9SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
2015f48b12bSBarry Smith   PC_CompositeLink next = jac->head;
20269d2c0f9SBarry Smith 
20369d2c0f9SBarry Smith   PetscFunctionBegin;
20469d2c0f9SBarry Smith   while (next) {
2059566063dSJacob Faibussowitsch     PetscCall(PCReset(next->pc));
20669d2c0f9SBarry Smith     next = next->next;
20769d2c0f9SBarry Smith   }
2089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->work1));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&jac->work2));
210*78a3110eSAlexander   PetscCall(MatDestroy(&jac->alpha_mat));
2113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
21269d2c0f9SBarry Smith }
21369d2c0f9SBarry Smith 
214d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_Composite(PC pc)
215d71ae5a4SJacob Faibussowitsch {
2164b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
217724c2c99SHong Zhang   PC_CompositeLink next = jac->head, next_tmp;
2184b9ad928SBarry Smith 
2194b9ad928SBarry Smith   PetscFunctionBegin;
2209566063dSJacob Faibussowitsch   PetscCall(PCReset_Composite(pc));
2214b9ad928SBarry Smith   while (next) {
2229566063dSJacob Faibussowitsch     PetscCall(PCDestroy(&next->pc));
223724c2c99SHong Zhang     next_tmp = next;
2244b9ad928SBarry Smith     next     = next->next;
2259566063dSJacob Faibussowitsch     PetscCall(PetscFree(next_tmp));
2264b9ad928SBarry Smith   }
2272e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", NULL));
2282e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", NULL));
2292e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", NULL));
2302e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", NULL));
2312e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", NULL));
2322e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", NULL));
2332e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", NULL));
234*78a3110eSAlexander   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlphaMat_C", NULL));
2359566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
2363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2374b9ad928SBarry Smith }
2384b9ad928SBarry Smith 
239d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_Composite(PC pc, PetscOptionItems *PetscOptionsObject)
240d71ae5a4SJacob Faibussowitsch {
2414b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
2429dcbbd2bSBarry Smith   PetscInt         nmax = 8, i;
2434b9ad928SBarry Smith   PC_CompositeLink next;
244e5999256SBarry Smith   char            *pcs[8];
245ace3abfcSBarry Smith   PetscBool        flg;
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith   PetscFunctionBegin;
248d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Composite preconditioner options");
2499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_composite_type", "Type of composition", "PCCompositeSetType", PCCompositeTypes, (PetscEnum)jac->type, (PetscEnum *)&jac->type, &flg));
2501baa6e33SBarry Smith   if (flg) PetscCall(PCCompositeSetType(pc, jac->type));
2519566063dSJacob Faibussowitsch   PetscCall(PetscOptionsStringArray("-pc_composite_pcs", "List of composite solvers", "PCCompositeAddPCType", pcs, &nmax, &flg));
2524b9ad928SBarry Smith   if (flg) {
2534b9ad928SBarry Smith     for (i = 0; i < nmax; i++) {
2549566063dSJacob Faibussowitsch       PetscCall(PCCompositeAddPCType(pc, pcs[i]));
2559566063dSJacob Faibussowitsch       PetscCall(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */
2564b9ad928SBarry Smith     }
2574b9ad928SBarry Smith   }
258d0609cedSBarry Smith   PetscOptionsHeadEnd();
2594b9ad928SBarry Smith 
2604b9ad928SBarry Smith   next = jac->head;
2614b9ad928SBarry Smith   while (next) {
2629566063dSJacob Faibussowitsch     PetscCall(PCSetFromOptions(next->pc));
2634b9ad928SBarry Smith     next = next->next;
2644b9ad928SBarry Smith   }
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2664b9ad928SBarry Smith }
2674b9ad928SBarry Smith 
268d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_Composite(PC pc, PetscViewer viewer)
269d71ae5a4SJacob Faibussowitsch {
2704b9ad928SBarry Smith   PC_Composite    *jac  = (PC_Composite *)pc->data;
2714b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
272ace3abfcSBarry Smith   PetscBool        iascii;
2734b9ad928SBarry Smith 
2744b9ad928SBarry Smith   PetscFunctionBegin;
2759566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
27632077d6dSBarry Smith   if (iascii) {
2779566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Composite PC type - %s\n", PCCompositeTypes[jac->type]));
2789566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "PCs on composite preconditioner follow\n"));
2799566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n"));
2804b9ad928SBarry Smith   }
2811baa6e33SBarry Smith   if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer));
2824b9ad928SBarry Smith   while (next) {
2839566063dSJacob Faibussowitsch     PetscCall(PCView(next->pc, viewer));
2844b9ad928SBarry Smith     next = next->next;
2854b9ad928SBarry Smith   }
28632077d6dSBarry Smith   if (iascii) {
2879566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPopTab(viewer));
2889566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "---------------------------------\n"));
2894b9ad928SBarry Smith   }
2903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2914b9ad928SBarry Smith }
2924b9ad928SBarry Smith 
293d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc, PetscScalar alpha)
294d71ae5a4SJacob Faibussowitsch {
2954b9ad928SBarry Smith   PC_Composite *jac = (PC_Composite *)pc->data;
2965fd66863SKarl Rupp 
2974b9ad928SBarry Smith   PetscFunctionBegin;
2984b9ad928SBarry Smith   jac->alpha = alpha;
2993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3004b9ad928SBarry Smith }
3014b9ad928SBarry Smith 
302*78a3110eSAlexander static PetscErrorCode PCCompositeSpecialSetAlphaMat_Composite(PC pc, Mat alpha_mat)
303*78a3110eSAlexander {
304*78a3110eSAlexander   PC_Composite *jac = (PC_Composite *)pc->data;
305*78a3110eSAlexander 
306*78a3110eSAlexander   PetscFunctionBegin;
307*78a3110eSAlexander   if (alpha_mat) {
308*78a3110eSAlexander     PetscValidHeaderSpecific(alpha_mat, MAT_CLASSID, 2);
309*78a3110eSAlexander     PetscCall(PetscObjectReference((PetscObject)alpha_mat));
310*78a3110eSAlexander   }
311*78a3110eSAlexander   PetscCall(MatDestroy(&jac->alpha_mat));
312*78a3110eSAlexander   jac->alpha_mat = alpha_mat;
313*78a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
314*78a3110eSAlexander }
315*78a3110eSAlexander 
316d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeSetType_Composite(PC pc, PCCompositeType type)
317d71ae5a4SJacob Faibussowitsch {
318fad69fbaSJed Brown   PC_Composite *jac = (PC_Composite *)pc->data;
319fad69fbaSJed Brown 
3204b9ad928SBarry Smith   PetscFunctionBegin;
3214b9ad928SBarry Smith   if (type == PC_COMPOSITE_ADDITIVE) {
3224b9ad928SBarry Smith     pc->ops->apply          = PCApply_Composite_Additive;
3232533e041SBarry Smith     pc->ops->applytranspose = PCApplyTranspose_Composite_Additive;
324421e10b8SBarry Smith   } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
3254b9ad928SBarry Smith     pc->ops->apply          = PCApply_Composite_Multiplicative;
3262533e041SBarry Smith     pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative;
3274b9ad928SBarry Smith   } else if (type == PC_COMPOSITE_SPECIAL) {
3284b9ad928SBarry Smith     pc->ops->apply          = PCApply_Composite_Special;
3290298fd71SBarry Smith     pc->ops->applytranspose = NULL;
330a7261c6bSprj-   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Unknown composite preconditioner type");
331fad69fbaSJed Brown   jac->type = type;
3323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3334b9ad928SBarry Smith }
3344b9ad928SBarry Smith 
335d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeGetType_Composite(PC pc, PCCompositeType *type)
336d71ae5a4SJacob Faibussowitsch {
337c60c7ad4SBarry Smith   PC_Composite *jac = (PC_Composite *)pc->data;
338c60c7ad4SBarry Smith 
339c60c7ad4SBarry Smith   PetscFunctionBegin;
340c60c7ad4SBarry Smith   *type = jac->type;
3413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
342c60c7ad4SBarry Smith }
343c60c7ad4SBarry Smith 
344d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc)
345d71ae5a4SJacob Faibussowitsch {
3464b9ad928SBarry Smith   PC_Composite    *jac;
3475a9f2f41SSatish Balay   PC_CompositeLink next, ilink;
34879416396SBarry Smith   PetscInt         cnt = 0;
3492dcb1b2aSMatthew Knepley   const char      *prefix;
350d726e3a5SJed Brown   char             newprefix[20];
3514b9ad928SBarry Smith 
3524b9ad928SBarry Smith   PetscFunctionBegin;
3534dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&ilink));
3540a545947SLisandro Dalcin   ilink->next = NULL;
3558aa07aa6SMatthew G. Knepley   ilink->pc   = subpc;
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith   jac  = (PC_Composite *)pc->data;
3584b9ad928SBarry Smith   next = jac->head;
3594b9ad928SBarry Smith   if (!next) {
3605a9f2f41SSatish Balay     jac->head       = ilink;
3610298fd71SBarry Smith     ilink->previous = NULL;
3624b9ad928SBarry Smith   } else {
3634b9ad928SBarry Smith     cnt++;
3644b9ad928SBarry Smith     while (next->next) {
3654b9ad928SBarry Smith       next = next->next;
3664b9ad928SBarry Smith       cnt++;
3674b9ad928SBarry Smith     }
3685a9f2f41SSatish Balay     next->next      = ilink;
369421e10b8SBarry Smith     ilink->previous = next;
3704b9ad928SBarry Smith   }
3719566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
3729566063dSJacob Faibussowitsch   PetscCall(PCSetOptionsPrefix(subpc, prefix));
3739566063dSJacob Faibussowitsch   PetscCall(PetscSNPrintf(newprefix, 20, "sub_%d_", (int)cnt));
3749566063dSJacob Faibussowitsch   PetscCall(PCAppendOptionsPrefix(subpc, newprefix));
3759566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)subpc));
3763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3778aa07aa6SMatthew G. Knepley }
3788aa07aa6SMatthew G. Knepley 
379d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type)
380d71ae5a4SJacob Faibussowitsch {
3818aa07aa6SMatthew G. Knepley   PC subpc;
3828aa07aa6SMatthew G. Knepley 
3838aa07aa6SMatthew G. Knepley   PetscFunctionBegin;
3849566063dSJacob Faibussowitsch   PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), &subpc));
3859566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)subpc, (PetscObject)pc, 1));
3869566063dSJacob Faibussowitsch   PetscCall(PCCompositeAddPC_Composite(pc, subpc));
3874b9ad928SBarry Smith   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
3889566063dSJacob Faibussowitsch   PetscCall(PCSetType(subpc, type));
3899566063dSJacob Faibussowitsch   PetscCall(PCDestroy(&subpc));
3903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3914b9ad928SBarry Smith }
3924b9ad928SBarry Smith 
393d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc, PetscInt *n)
394d71ae5a4SJacob Faibussowitsch {
3958e6eba06SBarry Smith   PC_Composite    *jac;
3968e6eba06SBarry Smith   PC_CompositeLink next;
3978e6eba06SBarry Smith 
3988e6eba06SBarry Smith   PetscFunctionBegin;
3998e6eba06SBarry Smith   jac  = (PC_Composite *)pc->data;
4008e6eba06SBarry Smith   next = jac->head;
4018e6eba06SBarry Smith   *n   = 0;
4028e6eba06SBarry Smith   while (next) {
4038e6eba06SBarry Smith     next = next->next;
4048e6eba06SBarry Smith     (*n)++;
4058e6eba06SBarry Smith   }
4063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4078e6eba06SBarry Smith }
4088e6eba06SBarry Smith 
409d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCCompositeGetPC_Composite(PC pc, PetscInt n, PC *subpc)
410d71ae5a4SJacob Faibussowitsch {
4114b9ad928SBarry Smith   PC_Composite    *jac;
4124b9ad928SBarry Smith   PC_CompositeLink next;
41379416396SBarry Smith   PetscInt         i;
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith   PetscFunctionBegin;
4164b9ad928SBarry Smith   jac  = (PC_Composite *)pc->data;
4174b9ad928SBarry Smith   next = jac->head;
4184b9ad928SBarry Smith   for (i = 0; i < n; i++) {
41928b400f6SJacob Faibussowitsch     PetscCheck(next->next, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_INCOMP, "Not enough PCs in composite preconditioner");
4204b9ad928SBarry Smith     next = next->next;
4214b9ad928SBarry Smith   }
4224b9ad928SBarry Smith   *subpc = next->pc;
4233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4244b9ad928SBarry Smith }
4254b9ad928SBarry Smith 
426f39d8e23SSatish Balay /*@
4274b9ad928SBarry Smith   PCCompositeSetType - Sets the type of composite preconditioner.
4284b9ad928SBarry Smith 
429c3339decSBarry Smith   Logically Collective
4304b9ad928SBarry Smith 
431c60c7ad4SBarry Smith   Input Parameters:
4322a6744ebSBarry Smith + pc   - the preconditioner context
433f1580f4eSBarry Smith - type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith   Options Database Key:
4364b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
4374b9ad928SBarry Smith 
438f1580f4eSBarry Smith   Level: advanced
4394b9ad928SBarry Smith 
440562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
441f1580f4eSBarry Smith           `PCCompositeGetType()`
4424b9ad928SBarry Smith @*/
443d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeSetType(PC pc, PCCompositeType type)
444d71ae5a4SJacob Faibussowitsch {
4454b9ad928SBarry Smith   PetscFunctionBegin;
4460700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
447c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
448cac4c232SBarry Smith   PetscTryMethod(pc, "PCCompositeSetType_C", (PC, PCCompositeType), (pc, type));
4493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4504b9ad928SBarry Smith }
4514b9ad928SBarry Smith 
452c60c7ad4SBarry Smith /*@
453721f67b5SBarry Smith   PCCompositeGetType - Gets the type of composite preconditioner.
454c60c7ad4SBarry Smith 
455c3339decSBarry Smith   Logically Collective
456c60c7ad4SBarry Smith 
457c60c7ad4SBarry Smith   Input Parameter:
458c60c7ad4SBarry Smith . pc - the preconditioner context
459c60c7ad4SBarry Smith 
460c60c7ad4SBarry Smith   Output Parameter:
461f1580f4eSBarry Smith . type - `PC_COMPOSITE_ADDITIVE` (default), `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`
462c60c7ad4SBarry Smith 
463f1580f4eSBarry Smith   Level: advanced
464c60c7ad4SBarry Smith 
465562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
466f1580f4eSBarry Smith           `PCCompositeSetType()`
467c60c7ad4SBarry Smith @*/
468d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeGetType(PC pc, PCCompositeType *type)
469d71ae5a4SJacob Faibussowitsch {
470c60c7ad4SBarry Smith   PetscFunctionBegin;
471c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
472cac4c232SBarry Smith   PetscUseMethod(pc, "PCCompositeGetType_C", (PC, PCCompositeType *), (pc, type));
4733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
474c60c7ad4SBarry Smith }
475c60c7ad4SBarry Smith 
476f39d8e23SSatish Balay /*@
477f1580f4eSBarry Smith   PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner, `PC_COMPOSITE_SPECIAL`,
478af27ebaaSBarry Smith   for $\alpha I + R + S$
4794b9ad928SBarry Smith 
480c3339decSBarry Smith   Logically Collective
4814b9ad928SBarry Smith 
482d8d19677SJose E. Roman   Input Parameters:
4834b9ad928SBarry Smith + pc    - the preconditioner context
4844b9ad928SBarry Smith - alpha - scale on identity
4854b9ad928SBarry Smith 
486feefa0e1SJacob Faibussowitsch   Level: developer
4874b9ad928SBarry Smith 
488562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PC_COMPOSITE_ADDITIVE`, `PC_COMPOSITE_MULTIPLICATIVE`, `PC_COMPOSITE_SPECIAL`, `PCCompositeType`,
489f1580f4eSBarry Smith           `PCCompositeSetType()`, `PCCompositeGetType()`
4904b9ad928SBarry Smith @*/
491d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeSpecialSetAlpha(PC pc, PetscScalar alpha)
492d71ae5a4SJacob Faibussowitsch {
4934b9ad928SBarry Smith   PetscFunctionBegin;
4940700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
495c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(pc, alpha, 2);
496cac4c232SBarry Smith   PetscTryMethod(pc, "PCCompositeSpecialSetAlpha_C", (PC, PetscScalar), (pc, alpha));
4973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4984b9ad928SBarry Smith }
4994b9ad928SBarry Smith 
500*78a3110eSAlexander PetscErrorCode PCCompositeSpecialSetAlphaMat(PC pc, Mat alpha_mat)
501*78a3110eSAlexander {
502*78a3110eSAlexander   PetscFunctionBegin;
503*78a3110eSAlexander   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
504*78a3110eSAlexander   PetscTryMethod(pc, "PCCompositeSpecialSetAlphaMat_C", (PC, Mat), (pc, alpha_mat));
505*78a3110eSAlexander   PetscFunctionReturn(PETSC_SUCCESS);
506*78a3110eSAlexander }
507*78a3110eSAlexander 
5084b9ad928SBarry Smith /*@C
509f1580f4eSBarry Smith   PCCompositeAddPCType - Adds another `PC` of the given type to the composite `PC`.
5104b9ad928SBarry Smith 
511c3339decSBarry Smith   Collective
5124b9ad928SBarry Smith 
5134b9ad928SBarry Smith   Input Parameters:
5142a6744ebSBarry Smith + pc   - the preconditioner context
5152a6744ebSBarry Smith - type - the type of the new preconditioner
5164b9ad928SBarry Smith 
517f1580f4eSBarry Smith   Level: intermediate
5184b9ad928SBarry Smith 
519562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()`
5204b9ad928SBarry Smith @*/
521d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeAddPCType(PC pc, PCType type)
522d71ae5a4SJacob Faibussowitsch {
5234b9ad928SBarry Smith   PetscFunctionBegin;
5240700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
525cac4c232SBarry Smith   PetscTryMethod(pc, "PCCompositeAddPCType_C", (PC, PCType), (pc, type));
5263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5278aa07aa6SMatthew G. Knepley }
5288aa07aa6SMatthew G. Knepley 
5298aa07aa6SMatthew G. Knepley /*@
530f1580f4eSBarry Smith   PCCompositeAddPC - Adds another `PC` to the composite `PC`.
5318aa07aa6SMatthew G. Knepley 
532c3339decSBarry Smith   Collective
5338aa07aa6SMatthew G. Knepley 
5348aa07aa6SMatthew G. Knepley   Input Parameters:
5358aa07aa6SMatthew G. Knepley + pc    - the preconditioner context
5368aa07aa6SMatthew G. Knepley - subpc - the new preconditioner
5378aa07aa6SMatthew G. Knepley 
538f1580f4eSBarry Smith   Level: intermediate
5398aa07aa6SMatthew G. Knepley 
540562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`
5418aa07aa6SMatthew G. Knepley @*/
542d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeAddPC(PC pc, PC subpc)
543d71ae5a4SJacob Faibussowitsch {
5448aa07aa6SMatthew G. Knepley   PetscFunctionBegin;
5458aa07aa6SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5468aa07aa6SMatthew G. Knepley   PetscValidHeaderSpecific(subpc, PC_CLASSID, 2);
547cac4c232SBarry Smith   PetscTryMethod(pc, "PCCompositeAddPC_C", (PC, PC), (pc, subpc));
5483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5494b9ad928SBarry Smith }
5504b9ad928SBarry Smith 
5518e6eba06SBarry Smith /*@
552f1580f4eSBarry Smith   PCCompositeGetNumberPC - Gets the number of `PC` objects in the composite `PC`.
5538e6eba06SBarry Smith 
5548e6eba06SBarry Smith   Not Collective
5558e6eba06SBarry Smith 
5568e6eba06SBarry Smith   Input Parameter:
5578e6eba06SBarry Smith . pc - the preconditioner context
5588e6eba06SBarry Smith 
5598e6eba06SBarry Smith   Output Parameter:
5608e6eba06SBarry Smith . num - the number of sub pcs
5618e6eba06SBarry Smith 
562feefa0e1SJacob Faibussowitsch   Level: developer
5638e6eba06SBarry Smith 
564562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeGetPC()`, `PCCompositeAddPC()`, `PCCompositeAddPCType()`
5658e6eba06SBarry Smith @*/
566d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeGetNumberPC(PC pc, PetscInt *num)
567d71ae5a4SJacob Faibussowitsch {
5688e6eba06SBarry Smith   PetscFunctionBegin;
5698e6eba06SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5704f572ea9SToby Isaac   PetscAssertPointer(num, 2);
571cac4c232SBarry Smith   PetscUseMethod(pc, "PCCompositeGetNumberPC_C", (PC, PetscInt *), (pc, num));
5723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5738e6eba06SBarry Smith }
5748e6eba06SBarry Smith 
575f39d8e23SSatish Balay /*@
576f1580f4eSBarry Smith   PCCompositeGetPC - Gets one of the `PC` objects in the composite `PC`.
5774b9ad928SBarry Smith 
5784b9ad928SBarry Smith   Not Collective
5794b9ad928SBarry Smith 
580d8d19677SJose E. Roman   Input Parameters:
5812a6744ebSBarry Smith + pc - the preconditioner context
5822a6744ebSBarry Smith - n  - the number of the pc requested
5834b9ad928SBarry Smith 
584f1580f4eSBarry Smith   Output Parameter:
5854b9ad928SBarry Smith . subpc - the PC requested
5864b9ad928SBarry Smith 
587f1580f4eSBarry Smith   Level: intermediate
5884b9ad928SBarry Smith 
589f1580f4eSBarry Smith   Note:
590f1580f4eSBarry Smith   To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then
591f1580f4eSBarry Smith   call `PCSetOperators()` on that `PC`.
5922b1d202aSBarry Smith 
593562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCOMPOSITE`, `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`, `PCSetOperators()`
5944b9ad928SBarry Smith @*/
595d71ae5a4SJacob Faibussowitsch PetscErrorCode PCCompositeGetPC(PC pc, PetscInt n, PC *subpc)
596d71ae5a4SJacob Faibussowitsch {
5974b9ad928SBarry Smith   PetscFunctionBegin;
5980700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
5994f572ea9SToby Isaac   PetscAssertPointer(subpc, 3);
600cac4c232SBarry Smith   PetscUseMethod(pc, "PCCompositeGetPC_C", (PC, PetscInt, PC *), (pc, n, subpc));
6013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6024b9ad928SBarry Smith }
6034b9ad928SBarry Smith 
6044b9ad928SBarry Smith /*MC
6054b9ad928SBarry Smith      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
6064b9ad928SBarry Smith 
6074b9ad928SBarry Smith    Options Database Keys:
6082eab2d5bSJungho Lee +  -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type
609f1580f4eSBarry Smith .  -pc_use_amat                                                                                  - activates `PCSetUseAmat()`
61051f519a2SBarry Smith -  -pc_composite_pcs                                                                             - <pc0,pc1,...> list of PCs to compose
6114b9ad928SBarry Smith 
6124b9ad928SBarry Smith    Level: intermediate
6134b9ad928SBarry Smith 
61495452b02SPatrick Sanan    Notes:
615f1580f4eSBarry Smith    To use a Krylov method inside the composite preconditioner, set the `PCType` of one or more
616f1580f4eSBarry Smith    inner `PC`s to be `PCKSP`. Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
617f1580f4eSBarry Smith    the incorrect answer) unless you use `KSPFGMRES` as the outer Krylov method
618f1580f4eSBarry Smith 
619f1580f4eSBarry Smith    To use a different operator to construct one of the inner preconditioners first call `PCCompositeGetPC()`, then
620f1580f4eSBarry Smith    call `PCSetOperators()` on that `PC`.
6214b9ad928SBarry Smith 
622562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
623db781477SPatrick Sanan           `PCSHELL`, `PCKSP`, `PCCompositeSetType()`, `PCCompositeSpecialSetAlpha()`, `PCCompositeAddPCType()`,
624f1580f4eSBarry Smith           `PCCompositeGetPC()`, `PCSetUseAmat()`, `PCCompositeAddPC()`, `PCCompositeGetNumberPC()`
6254b9ad928SBarry Smith M*/
6264b9ad928SBarry Smith 
627d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc)
628d71ae5a4SJacob Faibussowitsch {
6294b9ad928SBarry Smith   PC_Composite *jac;
6304b9ad928SBarry Smith 
6314b9ad928SBarry Smith   PetscFunctionBegin;
6324dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&jac));
6332fa5cd67SKarl Rupp 
6344b9ad928SBarry Smith   pc->ops->apply           = PCApply_Composite_Additive;
6352533e041SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_Composite_Additive;
6364b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Composite;
6374f853519SStefano Zampini   pc->ops->setuponblocks   = PCSetUpOnBlocks_Composite;
63869d2c0f9SBarry Smith   pc->ops->reset           = PCReset_Composite;
6394b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Composite;
6404b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Composite;
6414b9ad928SBarry Smith   pc->ops->view            = PCView_Composite;
6420a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
6434b9ad928SBarry Smith 
6444b9ad928SBarry Smith   pc->data       = (void *)jac;
6454b9ad928SBarry Smith   jac->type      = PC_COMPOSITE_ADDITIVE;
6460a545947SLisandro Dalcin   jac->work1     = NULL;
6470a545947SLisandro Dalcin   jac->work2     = NULL;
6480a545947SLisandro Dalcin   jac->head      = NULL;
649*78a3110eSAlexander   jac->alpha_mat = NULL;
6504b9ad928SBarry Smith 
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSetType_C", PCCompositeSetType_Composite));
6529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetType_C", PCCompositeGetType_Composite));
6539566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPCType_C", PCCompositeAddPCType_Composite));
6549566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeAddPC_C", PCCompositeAddPC_Composite));
6559566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetNumberPC_C", PCCompositeGetNumberPC_Composite));
6569566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeGetPC_C", PCCompositeGetPC_Composite));
6579566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlpha_C", PCCompositeSpecialSetAlpha_Composite));
658*78a3110eSAlexander   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCCompositeSpecialSetAlphaMat_C", PCCompositeSpecialSetAlphaMat_Composite));
6593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6604b9ad928SBarry Smith }
661