xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision c5eb91543d2ee8daf880d93389b892228ddada03)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith       Defines a preconditioner that can consist of a collection of PCs
54b9ad928SBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"   /*I "petscpc.h" I*/
74b9ad928SBarry Smith #include "petscksp.h"            /*I "petscksp.h" I*/
84b9ad928SBarry Smith 
94b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink;
104b9ad928SBarry Smith struct _PC_CompositeLink {
114b9ad928SBarry Smith   PC               pc;
124b9ad928SBarry Smith   PC_CompositeLink next;
13421e10b8SBarry Smith   PC_CompositeLink previous;
144b9ad928SBarry Smith };
154b9ad928SBarry Smith 
164b9ad928SBarry Smith typedef struct {
174b9ad928SBarry Smith   PC_CompositeLink head;
184b9ad928SBarry Smith   PCCompositeType  type;
194b9ad928SBarry Smith   Vec              work1;
204b9ad928SBarry Smith   Vec              work2;
214b9ad928SBarry Smith   PetscScalar      alpha;
224b9ad928SBarry Smith   PetscTruth       use_true_matrix;
234b9ad928SBarry Smith } PC_Composite;
244b9ad928SBarry Smith 
254b9ad928SBarry Smith #undef __FUNCT__
264b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative"
276849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
284b9ad928SBarry Smith {
29dfbe8321SBarry Smith   PetscErrorCode   ierr;
304b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
314b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
324b9ad928SBarry Smith   Mat              mat = pc->pmat;
334b9ad928SBarry Smith 
344b9ad928SBarry Smith   PetscFunctionBegin;
35e7e72b3dSBarry Smith   if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
364b9ad928SBarry Smith   if (next->next && !jac->work2) { /* allocate second work vector */
374b9ad928SBarry Smith     ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr);
384b9ad928SBarry Smith   }
39d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
404b9ad928SBarry Smith   if (jac->use_true_matrix) mat = pc->mat;
414b9ad928SBarry Smith   while (next->next) {
424b9ad928SBarry Smith     next = next->next;
434b9ad928SBarry Smith     ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
44efb30889SBarry Smith     ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr);
4551f519a2SBarry Smith     ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
46d8fd42c4SBarry Smith     ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr);
47efb30889SBarry Smith     ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
484b9ad928SBarry Smith   }
49421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
50421e10b8SBarry Smith     while (next->previous) {
51421e10b8SBarry Smith       next = next->previous;
52421e10b8SBarry Smith       ierr  = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
53421e10b8SBarry Smith       ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr);
54421e10b8SBarry Smith       ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
55421e10b8SBarry Smith       ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr);
56421e10b8SBarry Smith       ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
57421e10b8SBarry Smith     }
58421e10b8SBarry Smith   }
594b9ad928SBarry Smith   PetscFunctionReturn(0);
604b9ad928SBarry Smith }
614b9ad928SBarry Smith 
624b9ad928SBarry Smith /*
634b9ad928SBarry Smith     This is very special for a matrix of the form alpha I + R + S
644b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from
654b9ad928SBarry Smith alpha I + R
664b9ad928SBarry Smith */
674b9ad928SBarry Smith #undef __FUNCT__
684b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special"
696849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
704b9ad928SBarry Smith {
71dfbe8321SBarry Smith   PetscErrorCode   ierr;
724b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
734b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
744b9ad928SBarry Smith 
754b9ad928SBarry Smith   PetscFunctionBegin;
76e7e72b3dSBarry Smith   if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
77e7e72b3dSBarry Smith   if (!next->next || next->next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
784b9ad928SBarry Smith 
79d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
80d8fd42c4SBarry Smith   ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr);
814b9ad928SBarry Smith   PetscFunctionReturn(0);
824b9ad928SBarry Smith }
834b9ad928SBarry Smith 
844b9ad928SBarry Smith #undef __FUNCT__
854b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive"
866849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
874b9ad928SBarry Smith {
88dfbe8321SBarry Smith   PetscErrorCode   ierr;
894b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
904b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
914b9ad928SBarry Smith 
924b9ad928SBarry Smith   PetscFunctionBegin;
93e7e72b3dSBarry Smith   if (!next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
94d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
954b9ad928SBarry Smith   while (next->next) {
964b9ad928SBarry Smith     next = next->next;
9751f519a2SBarry Smith     ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
98d8fd42c4SBarry Smith     ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
99efb30889SBarry Smith     ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
1004b9ad928SBarry Smith   }
1014b9ad928SBarry Smith   PetscFunctionReturn(0);
1024b9ad928SBarry Smith }
1034b9ad928SBarry Smith 
1044b9ad928SBarry Smith #undef __FUNCT__
1054b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite"
1066849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc)
1074b9ad928SBarry Smith {
108dfbe8321SBarry Smith   PetscErrorCode   ierr;
1094b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
1104b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith   PetscFunctionBegin;
1134b9ad928SBarry Smith   if (!jac->work1) {
11423ce1328SBarry Smith    ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr);
1154b9ad928SBarry Smith   }
1164b9ad928SBarry Smith   while (next) {
1174b9ad928SBarry Smith     ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1184b9ad928SBarry Smith     next = next->next;
1194b9ad928SBarry Smith   }
1204b9ad928SBarry Smith   PetscFunctionReturn(0);
1214b9ad928SBarry Smith }
1224b9ad928SBarry Smith 
1234b9ad928SBarry Smith #undef __FUNCT__
1244b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite"
1256849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc)
1264b9ad928SBarry Smith {
1274b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
128dfbe8321SBarry Smith   PetscErrorCode   ierr;
1294b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
1324b9ad928SBarry Smith   while (next) {
1334b9ad928SBarry Smith     ierr = PCDestroy(next->pc);CHKERRQ(ierr);
1344b9ad928SBarry Smith     next = next->next;
1354b9ad928SBarry Smith   }
1364b9ad928SBarry Smith 
1374b9ad928SBarry Smith   if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);}
1384b9ad928SBarry Smith   if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);}
1394b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1404b9ad928SBarry Smith   PetscFunctionReturn(0);
1414b9ad928SBarry Smith }
1424b9ad928SBarry Smith 
1434b9ad928SBarry Smith #undef __FUNCT__
1444b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite"
1456849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc)
1464b9ad928SBarry Smith {
1474b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
148dfbe8321SBarry Smith   PetscErrorCode   ierr;
1499dcbbd2bSBarry Smith   PetscInt         nmax = 8,i;
1504b9ad928SBarry Smith   PC_CompositeLink next;
151e5999256SBarry Smith   char             *pcs[8];
1524b9ad928SBarry Smith   PetscTruth       flg;
1534b9ad928SBarry Smith 
1544b9ad928SBarry Smith   PetscFunctionBegin;
1554b9ad928SBarry Smith   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
1569dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
15751f519a2SBarry Smith     if (flg) {
15851f519a2SBarry Smith       ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr);
15951f519a2SBarry Smith     }
16090d69ab7SBarry Smith     flg  = PETSC_FALSE;
16190d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1624b9ad928SBarry Smith     if (flg) {
1634b9ad928SBarry Smith       ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr);
1644b9ad928SBarry Smith     }
1654b9ad928SBarry Smith     ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr);
1664b9ad928SBarry Smith     if (flg) {
1674b9ad928SBarry Smith       for (i=0; i<nmax; i++) {
1684b9ad928SBarry Smith         ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr);
1694b9ad928SBarry Smith       }
1704b9ad928SBarry Smith     }
1714b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1724b9ad928SBarry Smith 
1734b9ad928SBarry Smith   next = jac->head;
1744b9ad928SBarry Smith   while (next) {
1754b9ad928SBarry Smith     ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr);
1764b9ad928SBarry Smith     next = next->next;
1774b9ad928SBarry Smith   }
1784b9ad928SBarry Smith   PetscFunctionReturn(0);
1794b9ad928SBarry Smith }
1804b9ad928SBarry Smith 
1814b9ad928SBarry Smith #undef __FUNCT__
1824b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite"
1836849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
1844b9ad928SBarry Smith {
1854b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
186dfbe8321SBarry Smith   PetscErrorCode   ierr;
1874b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
18832077d6dSBarry Smith   PetscTruth       iascii;
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   PetscFunctionBegin;
1912692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19232077d6dSBarry Smith   if (iascii) {
1939dcbbd2bSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr);
1944b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr);
1954b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
1964b9ad928SBarry Smith   } else {
19765e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
1984b9ad928SBarry Smith   }
19932077d6dSBarry Smith   if (iascii) {
2004b9ad928SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2014b9ad928SBarry Smith   }
2024b9ad928SBarry Smith   while (next) {
2034b9ad928SBarry Smith     ierr = PCView(next->pc,viewer);CHKERRQ(ierr);
2044b9ad928SBarry Smith     next = next->next;
2054b9ad928SBarry Smith   }
20632077d6dSBarry Smith   if (iascii) {
2074b9ad928SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2084b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
2094b9ad928SBarry Smith   }
2104b9ad928SBarry Smith   PetscFunctionReturn(0);
2114b9ad928SBarry Smith }
2124b9ad928SBarry Smith 
2134b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
2144b9ad928SBarry Smith 
2154b9ad928SBarry Smith EXTERN_C_BEGIN
2164b9ad928SBarry Smith #undef __FUNCT__
2174b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite"
218dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
2194b9ad928SBarry Smith {
2204b9ad928SBarry Smith   PC_Composite *jac = (PC_Composite*)pc->data;
2214b9ad928SBarry Smith   PetscFunctionBegin;
2224b9ad928SBarry Smith   jac->alpha = alpha;
2234b9ad928SBarry Smith   PetscFunctionReturn(0);
2244b9ad928SBarry Smith }
2254b9ad928SBarry Smith EXTERN_C_END
2264b9ad928SBarry Smith 
2274b9ad928SBarry Smith EXTERN_C_BEGIN
2284b9ad928SBarry Smith #undef __FUNCT__
2294b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite"
230dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type)
2314b9ad928SBarry Smith {
2324b9ad928SBarry Smith   PetscFunctionBegin;
2334b9ad928SBarry Smith   if (type == PC_COMPOSITE_ADDITIVE) {
2344b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Additive;
235421e10b8SBarry Smith   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
2364b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Multiplicative;
2374b9ad928SBarry Smith   } else if (type ==  PC_COMPOSITE_SPECIAL) {
2384b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Special;
239e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
2404b9ad928SBarry Smith   PetscFunctionReturn(0);
2414b9ad928SBarry Smith }
2424b9ad928SBarry Smith EXTERN_C_END
2434b9ad928SBarry Smith 
2444b9ad928SBarry Smith EXTERN_C_BEGIN
2454b9ad928SBarry Smith #undef __FUNCT__
2464b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite"
247dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type)
2484b9ad928SBarry Smith {
2494b9ad928SBarry Smith   PC_Composite     *jac;
2505a9f2f41SSatish Balay   PC_CompositeLink next,ilink;
251dfbe8321SBarry Smith   PetscErrorCode   ierr;
25279416396SBarry Smith   PetscInt         cnt = 0;
2532dcb1b2aSMatthew Knepley   const char       *prefix;
2542dcb1b2aSMatthew Knepley   char             newprefix[8];
2554b9ad928SBarry Smith 
2564b9ad928SBarry Smith   PetscFunctionBegin;
25738f2d2fdSLisandro Dalcin   ierr        = PetscNewLog(pc,struct _PC_CompositeLink,&ilink);CHKERRQ(ierr);
2585a9f2f41SSatish Balay   ilink->next = 0;
2597adad957SLisandro Dalcin   ierr = PCCreate(((PetscObject)pc)->comm,&ilink->pc);CHKERRQ(ierr);
2604b9ad928SBarry Smith 
2614b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2624b9ad928SBarry Smith   next = jac->head;
2634b9ad928SBarry Smith   if (!next) {
2645a9f2f41SSatish Balay     jac->head       = ilink;
265421e10b8SBarry Smith     ilink->previous = PETSC_NULL;
2664b9ad928SBarry Smith   } else {
2674b9ad928SBarry Smith     cnt++;
2684b9ad928SBarry Smith     while (next->next) {
2694b9ad928SBarry Smith       next = next->next;
2704b9ad928SBarry Smith       cnt++;
2714b9ad928SBarry Smith     }
2725a9f2f41SSatish Balay     next->next      = ilink;
273421e10b8SBarry Smith     ilink->previous = next;
2744b9ad928SBarry Smith   }
2754b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2765a9f2f41SSatish Balay   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
27713f74950SBarry Smith   sprintf(newprefix,"sub_%d_",(int)cnt);
2785a9f2f41SSatish Balay   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
2794b9ad928SBarry Smith   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
2805a9f2f41SSatish Balay   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
2814b9ad928SBarry Smith   PetscFunctionReturn(0);
2824b9ad928SBarry Smith }
2834b9ad928SBarry Smith EXTERN_C_END
2844b9ad928SBarry Smith 
2854b9ad928SBarry Smith EXTERN_C_BEGIN
2864b9ad928SBarry Smith #undef __FUNCT__
2874b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite"
288dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
2894b9ad928SBarry Smith {
2904b9ad928SBarry Smith   PC_Composite     *jac;
2914b9ad928SBarry Smith   PC_CompositeLink next;
29279416396SBarry Smith   PetscInt         i;
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith   PetscFunctionBegin;
2954b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2964b9ad928SBarry Smith   next = jac->head;
2974b9ad928SBarry Smith   for (i=0; i<n; i++) {
298e7e72b3dSBarry Smith     if (!next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
2994b9ad928SBarry Smith     next = next->next;
3004b9ad928SBarry Smith   }
3014b9ad928SBarry Smith   *subpc = next->pc;
3024b9ad928SBarry Smith   PetscFunctionReturn(0);
3034b9ad928SBarry Smith }
3044b9ad928SBarry Smith EXTERN_C_END
3054b9ad928SBarry Smith 
3064b9ad928SBarry Smith EXTERN_C_BEGIN
3074b9ad928SBarry Smith #undef __FUNCT__
3084b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
309dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
3104b9ad928SBarry Smith {
3114b9ad928SBarry Smith   PC_Composite   *jac;
3124b9ad928SBarry Smith 
3134b9ad928SBarry Smith   PetscFunctionBegin;
3144b9ad928SBarry Smith   jac                  = (PC_Composite*)pc->data;
3154b9ad928SBarry Smith   jac->use_true_matrix = PETSC_TRUE;
3164b9ad928SBarry Smith   PetscFunctionReturn(0);
3174b9ad928SBarry Smith }
3184b9ad928SBarry Smith EXTERN_C_END
3194b9ad928SBarry Smith 
3204b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */
3214b9ad928SBarry Smith #undef __FUNCT__
3224b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType"
323f39d8e23SSatish Balay /*@
3244b9ad928SBarry Smith    PCCompositeSetType - Sets the type of composite preconditioner.
3254b9ad928SBarry Smith 
326ad4df100SBarry Smith    Logically Collective on PC
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith    Input Parameter:
3292a6744ebSBarry Smith +  pc - the preconditioner context
3302a6744ebSBarry Smith -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
3314b9ad928SBarry Smith 
3324b9ad928SBarry Smith    Options Database Key:
3334b9ad928SBarry Smith .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
3344b9ad928SBarry Smith 
3354b9ad928SBarry Smith    Level: Developer
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3384b9ad928SBarry Smith @*/
339dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
3404b9ad928SBarry Smith {
341dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
3424b9ad928SBarry Smith 
3434b9ad928SBarry Smith   PetscFunctionBegin;
3440700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
345*c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
346*c5eb9154SBarry Smith 
3474b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
3484b9ad928SBarry Smith   if (f) {
3494b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
3504b9ad928SBarry Smith   }
3514b9ad928SBarry Smith   PetscFunctionReturn(0);
3524b9ad928SBarry Smith }
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith #undef __FUNCT__
3554b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha"
356f39d8e23SSatish Balay /*@
3574b9ad928SBarry Smith    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
3584b9ad928SBarry Smith      for alphaI + R + S
3594b9ad928SBarry Smith 
360ad4df100SBarry Smith    Logically Collective on PC
3614b9ad928SBarry Smith 
3624b9ad928SBarry Smith    Input Parameter:
3634b9ad928SBarry Smith +  pc - the preconditioner context
3644b9ad928SBarry Smith -  alpha - scale on identity
3654b9ad928SBarry Smith 
3664b9ad928SBarry Smith    Level: Developer
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3694b9ad928SBarry Smith @*/
370dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
3714b9ad928SBarry Smith {
372dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscScalar);
3734b9ad928SBarry Smith 
3744b9ad928SBarry Smith   PetscFunctionBegin;
3750700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
376*c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(pc,alpha,2);
377*c5eb9154SBarry Smith 
3784b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
3794b9ad928SBarry Smith   if (f) {
3804b9ad928SBarry Smith     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
3814b9ad928SBarry Smith   }
3824b9ad928SBarry Smith   PetscFunctionReturn(0);
3834b9ad928SBarry Smith }
3844b9ad928SBarry Smith 
3854b9ad928SBarry Smith #undef __FUNCT__
3864b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC"
3874b9ad928SBarry Smith /*@C
3884b9ad928SBarry Smith    PCCompositeAddPC - Adds another PC to the composite PC.
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith    Collective on PC
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith    Input Parameters:
3932a6744ebSBarry Smith +  pc - the preconditioner context
3942a6744ebSBarry Smith -  type - the type of the new preconditioner
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith    Level: Developer
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith .keywords: PC, composite preconditioner, add
3994b9ad928SBarry Smith @*/
400dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
4014b9ad928SBarry Smith {
402dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCType);
4034b9ad928SBarry Smith 
4044b9ad928SBarry Smith   PetscFunctionBegin;
4050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4064b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4074b9ad928SBarry Smith   if (f) {
4084b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
4094b9ad928SBarry Smith   }
4104b9ad928SBarry Smith   PetscFunctionReturn(0);
4114b9ad928SBarry Smith }
4124b9ad928SBarry Smith 
4134b9ad928SBarry Smith #undef __FUNCT__
4144b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC"
415f39d8e23SSatish Balay /*@
4164b9ad928SBarry Smith    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith    Not Collective
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith    Input Parameter:
4212a6744ebSBarry Smith +  pc - the preconditioner context
4222a6744ebSBarry Smith -  n - the number of the pc requested
4234b9ad928SBarry Smith 
4244b9ad928SBarry Smith    Output Parameters:
4254b9ad928SBarry Smith .  subpc - the PC requested
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Level: Developer
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith .seealso: PCCompositeAddPC()
4324b9ad928SBarry Smith @*/
433dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
4344b9ad928SBarry Smith {
43513f74950SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
4364b9ad928SBarry Smith 
4374b9ad928SBarry Smith   PetscFunctionBegin;
4380700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4394482741eSBarry Smith   PetscValidPointer(subpc,3);
4404b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4414b9ad928SBarry Smith   if (f) {
4424b9ad928SBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
443e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
4444b9ad928SBarry Smith   PetscFunctionReturn(0);
4454b9ad928SBarry Smith }
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith #undef __FUNCT__
4484b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue"
4494b9ad928SBarry Smith /*@
4504b9ad928SBarry Smith    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
4514b9ad928SBarry Smith                       the matrix used to define the preconditioner) is used to compute
4524b9ad928SBarry Smith                       the residual when the multiplicative scheme is used.
4534b9ad928SBarry Smith 
454ad4df100SBarry Smith    Logically Collective on PC
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith    Input Parameters:
4574b9ad928SBarry Smith .  pc - the preconditioner context
4584b9ad928SBarry Smith 
4594b9ad928SBarry Smith    Options Database Key:
4604b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith    Note:
4634b9ad928SBarry Smith    For the common case in which the preconditioning and linear
4644b9ad928SBarry Smith    system matrices are identical, this routine is unnecessary.
4654b9ad928SBarry Smith 
4664b9ad928SBarry Smith    Level: Developer
4674b9ad928SBarry Smith 
4684b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag
4694b9ad928SBarry Smith 
4704b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
4714b9ad928SBarry Smith @*/
472dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
4734b9ad928SBarry Smith {
474dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC);
4754b9ad928SBarry Smith 
4764b9ad928SBarry Smith   PetscFunctionBegin;
4770700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4784b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
4794b9ad928SBarry Smith   if (f) {
4804b9ad928SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
4814b9ad928SBarry Smith   }
4824b9ad928SBarry Smith   PetscFunctionReturn(0);
4834b9ad928SBarry Smith }
4844b9ad928SBarry Smith 
4854b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/
4864b9ad928SBarry Smith 
4874b9ad928SBarry Smith /*MC
4884b9ad928SBarry Smith      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith    Options Database Keys:
49151f519a2SBarry Smith +  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
4924b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
49351f519a2SBarry Smith -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
4944b9ad928SBarry Smith 
4954b9ad928SBarry Smith    Level: intermediate
4964b9ad928SBarry Smith 
4974b9ad928SBarry Smith    Concepts: composing solvers
4984b9ad928SBarry Smith 
4994b9ad928SBarry Smith    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
5004b9ad928SBarry Smith           inner PCs to be PCKSP.
5014b9ad928SBarry Smith           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
50279416396SBarry Smith           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
5034b9ad928SBarry Smith 
5044b9ad928SBarry Smith 
5054b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
5064b9ad928SBarry Smith            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
5074b9ad928SBarry Smith            PCCompositeGetPC(), PCCompositeSetUseTrue()
5084b9ad928SBarry Smith 
5094b9ad928SBarry Smith M*/
5104b9ad928SBarry Smith 
5114b9ad928SBarry Smith EXTERN_C_BEGIN
5124b9ad928SBarry Smith #undef __FUNCT__
5134b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite"
514dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
5154b9ad928SBarry Smith {
516dfbe8321SBarry Smith   PetscErrorCode ierr;
5174b9ad928SBarry Smith   PC_Composite   *jac;
5184b9ad928SBarry Smith 
5194b9ad928SBarry Smith   PetscFunctionBegin;
52038f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr);
5214b9ad928SBarry Smith   pc->ops->apply              = PCApply_Composite_Additive;
5224b9ad928SBarry Smith   pc->ops->setup              = PCSetUp_Composite;
5234b9ad928SBarry Smith   pc->ops->destroy            = PCDestroy_Composite;
5244b9ad928SBarry Smith   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
5254b9ad928SBarry Smith   pc->ops->view               = PCView_Composite;
5264b9ad928SBarry Smith   pc->ops->applyrichardson    = 0;
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith   pc->data               = (void*)jac;
5294b9ad928SBarry Smith   jac->type              = PC_COMPOSITE_ADDITIVE;
5304b9ad928SBarry Smith   jac->work1             = 0;
5314b9ad928SBarry Smith   jac->work2             = 0;
5324b9ad928SBarry Smith   jac->head              = 0;
5334b9ad928SBarry Smith 
5344b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
5354b9ad928SBarry Smith                     PCCompositeSetType_Composite);CHKERRQ(ierr);
5364b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
5374b9ad928SBarry Smith                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
5384b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
5394b9ad928SBarry Smith                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
5404b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
5414b9ad928SBarry Smith                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
5424b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
5434b9ad928SBarry Smith                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
5444b9ad928SBarry Smith 
5454b9ad928SBarry Smith   PetscFunctionReturn(0);
5464b9ad928SBarry Smith }
5474b9ad928SBarry Smith EXTERN_C_END
5484b9ad928SBarry Smith 
549