xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 421e10b8f79bbed49dcc4e803c884835d979c6ea)
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;
13*421e10b8SBarry 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;
354b9ad928SBarry Smith   if (!next) {
36*421e10b8SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
374b9ad928SBarry Smith   }
384b9ad928SBarry Smith   if (next->next && !jac->work2) { /* allocate second work vector */
394b9ad928SBarry Smith     ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr);
404b9ad928SBarry Smith   }
41d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
424b9ad928SBarry Smith   if (jac->use_true_matrix) mat = pc->mat;
434b9ad928SBarry Smith   while (next->next) {
444b9ad928SBarry Smith     next = next->next;
454b9ad928SBarry Smith     ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
46efb30889SBarry Smith     ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr);
4751f519a2SBarry Smith     ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
48d8fd42c4SBarry Smith     ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr);
49efb30889SBarry Smith     ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
504b9ad928SBarry Smith   }
51*421e10b8SBarry Smith   if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
52*421e10b8SBarry Smith     while (next->previous) {
53*421e10b8SBarry Smith       next = next->previous;
54*421e10b8SBarry Smith       ierr  = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
55*421e10b8SBarry Smith       ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr);
56*421e10b8SBarry Smith       ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
57*421e10b8SBarry Smith       ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr);
58*421e10b8SBarry Smith       ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
59*421e10b8SBarry Smith     }
60*421e10b8SBarry Smith   }
614b9ad928SBarry Smith   PetscFunctionReturn(0);
624b9ad928SBarry Smith }
634b9ad928SBarry Smith 
644b9ad928SBarry Smith /*
654b9ad928SBarry Smith     This is very special for a matrix of the form alpha I + R + S
664b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from
674b9ad928SBarry Smith alpha I + R
684b9ad928SBarry Smith */
694b9ad928SBarry Smith #undef __FUNCT__
704b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special"
716849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
724b9ad928SBarry Smith {
73dfbe8321SBarry Smith   PetscErrorCode   ierr;
744b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
754b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
764b9ad928SBarry Smith 
774b9ad928SBarry Smith   PetscFunctionBegin;
784b9ad928SBarry Smith   if (!next) {
79*421e10b8SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
804b9ad928SBarry Smith   }
814b9ad928SBarry Smith   if (!next->next || next->next->next) {
821302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
834b9ad928SBarry Smith   }
844b9ad928SBarry Smith 
85d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
86d8fd42c4SBarry Smith   ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr);
874b9ad928SBarry Smith   PetscFunctionReturn(0);
884b9ad928SBarry Smith }
894b9ad928SBarry Smith 
904b9ad928SBarry Smith #undef __FUNCT__
914b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive"
926849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
934b9ad928SBarry Smith {
94dfbe8321SBarry Smith   PetscErrorCode   ierr;
954b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
964b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
974b9ad928SBarry Smith 
984b9ad928SBarry Smith   PetscFunctionBegin;
994b9ad928SBarry Smith   if (!next) {
100*421e10b8SBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs");
1014b9ad928SBarry Smith   }
102d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
1034b9ad928SBarry Smith   while (next->next) {
1044b9ad928SBarry Smith     next = next->next;
10551f519a2SBarry Smith     ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr);  /* zero since some PC's may not set all entries in the result */
106d8fd42c4SBarry Smith     ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
107efb30889SBarry Smith     ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
1084b9ad928SBarry Smith   }
1094b9ad928SBarry Smith   PetscFunctionReturn(0);
1104b9ad928SBarry Smith }
1114b9ad928SBarry Smith 
1124b9ad928SBarry Smith #undef __FUNCT__
1134b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite"
1146849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc)
1154b9ad928SBarry Smith {
116dfbe8321SBarry Smith   PetscErrorCode   ierr;
1174b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
1184b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith   PetscFunctionBegin;
1214b9ad928SBarry Smith   if (!jac->work1) {
12223ce1328SBarry Smith    ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr);
1234b9ad928SBarry Smith   }
1244b9ad928SBarry Smith   while (next) {
1254b9ad928SBarry Smith     ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1264b9ad928SBarry Smith     next = next->next;
1274b9ad928SBarry Smith   }
1284b9ad928SBarry Smith   PetscFunctionReturn(0);
1294b9ad928SBarry Smith }
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith #undef __FUNCT__
1324b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite"
1336849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc)
1344b9ad928SBarry Smith {
1354b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
136dfbe8321SBarry Smith   PetscErrorCode   ierr;
1374b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   PetscFunctionBegin;
1404b9ad928SBarry Smith   while (next) {
1414b9ad928SBarry Smith     ierr = PCDestroy(next->pc);CHKERRQ(ierr);
1424b9ad928SBarry Smith     next = next->next;
1434b9ad928SBarry Smith   }
1444b9ad928SBarry Smith 
1454b9ad928SBarry Smith   if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);}
1464b9ad928SBarry Smith   if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);}
1474b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1484b9ad928SBarry Smith   PetscFunctionReturn(0);
1494b9ad928SBarry Smith }
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith #undef __FUNCT__
1524b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite"
1536849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc)
1544b9ad928SBarry Smith {
1554b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
156dfbe8321SBarry Smith   PetscErrorCode   ierr;
1579dcbbd2bSBarry Smith   PetscInt         nmax = 8,i;
1584b9ad928SBarry Smith   PC_CompositeLink next;
159e5999256SBarry Smith   char             *pcs[8];
1604b9ad928SBarry Smith   PetscTruth       flg;
1614b9ad928SBarry Smith 
1624b9ad928SBarry Smith   PetscFunctionBegin;
1634b9ad928SBarry Smith   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
1649dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
16551f519a2SBarry Smith     if (flg) {
16651f519a2SBarry Smith       ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr);
16751f519a2SBarry Smith     }
1684b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr);
1694b9ad928SBarry Smith     if (flg) {
1704b9ad928SBarry Smith       ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr);
1714b9ad928SBarry Smith     }
1724b9ad928SBarry Smith     ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr);
1734b9ad928SBarry Smith     if (flg) {
1744b9ad928SBarry Smith       for (i=0; i<nmax; i++) {
1754b9ad928SBarry Smith         ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr);
1764b9ad928SBarry Smith       }
1774b9ad928SBarry Smith     }
1784b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1794b9ad928SBarry Smith 
1804b9ad928SBarry Smith   next = jac->head;
1814b9ad928SBarry Smith   while (next) {
1824b9ad928SBarry Smith     ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr);
1834b9ad928SBarry Smith     next = next->next;
1844b9ad928SBarry Smith   }
1854b9ad928SBarry Smith   PetscFunctionReturn(0);
1864b9ad928SBarry Smith }
1874b9ad928SBarry Smith 
1884b9ad928SBarry Smith #undef __FUNCT__
1894b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite"
1906849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
1914b9ad928SBarry Smith {
1924b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
193dfbe8321SBarry Smith   PetscErrorCode   ierr;
1944b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
19532077d6dSBarry Smith   PetscTruth       iascii;
1964b9ad928SBarry Smith 
1974b9ad928SBarry Smith   PetscFunctionBegin;
19832077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
19932077d6dSBarry Smith   if (iascii) {
2009dcbbd2bSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr);
2014b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr);
2024b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
2034b9ad928SBarry Smith   } else {
20479a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
2054b9ad928SBarry Smith   }
20632077d6dSBarry Smith   if (iascii) {
2074b9ad928SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2084b9ad928SBarry Smith   }
2094b9ad928SBarry Smith   while (next) {
2104b9ad928SBarry Smith     ierr = PCView(next->pc,viewer);CHKERRQ(ierr);
2114b9ad928SBarry Smith     next = next->next;
2124b9ad928SBarry Smith   }
21332077d6dSBarry Smith   if (iascii) {
2144b9ad928SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2154b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
2164b9ad928SBarry Smith   }
2174b9ad928SBarry Smith   PetscFunctionReturn(0);
2184b9ad928SBarry Smith }
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
2214b9ad928SBarry Smith 
2224b9ad928SBarry Smith EXTERN_C_BEGIN
2234b9ad928SBarry Smith #undef __FUNCT__
2244b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite"
225dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
2264b9ad928SBarry Smith {
2274b9ad928SBarry Smith   PC_Composite *jac = (PC_Composite*)pc->data;
2284b9ad928SBarry Smith   PetscFunctionBegin;
2294b9ad928SBarry Smith   jac->alpha = alpha;
2304b9ad928SBarry Smith   PetscFunctionReturn(0);
2314b9ad928SBarry Smith }
2324b9ad928SBarry Smith EXTERN_C_END
2334b9ad928SBarry Smith 
2344b9ad928SBarry Smith EXTERN_C_BEGIN
2354b9ad928SBarry Smith #undef __FUNCT__
2364b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite"
237dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type)
2384b9ad928SBarry Smith {
2394b9ad928SBarry Smith   PetscFunctionBegin;
2404b9ad928SBarry Smith   if (type == PC_COMPOSITE_ADDITIVE) {
2414b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Additive;
242*421e10b8SBarry Smith   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
2434b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Multiplicative;
2444b9ad928SBarry Smith   } else if (type ==  PC_COMPOSITE_SPECIAL) {
2454b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Special;
2464b9ad928SBarry Smith   } else {
2471302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
2484b9ad928SBarry Smith   }
2494b9ad928SBarry Smith   PetscFunctionReturn(0);
2504b9ad928SBarry Smith }
2514b9ad928SBarry Smith EXTERN_C_END
2524b9ad928SBarry Smith 
2534b9ad928SBarry Smith EXTERN_C_BEGIN
2544b9ad928SBarry Smith #undef __FUNCT__
2554b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite"
256dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type)
2574b9ad928SBarry Smith {
2584b9ad928SBarry Smith   PC_Composite     *jac;
2595a9f2f41SSatish Balay   PC_CompositeLink next,ilink;
260dfbe8321SBarry Smith   PetscErrorCode   ierr;
26179416396SBarry Smith   PetscInt         cnt = 0;
2622dcb1b2aSMatthew Knepley   const char       *prefix;
2632dcb1b2aSMatthew Knepley   char             newprefix[8];
2644b9ad928SBarry Smith 
2654b9ad928SBarry Smith   PetscFunctionBegin;
2665a9f2f41SSatish Balay   ierr       = PetscNew(struct _PC_CompositeLink,&ilink);CHKERRQ(ierr);
2675a9f2f41SSatish Balay   ilink->next = 0;
2685a9f2f41SSatish Balay   ierr = PCCreate(pc->comm,&ilink->pc);CHKERRQ(ierr);
2694b9ad928SBarry Smith 
2704b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2714b9ad928SBarry Smith   next = jac->head;
2724b9ad928SBarry Smith   if (!next) {
2735a9f2f41SSatish Balay     jac->head       = ilink;
274*421e10b8SBarry Smith     ilink->previous = PETSC_NULL;
2754b9ad928SBarry Smith   } else {
2764b9ad928SBarry Smith     cnt++;
2774b9ad928SBarry Smith     while (next->next) {
2784b9ad928SBarry Smith       next = next->next;
2794b9ad928SBarry Smith       cnt++;
2804b9ad928SBarry Smith     }
2815a9f2f41SSatish Balay     next->next      = ilink;
282*421e10b8SBarry Smith     ilink->previous = next;
2834b9ad928SBarry Smith   }
2844b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2855a9f2f41SSatish Balay   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
28613f74950SBarry Smith   sprintf(newprefix,"sub_%d_",(int)cnt);
2875a9f2f41SSatish Balay   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
2884b9ad928SBarry Smith   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
2895a9f2f41SSatish Balay   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
2904b9ad928SBarry Smith   PetscFunctionReturn(0);
2914b9ad928SBarry Smith }
2924b9ad928SBarry Smith EXTERN_C_END
2934b9ad928SBarry Smith 
2944b9ad928SBarry Smith EXTERN_C_BEGIN
2954b9ad928SBarry Smith #undef __FUNCT__
2964b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite"
297dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
2984b9ad928SBarry Smith {
2994b9ad928SBarry Smith   PC_Composite     *jac;
3004b9ad928SBarry Smith   PC_CompositeLink next;
30179416396SBarry Smith   PetscInt         i;
3024b9ad928SBarry Smith 
3034b9ad928SBarry Smith   PetscFunctionBegin;
3044b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
3054b9ad928SBarry Smith   next = jac->head;
3064b9ad928SBarry Smith   for (i=0; i<n; i++) {
3074b9ad928SBarry Smith     if (!next->next) {
3081302d50aSBarry Smith       SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
3094b9ad928SBarry Smith     }
3104b9ad928SBarry Smith     next = next->next;
3114b9ad928SBarry Smith   }
3124b9ad928SBarry Smith   *subpc = next->pc;
3134b9ad928SBarry Smith   PetscFunctionReturn(0);
3144b9ad928SBarry Smith }
3154b9ad928SBarry Smith EXTERN_C_END
3164b9ad928SBarry Smith 
3174b9ad928SBarry Smith EXTERN_C_BEGIN
3184b9ad928SBarry Smith #undef __FUNCT__
3194b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
320dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
3214b9ad928SBarry Smith {
3224b9ad928SBarry Smith   PC_Composite   *jac;
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith   PetscFunctionBegin;
3254b9ad928SBarry Smith   jac                  = (PC_Composite*)pc->data;
3264b9ad928SBarry Smith   jac->use_true_matrix = PETSC_TRUE;
3274b9ad928SBarry Smith   PetscFunctionReturn(0);
3284b9ad928SBarry Smith }
3294b9ad928SBarry Smith EXTERN_C_END
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */
3324b9ad928SBarry Smith #undef __FUNCT__
3334b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType"
334f39d8e23SSatish Balay /*@
3354b9ad928SBarry Smith    PCCompositeSetType - Sets the type of composite preconditioner.
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith    Collective on PC
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith    Input Parameter:
3402a6744ebSBarry Smith +  pc - the preconditioner context
3412a6744ebSBarry Smith -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
3424b9ad928SBarry Smith 
3434b9ad928SBarry Smith    Options Database Key:
3444b9ad928SBarry Smith .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
3454b9ad928SBarry Smith 
3464b9ad928SBarry Smith    Level: Developer
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3494b9ad928SBarry Smith @*/
350dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
3514b9ad928SBarry Smith {
352dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith   PetscFunctionBegin;
3554482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3564b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
3574b9ad928SBarry Smith   if (f) {
3584b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
3594b9ad928SBarry Smith   }
3604b9ad928SBarry Smith   PetscFunctionReturn(0);
3614b9ad928SBarry Smith }
3624b9ad928SBarry Smith 
3634b9ad928SBarry Smith #undef __FUNCT__
3644b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha"
365f39d8e23SSatish Balay /*@
3664b9ad928SBarry Smith    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
3674b9ad928SBarry Smith      for alphaI + R + S
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith    Collective on PC
3704b9ad928SBarry Smith 
3714b9ad928SBarry Smith    Input Parameter:
3724b9ad928SBarry Smith +  pc - the preconditioner context
3734b9ad928SBarry Smith -  alpha - scale on identity
3744b9ad928SBarry Smith 
3754b9ad928SBarry Smith    Level: Developer
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3784b9ad928SBarry Smith @*/
379dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
3804b9ad928SBarry Smith {
381dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscScalar);
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith   PetscFunctionBegin;
3844482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3854b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
3864b9ad928SBarry Smith   if (f) {
3874b9ad928SBarry Smith     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
3884b9ad928SBarry Smith   }
3894b9ad928SBarry Smith   PetscFunctionReturn(0);
3904b9ad928SBarry Smith }
3914b9ad928SBarry Smith 
3924b9ad928SBarry Smith #undef __FUNCT__
3934b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC"
3944b9ad928SBarry Smith /*@C
3954b9ad928SBarry Smith    PCCompositeAddPC - Adds another PC to the composite PC.
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith    Collective on PC
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith    Input Parameters:
4002a6744ebSBarry Smith +  pc - the preconditioner context
4012a6744ebSBarry Smith -  type - the type of the new preconditioner
4024b9ad928SBarry Smith 
4034b9ad928SBarry Smith    Level: Developer
4044b9ad928SBarry Smith 
4054b9ad928SBarry Smith .keywords: PC, composite preconditioner, add
4064b9ad928SBarry Smith @*/
407dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
4084b9ad928SBarry Smith {
409dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCType);
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith   PetscFunctionBegin;
4124482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4134b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4144b9ad928SBarry Smith   if (f) {
4154b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
4164b9ad928SBarry Smith   }
4174b9ad928SBarry Smith   PetscFunctionReturn(0);
4184b9ad928SBarry Smith }
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith #undef __FUNCT__
4214b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC"
422f39d8e23SSatish Balay /*@
4234b9ad928SBarry Smith    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith    Not Collective
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Input Parameter:
4282a6744ebSBarry Smith +  pc - the preconditioner context
4292a6744ebSBarry Smith -  n - the number of the pc requested
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith    Output Parameters:
4324b9ad928SBarry Smith .  subpc - the PC requested
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith    Level: Developer
4354b9ad928SBarry Smith 
4364b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith .seealso: PCCompositeAddPC()
4394b9ad928SBarry Smith @*/
440dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
4414b9ad928SBarry Smith {
44213f74950SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith   PetscFunctionBegin;
4454482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4464482741eSBarry Smith   PetscValidPointer(subpc,3);
4474b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4484b9ad928SBarry Smith   if (f) {
4494b9ad928SBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
4504b9ad928SBarry Smith   } else {
4511302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
4524b9ad928SBarry Smith   }
4534b9ad928SBarry Smith   PetscFunctionReturn(0);
4544b9ad928SBarry Smith }
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith #undef __FUNCT__
4574b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue"
4584b9ad928SBarry Smith /*@
4594b9ad928SBarry Smith    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
4604b9ad928SBarry Smith                       the matrix used to define the preconditioner) is used to compute
4614b9ad928SBarry Smith                       the residual when the multiplicative scheme is used.
4624b9ad928SBarry Smith 
4634b9ad928SBarry Smith    Collective on PC
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith    Input Parameters:
4664b9ad928SBarry Smith .  pc - the preconditioner context
4674b9ad928SBarry Smith 
4684b9ad928SBarry Smith    Options Database Key:
4694b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
4704b9ad928SBarry Smith 
4714b9ad928SBarry Smith    Note:
4724b9ad928SBarry Smith    For the common case in which the preconditioning and linear
4734b9ad928SBarry Smith    system matrices are identical, this routine is unnecessary.
4744b9ad928SBarry Smith 
4754b9ad928SBarry Smith    Level: Developer
4764b9ad928SBarry Smith 
4774b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag
4784b9ad928SBarry Smith 
4794b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
4804b9ad928SBarry Smith @*/
481dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
4824b9ad928SBarry Smith {
483dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC);
4844b9ad928SBarry Smith 
4854b9ad928SBarry Smith   PetscFunctionBegin;
4864482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4874b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
4884b9ad928SBarry Smith   if (f) {
4894b9ad928SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
4904b9ad928SBarry Smith   }
4914b9ad928SBarry Smith   PetscFunctionReturn(0);
4924b9ad928SBarry Smith }
4934b9ad928SBarry Smith 
4944b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/
4954b9ad928SBarry Smith 
4964b9ad928SBarry Smith /*MC
4974b9ad928SBarry Smith      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
4984b9ad928SBarry Smith 
4994b9ad928SBarry Smith    Options Database Keys:
50051f519a2SBarry Smith +  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
5014b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
50251f519a2SBarry Smith -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
5034b9ad928SBarry Smith 
5044b9ad928SBarry Smith    Level: intermediate
5054b9ad928SBarry Smith 
5064b9ad928SBarry Smith    Concepts: composing solvers
5074b9ad928SBarry Smith 
5084b9ad928SBarry Smith    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
5094b9ad928SBarry Smith           inner PCs to be PCKSP.
5104b9ad928SBarry Smith           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
51179416396SBarry Smith           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
5124b9ad928SBarry Smith 
5134b9ad928SBarry Smith 
5144b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
5154b9ad928SBarry Smith            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
5164b9ad928SBarry Smith            PCCompositeGetPC(), PCCompositeSetUseTrue()
5174b9ad928SBarry Smith 
5184b9ad928SBarry Smith M*/
5194b9ad928SBarry Smith 
5204b9ad928SBarry Smith EXTERN_C_BEGIN
5214b9ad928SBarry Smith #undef __FUNCT__
5224b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite"
523dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
5244b9ad928SBarry Smith {
525dfbe8321SBarry Smith   PetscErrorCode ierr;
5264b9ad928SBarry Smith   PC_Composite   *jac;
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith   PetscFunctionBegin;
5294b9ad928SBarry Smith   ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr);
53052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_Composite));CHKERRQ(ierr);
5314b9ad928SBarry Smith   pc->ops->apply              = PCApply_Composite_Additive;
5324b9ad928SBarry Smith   pc->ops->setup              = PCSetUp_Composite;
5334b9ad928SBarry Smith   pc->ops->destroy            = PCDestroy_Composite;
5344b9ad928SBarry Smith   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
5354b9ad928SBarry Smith   pc->ops->view               = PCView_Composite;
5364b9ad928SBarry Smith   pc->ops->applyrichardson    = 0;
5374b9ad928SBarry Smith 
5384b9ad928SBarry Smith   pc->data               = (void*)jac;
5394b9ad928SBarry Smith   jac->type              = PC_COMPOSITE_ADDITIVE;
5404b9ad928SBarry Smith   jac->work1             = 0;
5414b9ad928SBarry Smith   jac->work2             = 0;
5424b9ad928SBarry Smith   jac->head              = 0;
5434b9ad928SBarry Smith 
5444b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
5454b9ad928SBarry Smith                     PCCompositeSetType_Composite);CHKERRQ(ierr);
5464b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
5474b9ad928SBarry Smith                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
5484b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
5494b9ad928SBarry Smith                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
5504b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
5514b9ad928SBarry Smith                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
5524b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
5534b9ad928SBarry Smith                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
5544b9ad928SBarry Smith 
5554b9ad928SBarry Smith   PetscFunctionReturn(0);
5564b9ad928SBarry Smith }
5574b9ad928SBarry Smith EXTERN_C_END
5584b9ad928SBarry Smith 
559