xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 724c2c9915ab311ad7a70ce94136629dd9f613ae)
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;
129*724c2c99SHong Zhang   PC_CompositeLink next = jac->head,next_tmp;
1304b9ad928SBarry Smith 
1314b9ad928SBarry Smith   PetscFunctionBegin;
1324b9ad928SBarry Smith   while (next) {
1334b9ad928SBarry Smith     ierr = PCDestroy(next->pc);CHKERRQ(ierr);
134*724c2c99SHong Zhang     next_tmp = next;
1354b9ad928SBarry Smith     next     = next->next;
136*724c2c99SHong Zhang     ierr = PetscFree(next_tmp);CHKERRQ(ierr);
1374b9ad928SBarry Smith   }
1384b9ad928SBarry Smith 
1394b9ad928SBarry Smith   if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);}
1404b9ad928SBarry Smith   if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);}
1414b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1424b9ad928SBarry Smith   PetscFunctionReturn(0);
1434b9ad928SBarry Smith }
1444b9ad928SBarry Smith 
1454b9ad928SBarry Smith #undef __FUNCT__
1464b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite"
1476849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc)
1484b9ad928SBarry Smith {
1494b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
150dfbe8321SBarry Smith   PetscErrorCode   ierr;
1519dcbbd2bSBarry Smith   PetscInt         nmax = 8,i;
1524b9ad928SBarry Smith   PC_CompositeLink next;
153e5999256SBarry Smith   char             *pcs[8];
1544b9ad928SBarry Smith   PetscTruth       flg;
1554b9ad928SBarry Smith 
1564b9ad928SBarry Smith   PetscFunctionBegin;
1574b9ad928SBarry Smith   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
1589dcbbd2bSBarry Smith     ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
15951f519a2SBarry Smith     if (flg) {
16051f519a2SBarry Smith       ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr);
16151f519a2SBarry Smith     }
16290d69ab7SBarry Smith     flg  = PETSC_FALSE;
16390d69ab7SBarry Smith     ierr = PetscOptionsTruth("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);CHKERRQ(ierr);
1644b9ad928SBarry Smith     if (flg) {
1654b9ad928SBarry Smith       ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr);
1664b9ad928SBarry Smith     }
1674b9ad928SBarry Smith     ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr);
1684b9ad928SBarry Smith     if (flg) {
1694b9ad928SBarry Smith       for (i=0; i<nmax; i++) {
1704b9ad928SBarry Smith         ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr);
171*724c2c99SHong Zhang         ierr = PetscFree(pcs[i]);CHKERRQ(ierr); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */
1724b9ad928SBarry Smith       }
1734b9ad928SBarry Smith     }
1744b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith   next = jac->head;
1774b9ad928SBarry Smith   while (next) {
1784b9ad928SBarry Smith     ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr);
1794b9ad928SBarry Smith     next = next->next;
1804b9ad928SBarry Smith   }
1814b9ad928SBarry Smith   PetscFunctionReturn(0);
1824b9ad928SBarry Smith }
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith #undef __FUNCT__
1854b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite"
1866849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
1874b9ad928SBarry Smith {
1884b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
189dfbe8321SBarry Smith   PetscErrorCode   ierr;
1904b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
19132077d6dSBarry Smith   PetscTruth       iascii;
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith   PetscFunctionBegin;
1942692d6eeSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
19532077d6dSBarry Smith   if (iascii) {
1969dcbbd2bSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr);
1974b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr);
1984b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
1994b9ad928SBarry Smith   } else {
20065e19b50SBarry Smith     SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
2014b9ad928SBarry Smith   }
20232077d6dSBarry Smith   if (iascii) {
2034b9ad928SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
2044b9ad928SBarry Smith   }
2054b9ad928SBarry Smith   while (next) {
2064b9ad928SBarry Smith     ierr = PCView(next->pc,viewer);CHKERRQ(ierr);
2074b9ad928SBarry Smith     next = next->next;
2084b9ad928SBarry Smith   }
20932077d6dSBarry Smith   if (iascii) {
2104b9ad928SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2114b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
2124b9ad928SBarry Smith   }
2134b9ad928SBarry Smith   PetscFunctionReturn(0);
2144b9ad928SBarry Smith }
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
2174b9ad928SBarry Smith 
2184b9ad928SBarry Smith EXTERN_C_BEGIN
2194b9ad928SBarry Smith #undef __FUNCT__
2204b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite"
221dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
2224b9ad928SBarry Smith {
2234b9ad928SBarry Smith   PC_Composite *jac = (PC_Composite*)pc->data;
2244b9ad928SBarry Smith   PetscFunctionBegin;
2254b9ad928SBarry Smith   jac->alpha = alpha;
2264b9ad928SBarry Smith   PetscFunctionReturn(0);
2274b9ad928SBarry Smith }
2284b9ad928SBarry Smith EXTERN_C_END
2294b9ad928SBarry Smith 
2304b9ad928SBarry Smith EXTERN_C_BEGIN
2314b9ad928SBarry Smith #undef __FUNCT__
2324b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite"
233dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type)
2344b9ad928SBarry Smith {
2354b9ad928SBarry Smith   PetscFunctionBegin;
2364b9ad928SBarry Smith   if (type == PC_COMPOSITE_ADDITIVE) {
2374b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Additive;
238421e10b8SBarry Smith   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) {
2394b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Multiplicative;
2404b9ad928SBarry Smith   } else if (type ==  PC_COMPOSITE_SPECIAL) {
2414b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Special;
242e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
2434b9ad928SBarry Smith   PetscFunctionReturn(0);
2444b9ad928SBarry Smith }
2454b9ad928SBarry Smith EXTERN_C_END
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith EXTERN_C_BEGIN
2484b9ad928SBarry Smith #undef __FUNCT__
2494b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite"
250dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type)
2514b9ad928SBarry Smith {
2524b9ad928SBarry Smith   PC_Composite     *jac;
2535a9f2f41SSatish Balay   PC_CompositeLink next,ilink;
254dfbe8321SBarry Smith   PetscErrorCode   ierr;
25579416396SBarry Smith   PetscInt         cnt = 0;
2562dcb1b2aSMatthew Knepley   const char       *prefix;
2572dcb1b2aSMatthew Knepley   char             newprefix[8];
2584b9ad928SBarry Smith 
2594b9ad928SBarry Smith   PetscFunctionBegin;
26038f2d2fdSLisandro Dalcin   ierr        = PetscNewLog(pc,struct _PC_CompositeLink,&ilink);CHKERRQ(ierr);
2615a9f2f41SSatish Balay   ilink->next = 0;
2627adad957SLisandro Dalcin   ierr = PCCreate(((PetscObject)pc)->comm,&ilink->pc);CHKERRQ(ierr);
2634b9ad928SBarry Smith 
2644b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2654b9ad928SBarry Smith   next = jac->head;
2664b9ad928SBarry Smith   if (!next) {
2675a9f2f41SSatish Balay     jac->head       = ilink;
268421e10b8SBarry Smith     ilink->previous = PETSC_NULL;
2694b9ad928SBarry Smith   } else {
2704b9ad928SBarry Smith     cnt++;
2714b9ad928SBarry Smith     while (next->next) {
2724b9ad928SBarry Smith       next = next->next;
2734b9ad928SBarry Smith       cnt++;
2744b9ad928SBarry Smith     }
2755a9f2f41SSatish Balay     next->next      = ilink;
276421e10b8SBarry Smith     ilink->previous = next;
2774b9ad928SBarry Smith   }
2784b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2795a9f2f41SSatish Balay   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
28013f74950SBarry Smith   sprintf(newprefix,"sub_%d_",(int)cnt);
2815a9f2f41SSatish Balay   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
2824b9ad928SBarry Smith   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
2835a9f2f41SSatish Balay   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
2844b9ad928SBarry Smith   PetscFunctionReturn(0);
2854b9ad928SBarry Smith }
2864b9ad928SBarry Smith EXTERN_C_END
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith EXTERN_C_BEGIN
2894b9ad928SBarry Smith #undef __FUNCT__
2904b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite"
291dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
2924b9ad928SBarry Smith {
2934b9ad928SBarry Smith   PC_Composite     *jac;
2944b9ad928SBarry Smith   PC_CompositeLink next;
29579416396SBarry Smith   PetscInt         i;
2964b9ad928SBarry Smith 
2974b9ad928SBarry Smith   PetscFunctionBegin;
2984b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2994b9ad928SBarry Smith   next = jac->head;
3004b9ad928SBarry Smith   for (i=0; i<n; i++) {
301e7e72b3dSBarry Smith     if (!next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
3024b9ad928SBarry Smith     next = next->next;
3034b9ad928SBarry Smith   }
3044b9ad928SBarry Smith   *subpc = next->pc;
3054b9ad928SBarry Smith   PetscFunctionReturn(0);
3064b9ad928SBarry Smith }
3074b9ad928SBarry Smith EXTERN_C_END
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith EXTERN_C_BEGIN
3104b9ad928SBarry Smith #undef __FUNCT__
3114b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
312dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
3134b9ad928SBarry Smith {
3144b9ad928SBarry Smith   PC_Composite   *jac;
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith   PetscFunctionBegin;
3174b9ad928SBarry Smith   jac                  = (PC_Composite*)pc->data;
3184b9ad928SBarry Smith   jac->use_true_matrix = PETSC_TRUE;
3194b9ad928SBarry Smith   PetscFunctionReturn(0);
3204b9ad928SBarry Smith }
3214b9ad928SBarry Smith EXTERN_C_END
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */
3244b9ad928SBarry Smith #undef __FUNCT__
3254b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType"
326f39d8e23SSatish Balay /*@
3274b9ad928SBarry Smith    PCCompositeSetType - Sets the type of composite preconditioner.
3284b9ad928SBarry Smith 
329ad4df100SBarry Smith    Logically Collective on PC
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith    Input Parameter:
3322a6744ebSBarry Smith +  pc - the preconditioner context
3332a6744ebSBarry Smith -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
3344b9ad928SBarry Smith 
3354b9ad928SBarry Smith    Options Database Key:
3364b9ad928SBarry Smith .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith    Level: Developer
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3414b9ad928SBarry Smith @*/
342dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
3434b9ad928SBarry Smith {
344dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
3454b9ad928SBarry Smith 
3464b9ad928SBarry Smith   PetscFunctionBegin;
3470700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
348c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc,type,2);
349c5eb9154SBarry Smith 
3504b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
3514b9ad928SBarry Smith   if (f) {
3524b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
3534b9ad928SBarry Smith   }
3544b9ad928SBarry Smith   PetscFunctionReturn(0);
3554b9ad928SBarry Smith }
3564b9ad928SBarry Smith 
3574b9ad928SBarry Smith #undef __FUNCT__
3584b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha"
359f39d8e23SSatish Balay /*@
3604b9ad928SBarry Smith    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
3614b9ad928SBarry Smith      for alphaI + R + S
3624b9ad928SBarry Smith 
363ad4df100SBarry Smith    Logically Collective on PC
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Input Parameter:
3664b9ad928SBarry Smith +  pc - the preconditioner context
3674b9ad928SBarry Smith -  alpha - scale on identity
3684b9ad928SBarry Smith 
3694b9ad928SBarry Smith    Level: Developer
3704b9ad928SBarry Smith 
3714b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3724b9ad928SBarry Smith @*/
373dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
3744b9ad928SBarry Smith {
375dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscScalar);
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith   PetscFunctionBegin;
3780700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
379c5eb9154SBarry Smith   PetscValidLogicalCollectiveScalar(pc,alpha,2);
380c5eb9154SBarry Smith 
3814b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
3824b9ad928SBarry Smith   if (f) {
3834b9ad928SBarry Smith     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
3844b9ad928SBarry Smith   }
3854b9ad928SBarry Smith   PetscFunctionReturn(0);
3864b9ad928SBarry Smith }
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith #undef __FUNCT__
3894b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC"
3904b9ad928SBarry Smith /*@C
3914b9ad928SBarry Smith    PCCompositeAddPC - Adds another PC to the composite PC.
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Collective on PC
3944b9ad928SBarry Smith 
3954b9ad928SBarry Smith    Input Parameters:
3962a6744ebSBarry Smith +  pc - the preconditioner context
3972a6744ebSBarry Smith -  type - the type of the new preconditioner
3984b9ad928SBarry Smith 
3994b9ad928SBarry Smith    Level: Developer
4004b9ad928SBarry Smith 
4014b9ad928SBarry Smith .keywords: PC, composite preconditioner, add
4024b9ad928SBarry Smith @*/
403dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
4044b9ad928SBarry Smith {
405dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCType);
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith   PetscFunctionBegin;
4080700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4094b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4104b9ad928SBarry Smith   if (f) {
4114b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
4124b9ad928SBarry Smith   }
4134b9ad928SBarry Smith   PetscFunctionReturn(0);
4144b9ad928SBarry Smith }
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith #undef __FUNCT__
4174b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC"
418f39d8e23SSatish Balay /*@
4194b9ad928SBarry Smith    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith    Not Collective
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith    Input Parameter:
4242a6744ebSBarry Smith +  pc - the preconditioner context
4252a6744ebSBarry Smith -  n - the number of the pc requested
4264b9ad928SBarry Smith 
4274b9ad928SBarry Smith    Output Parameters:
4284b9ad928SBarry Smith .  subpc - the PC requested
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith    Level: Developer
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith .seealso: PCCompositeAddPC()
4354b9ad928SBarry Smith @*/
436dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
4374b9ad928SBarry Smith {
43813f74950SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
4394b9ad928SBarry Smith 
4404b9ad928SBarry Smith   PetscFunctionBegin;
4410700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4424482741eSBarry Smith   PetscValidPointer(subpc,3);
4434b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4444b9ad928SBarry Smith   if (f) {
4454b9ad928SBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
446e7e72b3dSBarry Smith   } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
4474b9ad928SBarry Smith   PetscFunctionReturn(0);
4484b9ad928SBarry Smith }
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith #undef __FUNCT__
4514b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue"
4524b9ad928SBarry Smith /*@
4534b9ad928SBarry Smith    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
4544b9ad928SBarry Smith                       the matrix used to define the preconditioner) is used to compute
4554b9ad928SBarry Smith                       the residual when the multiplicative scheme is used.
4564b9ad928SBarry Smith 
457ad4df100SBarry Smith    Logically Collective on PC
4584b9ad928SBarry Smith 
4594b9ad928SBarry Smith    Input Parameters:
4604b9ad928SBarry Smith .  pc - the preconditioner context
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith    Options Database Key:
4634b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
4644b9ad928SBarry Smith 
4654b9ad928SBarry Smith    Note:
4664b9ad928SBarry Smith    For the common case in which the preconditioning and linear
4674b9ad928SBarry Smith    system matrices are identical, this routine is unnecessary.
4684b9ad928SBarry Smith 
4694b9ad928SBarry Smith    Level: Developer
4704b9ad928SBarry Smith 
4714b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag
4724b9ad928SBarry Smith 
4734b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
4744b9ad928SBarry Smith @*/
475dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
4764b9ad928SBarry Smith {
477dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC);
4784b9ad928SBarry Smith 
4794b9ad928SBarry Smith   PetscFunctionBegin;
4800700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4814b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
4824b9ad928SBarry Smith   if (f) {
4834b9ad928SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
4844b9ad928SBarry Smith   }
4854b9ad928SBarry Smith   PetscFunctionReturn(0);
4864b9ad928SBarry Smith }
4874b9ad928SBarry Smith 
4884b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith /*MC
4914b9ad928SBarry Smith      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
4924b9ad928SBarry Smith 
4934b9ad928SBarry Smith    Options Database Keys:
49451f519a2SBarry Smith +  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
4954b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
49651f519a2SBarry Smith -  -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose
4974b9ad928SBarry Smith 
4984b9ad928SBarry Smith    Level: intermediate
4994b9ad928SBarry Smith 
5004b9ad928SBarry Smith    Concepts: composing solvers
5014b9ad928SBarry Smith 
5024b9ad928SBarry Smith    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
5034b9ad928SBarry Smith           inner PCs to be PCKSP.
5044b9ad928SBarry Smith           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
50579416396SBarry Smith           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
5064b9ad928SBarry Smith 
5074b9ad928SBarry Smith 
5084b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
5094b9ad928SBarry Smith            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
5104b9ad928SBarry Smith            PCCompositeGetPC(), PCCompositeSetUseTrue()
5114b9ad928SBarry Smith 
5124b9ad928SBarry Smith M*/
5134b9ad928SBarry Smith 
5144b9ad928SBarry Smith EXTERN_C_BEGIN
5154b9ad928SBarry Smith #undef __FUNCT__
5164b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite"
517dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
5184b9ad928SBarry Smith {
519dfbe8321SBarry Smith   PetscErrorCode ierr;
5204b9ad928SBarry Smith   PC_Composite   *jac;
5214b9ad928SBarry Smith 
5224b9ad928SBarry Smith   PetscFunctionBegin;
52338f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr);
5244b9ad928SBarry Smith   pc->ops->apply              = PCApply_Composite_Additive;
5254b9ad928SBarry Smith   pc->ops->setup              = PCSetUp_Composite;
5264b9ad928SBarry Smith   pc->ops->destroy            = PCDestroy_Composite;
5274b9ad928SBarry Smith   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
5284b9ad928SBarry Smith   pc->ops->view               = PCView_Composite;
5294b9ad928SBarry Smith   pc->ops->applyrichardson    = 0;
5304b9ad928SBarry Smith 
5314b9ad928SBarry Smith   pc->data               = (void*)jac;
5324b9ad928SBarry Smith   jac->type              = PC_COMPOSITE_ADDITIVE;
5334b9ad928SBarry Smith   jac->work1             = 0;
5344b9ad928SBarry Smith   jac->work2             = 0;
5354b9ad928SBarry Smith   jac->head              = 0;
5364b9ad928SBarry Smith 
5374b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
5384b9ad928SBarry Smith                     PCCompositeSetType_Composite);CHKERRQ(ierr);
5394b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
5404b9ad928SBarry Smith                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
5414b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
5424b9ad928SBarry Smith                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
5434b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
5444b9ad928SBarry Smith                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
5454b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
5464b9ad928SBarry Smith                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
5474b9ad928SBarry Smith 
5484b9ad928SBarry Smith   PetscFunctionReturn(0);
5494b9ad928SBarry Smith }
5504b9ad928SBarry Smith EXTERN_C_END
5514b9ad928SBarry Smith 
552