xref: /petsc/src/ksp/pc/impls/composite/composite.c (revision 2a6744eb01855f5aa328eb8fdf4b0d01e72ad151)
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;
134b9ad928SBarry Smith };
144b9ad928SBarry Smith 
154b9ad928SBarry Smith typedef struct {
164b9ad928SBarry Smith   PC_CompositeLink head;
174b9ad928SBarry Smith   PCCompositeType  type;
184b9ad928SBarry Smith   Vec              work1;
194b9ad928SBarry Smith   Vec              work2;
204b9ad928SBarry Smith   PetscScalar      alpha;
214b9ad928SBarry Smith   PetscTruth       use_true_matrix;
224b9ad928SBarry Smith } PC_Composite;
234b9ad928SBarry Smith 
244b9ad928SBarry Smith #undef __FUNCT__
254b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative"
266849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y)
274b9ad928SBarry Smith {
28dfbe8321SBarry Smith   PetscErrorCode   ierr;
294b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
304b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
314b9ad928SBarry Smith   Mat              mat = pc->pmat;
324b9ad928SBarry Smith 
334b9ad928SBarry Smith   PetscFunctionBegin;
344b9ad928SBarry Smith   if (!next) {
351302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
364b9ad928SBarry Smith   }
374b9ad928SBarry Smith   if (next->next && !jac->work2) { /* allocate second work vector */
384b9ad928SBarry Smith     ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr);
394b9ad928SBarry Smith   }
40d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
414b9ad928SBarry Smith   if (jac->use_true_matrix) mat = pc->mat;
424b9ad928SBarry Smith   while (next->next) {
434b9ad928SBarry Smith     next = next->next;
444b9ad928SBarry Smith     ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr);
45efb30889SBarry Smith     ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr);
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   }
494b9ad928SBarry Smith 
504b9ad928SBarry Smith   PetscFunctionReturn(0);
514b9ad928SBarry Smith }
524b9ad928SBarry Smith 
534b9ad928SBarry Smith /*
544b9ad928SBarry Smith     This is very special for a matrix of the form alpha I + R + S
554b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from
564b9ad928SBarry Smith alpha I + R
574b9ad928SBarry Smith */
584b9ad928SBarry Smith #undef __FUNCT__
594b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special"
606849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y)
614b9ad928SBarry Smith {
62dfbe8321SBarry Smith   PetscErrorCode   ierr;
634b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
644b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
654b9ad928SBarry Smith 
664b9ad928SBarry Smith   PetscFunctionBegin;
674b9ad928SBarry Smith   if (!next) {
681302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
694b9ad928SBarry Smith   }
704b9ad928SBarry Smith   if (!next->next || next->next->next) {
711302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs");
724b9ad928SBarry Smith   }
734b9ad928SBarry Smith 
74d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
75d8fd42c4SBarry Smith   ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr);
764b9ad928SBarry Smith   PetscFunctionReturn(0);
774b9ad928SBarry Smith }
784b9ad928SBarry Smith 
794b9ad928SBarry Smith #undef __FUNCT__
804b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive"
816849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y)
824b9ad928SBarry Smith {
83dfbe8321SBarry Smith   PetscErrorCode   ierr;
844b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
854b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
864b9ad928SBarry Smith 
874b9ad928SBarry Smith   PetscFunctionBegin;
884b9ad928SBarry Smith   if (!next) {
891302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()");
904b9ad928SBarry Smith   }
91d8fd42c4SBarry Smith   ierr = PCApply(next->pc,x,y);CHKERRQ(ierr);
924b9ad928SBarry Smith   while (next->next) {
934b9ad928SBarry Smith     next = next->next;
94d8fd42c4SBarry Smith     ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr);
95efb30889SBarry Smith     ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr);
964b9ad928SBarry Smith   }
974b9ad928SBarry Smith   PetscFunctionReturn(0);
984b9ad928SBarry Smith }
994b9ad928SBarry Smith 
1004b9ad928SBarry Smith #undef __FUNCT__
1014b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite"
1026849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc)
1034b9ad928SBarry Smith {
104dfbe8321SBarry Smith   PetscErrorCode   ierr;
1054b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
1064b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1074b9ad928SBarry Smith 
1084b9ad928SBarry Smith   PetscFunctionBegin;
1094b9ad928SBarry Smith   if (!jac->work1) {
11023ce1328SBarry Smith    ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr);
1114b9ad928SBarry Smith   }
1124b9ad928SBarry Smith   while (next) {
1134b9ad928SBarry Smith     ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1144b9ad928SBarry Smith     next = next->next;
1154b9ad928SBarry Smith   }
1164b9ad928SBarry Smith 
1174b9ad928SBarry Smith   PetscFunctionReturn(0);
1184b9ad928SBarry Smith }
1194b9ad928SBarry Smith 
1204b9ad928SBarry Smith #undef __FUNCT__
1214b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite"
1226849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc)
1234b9ad928SBarry Smith {
1244b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
125dfbe8321SBarry Smith   PetscErrorCode   ierr;
1264b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
1274b9ad928SBarry Smith 
1284b9ad928SBarry Smith   PetscFunctionBegin;
1294b9ad928SBarry Smith   while (next) {
1304b9ad928SBarry Smith     ierr = PCDestroy(next->pc);CHKERRQ(ierr);
1314b9ad928SBarry Smith     next = next->next;
1324b9ad928SBarry Smith   }
1334b9ad928SBarry Smith 
1344b9ad928SBarry Smith   if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);}
1354b9ad928SBarry Smith   if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);}
1364b9ad928SBarry Smith   ierr = PetscFree(jac);CHKERRQ(ierr);
1374b9ad928SBarry Smith   PetscFunctionReturn(0);
1384b9ad928SBarry Smith }
1394b9ad928SBarry Smith 
1404b9ad928SBarry Smith #undef __FUNCT__
1414b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite"
1426849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc)
1434b9ad928SBarry Smith {
1444b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
145dfbe8321SBarry Smith   PetscErrorCode   ierr;
1469dcbbd2bSBarry Smith   PetscInt         nmax = 8,i;
1474b9ad928SBarry Smith   PC_CompositeLink next;
148e5999256SBarry Smith   char             *pcs[8];
1494b9ad928SBarry Smith   PetscTruth       flg;
1504b9ad928SBarry Smith 
1514b9ad928SBarry Smith   PetscFunctionBegin;
1524b9ad928SBarry Smith   ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr);
1539dcbbd2bSBarry Smith   ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr);
1544b9ad928SBarry Smith     ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr);
1554b9ad928SBarry Smith     if (flg) {
1564b9ad928SBarry Smith       ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr);
1574b9ad928SBarry Smith     }
1584b9ad928SBarry Smith     ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr);
1594b9ad928SBarry Smith     if (flg) {
1604b9ad928SBarry Smith       for (i=0; i<nmax; i++) {
1614b9ad928SBarry Smith         ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr);
1624b9ad928SBarry Smith       }
1634b9ad928SBarry Smith     }
1644b9ad928SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
1654b9ad928SBarry Smith 
1664b9ad928SBarry Smith   next = jac->head;
1674b9ad928SBarry Smith   while (next) {
1684b9ad928SBarry Smith     ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr);
1694b9ad928SBarry Smith     next = next->next;
1704b9ad928SBarry Smith   }
1714b9ad928SBarry Smith   PetscFunctionReturn(0);
1724b9ad928SBarry Smith }
1734b9ad928SBarry Smith 
1744b9ad928SBarry Smith #undef __FUNCT__
1754b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite"
1766849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer)
1774b9ad928SBarry Smith {
1784b9ad928SBarry Smith   PC_Composite     *jac = (PC_Composite*)pc->data;
179dfbe8321SBarry Smith   PetscErrorCode   ierr;
1804b9ad928SBarry Smith   PC_CompositeLink next = jac->head;
18132077d6dSBarry Smith   PetscTruth       iascii;
1824b9ad928SBarry Smith 
1834b9ad928SBarry Smith   PetscFunctionBegin;
18432077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
18532077d6dSBarry Smith   if (iascii) {
1869dcbbd2bSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr);
1874b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr);
1884b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
1894b9ad928SBarry Smith   } else {
19079a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name);
1914b9ad928SBarry Smith   }
19232077d6dSBarry Smith   if (iascii) {
1934b9ad928SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
1944b9ad928SBarry Smith   }
1954b9ad928SBarry Smith   while (next) {
1964b9ad928SBarry Smith     ierr = PCView(next->pc,viewer);CHKERRQ(ierr);
1974b9ad928SBarry Smith     next = next->next;
1984b9ad928SBarry Smith   }
19932077d6dSBarry Smith   if (iascii) {
2004b9ad928SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
2014b9ad928SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr);
2024b9ad928SBarry Smith   }
2034b9ad928SBarry Smith   PetscFunctionReturn(0);
2044b9ad928SBarry Smith }
2054b9ad928SBarry Smith 
2064b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/
2074b9ad928SBarry Smith 
2084b9ad928SBarry Smith EXTERN_C_BEGIN
2094b9ad928SBarry Smith #undef __FUNCT__
2104b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite"
211dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha)
2124b9ad928SBarry Smith {
2134b9ad928SBarry Smith   PC_Composite *jac = (PC_Composite*)pc->data;
2144b9ad928SBarry Smith   PetscFunctionBegin;
2154b9ad928SBarry Smith   jac->alpha = alpha;
2164b9ad928SBarry Smith   PetscFunctionReturn(0);
2174b9ad928SBarry Smith }
2184b9ad928SBarry Smith EXTERN_C_END
2194b9ad928SBarry Smith 
2204b9ad928SBarry Smith EXTERN_C_BEGIN
2214b9ad928SBarry Smith #undef __FUNCT__
2224b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite"
223dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type)
2244b9ad928SBarry Smith {
2254b9ad928SBarry Smith   PetscFunctionBegin;
2264b9ad928SBarry Smith   if (type == PC_COMPOSITE_ADDITIVE) {
2274b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Additive;
2284b9ad928SBarry Smith   } else if (type ==  PC_COMPOSITE_MULTIPLICATIVE) {
2294b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Multiplicative;
2304b9ad928SBarry Smith   } else if (type ==  PC_COMPOSITE_SPECIAL) {
2314b9ad928SBarry Smith     pc->ops->apply = PCApply_Composite_Special;
2324b9ad928SBarry Smith   } else {
2331302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type");
2344b9ad928SBarry Smith   }
2354b9ad928SBarry Smith   PetscFunctionReturn(0);
2364b9ad928SBarry Smith }
2374b9ad928SBarry Smith EXTERN_C_END
2384b9ad928SBarry Smith 
2394b9ad928SBarry Smith EXTERN_C_BEGIN
2404b9ad928SBarry Smith #undef __FUNCT__
2414b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite"
242dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type)
2434b9ad928SBarry Smith {
2444b9ad928SBarry Smith   PC_Composite     *jac;
2455a9f2f41SSatish Balay   PC_CompositeLink next,ilink;
246dfbe8321SBarry Smith   PetscErrorCode   ierr;
24779416396SBarry Smith   PetscInt         cnt = 0;
2482dcb1b2aSMatthew Knepley   const char       *prefix;
2492dcb1b2aSMatthew Knepley   char             newprefix[8];
2504b9ad928SBarry Smith 
2514b9ad928SBarry Smith   PetscFunctionBegin;
2525a9f2f41SSatish Balay   ierr       = PetscNew(struct _PC_CompositeLink,&ilink);CHKERRQ(ierr);
2535a9f2f41SSatish Balay   ilink->next = 0;
2545a9f2f41SSatish Balay   ierr = PCCreate(pc->comm,&ilink->pc);CHKERRQ(ierr);
2554b9ad928SBarry Smith 
2564b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2574b9ad928SBarry Smith   next = jac->head;
2584b9ad928SBarry Smith   if (!next) {
2595a9f2f41SSatish Balay     jac->head = ilink;
2604b9ad928SBarry Smith   } else {
2614b9ad928SBarry Smith     cnt++;
2624b9ad928SBarry Smith     while (next->next) {
2634b9ad928SBarry Smith       next = next->next;
2644b9ad928SBarry Smith       cnt++;
2654b9ad928SBarry Smith     }
2665a9f2f41SSatish Balay     next->next = ilink;
2674b9ad928SBarry Smith   }
2684b9ad928SBarry Smith   ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
2695a9f2f41SSatish Balay   ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr);
27013f74950SBarry Smith   sprintf(newprefix,"sub_%d_",(int)cnt);
2715a9f2f41SSatish Balay   ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr);
2724b9ad928SBarry Smith   /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */
2735a9f2f41SSatish Balay   ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr);
2744b9ad928SBarry Smith 
2754b9ad928SBarry Smith   PetscFunctionReturn(0);
2764b9ad928SBarry Smith }
2774b9ad928SBarry Smith EXTERN_C_END
2784b9ad928SBarry Smith 
2794b9ad928SBarry Smith EXTERN_C_BEGIN
2804b9ad928SBarry Smith #undef __FUNCT__
2814b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite"
282dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc)
2834b9ad928SBarry Smith {
2844b9ad928SBarry Smith   PC_Composite     *jac;
2854b9ad928SBarry Smith   PC_CompositeLink next;
28679416396SBarry Smith   PetscInt         i;
2874b9ad928SBarry Smith 
2884b9ad928SBarry Smith   PetscFunctionBegin;
2894b9ad928SBarry Smith   jac  = (PC_Composite*)pc->data;
2904b9ad928SBarry Smith   next = jac->head;
2914b9ad928SBarry Smith   for (i=0; i<n; i++) {
2924b9ad928SBarry Smith     if (!next->next) {
2931302d50aSBarry Smith       SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner");
2944b9ad928SBarry Smith     }
2954b9ad928SBarry Smith     next = next->next;
2964b9ad928SBarry Smith   }
2974b9ad928SBarry Smith   *subpc = next->pc;
2984b9ad928SBarry Smith   PetscFunctionReturn(0);
2994b9ad928SBarry Smith }
3004b9ad928SBarry Smith EXTERN_C_END
3014b9ad928SBarry Smith 
3024b9ad928SBarry Smith EXTERN_C_BEGIN
3034b9ad928SBarry Smith #undef __FUNCT__
3044b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite"
305dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc)
3064b9ad928SBarry Smith {
3074b9ad928SBarry Smith   PC_Composite   *jac;
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith   PetscFunctionBegin;
3104b9ad928SBarry Smith   jac                  = (PC_Composite*)pc->data;
3114b9ad928SBarry Smith   jac->use_true_matrix = PETSC_TRUE;
3124b9ad928SBarry Smith   PetscFunctionReturn(0);
3134b9ad928SBarry Smith }
3144b9ad928SBarry Smith EXTERN_C_END
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */
3174b9ad928SBarry Smith #undef __FUNCT__
3184b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType"
319f39d8e23SSatish Balay /*@
3204b9ad928SBarry Smith    PCCompositeSetType - Sets the type of composite preconditioner.
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Collective on PC
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith    Input Parameter:
325*2a6744ebSBarry Smith +  pc - the preconditioner context
326*2a6744ebSBarry Smith -  type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith    Options Database Key:
3294b9ad928SBarry Smith .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith    Level: Developer
3324b9ad928SBarry Smith 
3334b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3344b9ad928SBarry Smith @*/
335dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type)
3364b9ad928SBarry Smith {
337dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCCompositeType);
3384b9ad928SBarry Smith 
3394b9ad928SBarry Smith   PetscFunctionBegin;
3404482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3414b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr);
3424b9ad928SBarry Smith   if (f) {
3434b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
3444b9ad928SBarry Smith   }
3454b9ad928SBarry Smith   PetscFunctionReturn(0);
3464b9ad928SBarry Smith }
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith #undef __FUNCT__
3494b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha"
350f39d8e23SSatish Balay /*@
3514b9ad928SBarry Smith    PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner
3524b9ad928SBarry Smith      for alphaI + R + S
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith    Collective on PC
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith    Input Parameter:
3574b9ad928SBarry Smith +  pc - the preconditioner context
3584b9ad928SBarry Smith -  alpha - scale on identity
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith    Level: Developer
3614b9ad928SBarry Smith 
3624b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative
3634b9ad928SBarry Smith @*/
364dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha)
3654b9ad928SBarry Smith {
366dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscScalar);
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith   PetscFunctionBegin;
3694482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3704b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr);
3714b9ad928SBarry Smith   if (f) {
3724b9ad928SBarry Smith     ierr = (*f)(pc,alpha);CHKERRQ(ierr);
3734b9ad928SBarry Smith   }
3744b9ad928SBarry Smith   PetscFunctionReturn(0);
3754b9ad928SBarry Smith }
3764b9ad928SBarry Smith 
3774b9ad928SBarry Smith #undef __FUNCT__
3784b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC"
3794b9ad928SBarry Smith /*@C
3804b9ad928SBarry Smith    PCCompositeAddPC - Adds another PC to the composite PC.
3814b9ad928SBarry Smith 
3824b9ad928SBarry Smith    Collective on PC
3834b9ad928SBarry Smith 
3844b9ad928SBarry Smith    Input Parameters:
385*2a6744ebSBarry Smith +  pc - the preconditioner context
386*2a6744ebSBarry Smith -  type - the type of the new preconditioner
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith    Level: Developer
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith .keywords: PC, composite preconditioner, add
3914b9ad928SBarry Smith @*/
392dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type)
3934b9ad928SBarry Smith {
394dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PCType);
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith   PetscFunctionBegin;
3974482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3984b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr);
3994b9ad928SBarry Smith   if (f) {
4004b9ad928SBarry Smith     ierr = (*f)(pc,type);CHKERRQ(ierr);
4014b9ad928SBarry Smith   }
4024b9ad928SBarry Smith   PetscFunctionReturn(0);
4034b9ad928SBarry Smith }
4044b9ad928SBarry Smith 
4054b9ad928SBarry Smith #undef __FUNCT__
4064b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC"
407f39d8e23SSatish Balay /*@
4084b9ad928SBarry Smith    PCCompositeGetPC - Gets one of the PC objects in the composite PC.
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith    Not Collective
4114b9ad928SBarry Smith 
4124b9ad928SBarry Smith    Input Parameter:
413*2a6744ebSBarry Smith +  pc - the preconditioner context
414*2a6744ebSBarry Smith -  n - the number of the pc requested
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith    Output Parameters:
4174b9ad928SBarry Smith .  subpc - the PC requested
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith    Level: Developer
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith .seealso: PCCompositeAddPC()
4244b9ad928SBarry Smith @*/
425dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc)
4264b9ad928SBarry Smith {
42713f74950SBarry Smith   PetscErrorCode ierr,(*f)(PC,PetscInt,PC *);
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith   PetscFunctionBegin;
4304482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4314482741eSBarry Smith   PetscValidPointer(subpc,3);
4324b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4334b9ad928SBarry Smith   if (f) {
4344b9ad928SBarry Smith     ierr = (*f)(pc,n,subpc);CHKERRQ(ierr);
4354b9ad928SBarry Smith   } else {
4361302d50aSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type");
4374b9ad928SBarry Smith   }
4384b9ad928SBarry Smith   PetscFunctionReturn(0);
4394b9ad928SBarry Smith }
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith #undef __FUNCT__
4424b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue"
4434b9ad928SBarry Smith /*@
4444b9ad928SBarry Smith    PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than
4454b9ad928SBarry Smith                       the matrix used to define the preconditioner) is used to compute
4464b9ad928SBarry Smith                       the residual when the multiplicative scheme is used.
4474b9ad928SBarry Smith 
4484b9ad928SBarry Smith    Collective on PC
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith    Input Parameters:
4514b9ad928SBarry Smith .  pc - the preconditioner context
4524b9ad928SBarry Smith 
4534b9ad928SBarry Smith    Options Database Key:
4544b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
4554b9ad928SBarry Smith 
4564b9ad928SBarry Smith    Note:
4574b9ad928SBarry Smith    For the common case in which the preconditioning and linear
4584b9ad928SBarry Smith    system matrices are identical, this routine is unnecessary.
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith    Level: Developer
4614b9ad928SBarry Smith 
4624b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue()
4654b9ad928SBarry Smith @*/
466dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc)
4674b9ad928SBarry Smith {
468dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC);
4694b9ad928SBarry Smith 
4704b9ad928SBarry Smith   PetscFunctionBegin;
4714482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4724b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr);
4734b9ad928SBarry Smith   if (f) {
4744b9ad928SBarry Smith     ierr = (*f)(pc);CHKERRQ(ierr);
4754b9ad928SBarry Smith   }
4764b9ad928SBarry Smith   PetscFunctionReturn(0);
4774b9ad928SBarry Smith }
4784b9ad928SBarry Smith 
4794b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/
4804b9ad928SBarry Smith 
4814b9ad928SBarry Smith /*MC
4824b9ad928SBarry Smith      PCCOMPOSITE - Build a preconditioner by composing together several preconditioners
4834b9ad928SBarry Smith 
4844b9ad928SBarry Smith    Options Database Keys:
4854b9ad928SBarry Smith .  -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type
4864b9ad928SBarry Smith .  -pc_composite_true - Activates PCCompositeSetUseTrue()
4874b9ad928SBarry Smith 
4884b9ad928SBarry Smith    Level: intermediate
4894b9ad928SBarry Smith 
4904b9ad928SBarry Smith    Concepts: composing solvers
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith    Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more
4934b9ad928SBarry Smith           inner PCs to be PCKSP.
4944b9ad928SBarry Smith           Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
49579416396SBarry Smith           the incorrect answer) unless you use KSPFGMRES as the outter Krylov method
4964b9ad928SBarry Smith 
4974b9ad928SBarry Smith 
4984b9ad928SBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
4994b9ad928SBarry Smith            PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(),
5004b9ad928SBarry Smith            PCCompositeGetPC(), PCCompositeSetUseTrue()
5014b9ad928SBarry Smith 
5024b9ad928SBarry Smith M*/
5034b9ad928SBarry Smith 
5044b9ad928SBarry Smith EXTERN_C_BEGIN
5054b9ad928SBarry Smith #undef __FUNCT__
5064b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite"
507dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc)
5084b9ad928SBarry Smith {
509dfbe8321SBarry Smith   PetscErrorCode ierr;
5104b9ad928SBarry Smith   PC_Composite   *jac;
5114b9ad928SBarry Smith 
5124b9ad928SBarry Smith   PetscFunctionBegin;
5134b9ad928SBarry Smith   ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr);
51452e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_Composite));CHKERRQ(ierr);
5154b9ad928SBarry Smith   pc->ops->apply              = PCApply_Composite_Additive;
5164b9ad928SBarry Smith   pc->ops->setup              = PCSetUp_Composite;
5174b9ad928SBarry Smith   pc->ops->destroy            = PCDestroy_Composite;
5184b9ad928SBarry Smith   pc->ops->setfromoptions     = PCSetFromOptions_Composite;
5194b9ad928SBarry Smith   pc->ops->view               = PCView_Composite;
5204b9ad928SBarry Smith   pc->ops->applyrichardson    = 0;
5214b9ad928SBarry Smith 
5224b9ad928SBarry Smith   pc->data               = (void*)jac;
5234b9ad928SBarry Smith   jac->type              = PC_COMPOSITE_ADDITIVE;
5244b9ad928SBarry Smith   jac->work1             = 0;
5254b9ad928SBarry Smith   jac->work2             = 0;
5264b9ad928SBarry Smith   jac->head              = 0;
5274b9ad928SBarry Smith 
5284b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite",
5294b9ad928SBarry Smith                     PCCompositeSetType_Composite);CHKERRQ(ierr);
5304b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite",
5314b9ad928SBarry Smith                     PCCompositeAddPC_Composite);CHKERRQ(ierr);
5324b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite",
5334b9ad928SBarry Smith                     PCCompositeGetPC_Composite);CHKERRQ(ierr);
5344b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite",
5354b9ad928SBarry Smith                     PCCompositeSetUseTrue_Composite);CHKERRQ(ierr);
5364b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite",
5374b9ad928SBarry Smith                     PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr);
5384b9ad928SBarry Smith 
5394b9ad928SBarry Smith   PetscFunctionReturn(0);
5404b9ad928SBarry Smith }
5414b9ad928SBarry Smith EXTERN_C_END
5424b9ad928SBarry Smith 
543