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; 22ace3abfcSBarry Smith PetscBool 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__ 124*69d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Composite" 125*69d2c0f9SBarry Smith static PetscErrorCode PCReset_Composite(PC pc) 126*69d2c0f9SBarry Smith { 127*69d2c0f9SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 128*69d2c0f9SBarry Smith PetscErrorCode ierr; 129*69d2c0f9SBarry Smith PC_CompositeLink next = jac->head,next_tmp; 130*69d2c0f9SBarry Smith 131*69d2c0f9SBarry Smith PetscFunctionBegin; 132*69d2c0f9SBarry Smith while (next) { 133*69d2c0f9SBarry Smith ierr = PCReset(next->pc);CHKERRQ(ierr); 134*69d2c0f9SBarry Smith next_tmp = next; 135*69d2c0f9SBarry Smith next = next->next; 136*69d2c0f9SBarry Smith } 137*69d2c0f9SBarry Smith if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 138*69d2c0f9SBarry Smith if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 139*69d2c0f9SBarry Smith PetscFunctionReturn(0); 140*69d2c0f9SBarry Smith } 141*69d2c0f9SBarry Smith 142*69d2c0f9SBarry Smith #undef __FUNCT__ 1434b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 1446849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc) 1454b9ad928SBarry Smith { 1464b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 147dfbe8321SBarry Smith PetscErrorCode ierr; 148724c2c99SHong Zhang PC_CompositeLink next = jac->head,next_tmp; 1494b9ad928SBarry Smith 1504b9ad928SBarry Smith PetscFunctionBegin; 151*69d2c0f9SBarry Smith ierr = PCReset_Composite(pc);CHKERRQ(ierr); 1524b9ad928SBarry Smith while (next) { 1534b9ad928SBarry Smith ierr = PCDestroy(next->pc);CHKERRQ(ierr); 154724c2c99SHong Zhang next_tmp = next; 1554b9ad928SBarry Smith next = next->next; 156724c2c99SHong Zhang ierr = PetscFree(next_tmp);CHKERRQ(ierr); 1574b9ad928SBarry Smith } 158c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1594b9ad928SBarry Smith PetscFunctionReturn(0); 1604b9ad928SBarry Smith } 1614b9ad928SBarry Smith 1624b9ad928SBarry Smith #undef __FUNCT__ 1634b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 1646849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc) 1654b9ad928SBarry Smith { 1664b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 167dfbe8321SBarry Smith PetscErrorCode ierr; 1689dcbbd2bSBarry Smith PetscInt nmax = 8,i; 1694b9ad928SBarry Smith PC_CompositeLink next; 170e5999256SBarry Smith char *pcs[8]; 171ace3abfcSBarry Smith PetscBool flg; 1724b9ad928SBarry Smith 1734b9ad928SBarry Smith PetscFunctionBegin; 1744b9ad928SBarry Smith ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 1759dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 17651f519a2SBarry Smith if (flg) { 17751f519a2SBarry Smith ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr); 17851f519a2SBarry Smith } 17990d69ab7SBarry Smith flg = PETSC_FALSE; 180acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 1814b9ad928SBarry Smith if (flg) { 1824b9ad928SBarry Smith ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 1854b9ad928SBarry Smith if (flg) { 1864b9ad928SBarry Smith for (i=0; i<nmax; i++) { 1874b9ad928SBarry Smith ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 188724c2c99SHong Zhang ierr = PetscFree(pcs[i]);CHKERRQ(ierr); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 1894b9ad928SBarry Smith } 1904b9ad928SBarry Smith } 1914b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith next = jac->head; 1944b9ad928SBarry Smith while (next) { 1954b9ad928SBarry Smith ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 1964b9ad928SBarry Smith next = next->next; 1974b9ad928SBarry Smith } 1984b9ad928SBarry Smith PetscFunctionReturn(0); 1994b9ad928SBarry Smith } 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith #undef __FUNCT__ 2024b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite" 2036849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 2044b9ad928SBarry Smith { 2054b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 206dfbe8321SBarry Smith PetscErrorCode ierr; 2074b9ad928SBarry Smith PC_CompositeLink next = jac->head; 208ace3abfcSBarry Smith PetscBool iascii; 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith PetscFunctionBegin; 2112692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 21232077d6dSBarry Smith if (iascii) { 2139dcbbd2bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 2144b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 2154b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2164b9ad928SBarry Smith } else { 21765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 2184b9ad928SBarry Smith } 21932077d6dSBarry Smith if (iascii) { 2204b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2214b9ad928SBarry Smith } 2224b9ad928SBarry Smith while (next) { 2234b9ad928SBarry Smith ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 2244b9ad928SBarry Smith next = next->next; 2254b9ad928SBarry Smith } 22632077d6dSBarry Smith if (iascii) { 2274b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2284b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2294b9ad928SBarry Smith } 2304b9ad928SBarry Smith PetscFunctionReturn(0); 2314b9ad928SBarry Smith } 2324b9ad928SBarry Smith 2334b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith EXTERN_C_BEGIN 2364b9ad928SBarry Smith #undef __FUNCT__ 2374b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 2387087cfbeSBarry Smith PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2394b9ad928SBarry Smith { 2404b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2414b9ad928SBarry Smith PetscFunctionBegin; 2424b9ad928SBarry Smith jac->alpha = alpha; 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__ "PCCompositeSetType_Composite" 2507087cfbeSBarry Smith PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 2514b9ad928SBarry Smith { 2524b9ad928SBarry Smith PetscFunctionBegin; 2534b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 2544b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 255421e10b8SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 2564b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 2574b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 2584b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 259e7e72b3dSBarry Smith } else SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 2604b9ad928SBarry Smith PetscFunctionReturn(0); 2614b9ad928SBarry Smith } 2624b9ad928SBarry Smith EXTERN_C_END 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith EXTERN_C_BEGIN 2654b9ad928SBarry Smith #undef __FUNCT__ 2664b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite" 2677087cfbeSBarry Smith PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type) 2684b9ad928SBarry Smith { 2694b9ad928SBarry Smith PC_Composite *jac; 2705a9f2f41SSatish Balay PC_CompositeLink next,ilink; 271dfbe8321SBarry Smith PetscErrorCode ierr; 27279416396SBarry Smith PetscInt cnt = 0; 2732dcb1b2aSMatthew Knepley const char *prefix; 2742dcb1b2aSMatthew Knepley char newprefix[8]; 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith PetscFunctionBegin; 27738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,struct _PC_CompositeLink,&ilink);CHKERRQ(ierr); 2785a9f2f41SSatish Balay ilink->next = 0; 2797adad957SLisandro Dalcin ierr = PCCreate(((PetscObject)pc)->comm,&ilink->pc);CHKERRQ(ierr); 280ace3abfcSBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);CHKERRQ(ierr); 2814b9ad928SBarry Smith 2824b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2834b9ad928SBarry Smith next = jac->head; 2844b9ad928SBarry Smith if (!next) { 2855a9f2f41SSatish Balay jac->head = ilink; 286421e10b8SBarry Smith ilink->previous = PETSC_NULL; 2874b9ad928SBarry Smith } else { 2884b9ad928SBarry Smith cnt++; 2894b9ad928SBarry Smith while (next->next) { 2904b9ad928SBarry Smith next = next->next; 2914b9ad928SBarry Smith cnt++; 2924b9ad928SBarry Smith } 2935a9f2f41SSatish Balay next->next = ilink; 294421e10b8SBarry Smith ilink->previous = next; 2954b9ad928SBarry Smith } 2964b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2975a9f2f41SSatish Balay ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 29813f74950SBarry Smith sprintf(newprefix,"sub_%d_",(int)cnt); 2995a9f2f41SSatish Balay ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 3004b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 3015a9f2f41SSatish Balay ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 3024b9ad928SBarry Smith PetscFunctionReturn(0); 3034b9ad928SBarry Smith } 3044b9ad928SBarry Smith EXTERN_C_END 3054b9ad928SBarry Smith 3064b9ad928SBarry Smith EXTERN_C_BEGIN 3074b9ad928SBarry Smith #undef __FUNCT__ 3084b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite" 3097087cfbeSBarry Smith PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 3104b9ad928SBarry Smith { 3114b9ad928SBarry Smith PC_Composite *jac; 3124b9ad928SBarry Smith PC_CompositeLink next; 31379416396SBarry Smith PetscInt i; 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith PetscFunctionBegin; 3164b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3174b9ad928SBarry Smith next = jac->head; 3184b9ad928SBarry Smith for (i=0; i<n; i++) { 319e7e72b3dSBarry Smith if (!next->next) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 3204b9ad928SBarry Smith next = next->next; 3214b9ad928SBarry Smith } 3224b9ad928SBarry Smith *subpc = next->pc; 3234b9ad928SBarry Smith PetscFunctionReturn(0); 3244b9ad928SBarry Smith } 3254b9ad928SBarry Smith EXTERN_C_END 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith EXTERN_C_BEGIN 3284b9ad928SBarry Smith #undef __FUNCT__ 3294b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite" 3307087cfbeSBarry Smith PetscErrorCode PCCompositeSetUseTrue_Composite(PC pc) 3314b9ad928SBarry Smith { 3324b9ad928SBarry Smith PC_Composite *jac; 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith PetscFunctionBegin; 3354b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3364b9ad928SBarry Smith jac->use_true_matrix = PETSC_TRUE; 3374b9ad928SBarry Smith PetscFunctionReturn(0); 3384b9ad928SBarry Smith } 3394b9ad928SBarry Smith EXTERN_C_END 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 3424b9ad928SBarry Smith #undef __FUNCT__ 3434b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType" 344f39d8e23SSatish Balay /*@ 3454b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 3464b9ad928SBarry Smith 347ad4df100SBarry Smith Logically Collective on PC 3484b9ad928SBarry Smith 3494b9ad928SBarry Smith Input Parameter: 3502a6744ebSBarry Smith + pc - the preconditioner context 3512a6744ebSBarry Smith - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith Options Database Key: 3544b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 3554b9ad928SBarry Smith 3564b9ad928SBarry Smith Level: Developer 3574b9ad928SBarry Smith 3584b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3594b9ad928SBarry Smith @*/ 3607087cfbeSBarry Smith PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 3614b9ad928SBarry Smith { 3627bb14e67SBarry Smith PetscErrorCode ierr; 3634b9ad928SBarry Smith 3644b9ad928SBarry Smith PetscFunctionBegin; 3650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 366c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 3677bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 3684b9ad928SBarry Smith PetscFunctionReturn(0); 3694b9ad928SBarry Smith } 3704b9ad928SBarry Smith 3714b9ad928SBarry Smith #undef __FUNCT__ 3724b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha" 373f39d8e23SSatish Balay /*@ 3744b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 3754b9ad928SBarry Smith for alphaI + R + S 3764b9ad928SBarry Smith 377ad4df100SBarry Smith Logically Collective on PC 3784b9ad928SBarry Smith 3794b9ad928SBarry Smith Input Parameter: 3804b9ad928SBarry Smith + pc - the preconditioner context 3814b9ad928SBarry Smith - alpha - scale on identity 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Level: Developer 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3864b9ad928SBarry Smith @*/ 3877087cfbeSBarry Smith PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 3884b9ad928SBarry Smith { 3894ac538c5SBarry Smith PetscErrorCode ierr; 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith PetscFunctionBegin; 3920700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 393c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(pc,alpha,2); 3944ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha));CHKERRQ(ierr); 3954b9ad928SBarry Smith PetscFunctionReturn(0); 3964b9ad928SBarry Smith } 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith #undef __FUNCT__ 3994b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC" 4004b9ad928SBarry Smith /*@C 4014b9ad928SBarry Smith PCCompositeAddPC - Adds another PC to the composite PC. 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith Collective on PC 4044b9ad928SBarry Smith 4054b9ad928SBarry Smith Input Parameters: 4062a6744ebSBarry Smith + pc - the preconditioner context 4072a6744ebSBarry Smith - type - the type of the new preconditioner 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith Level: Developer 4104b9ad928SBarry Smith 4114b9ad928SBarry Smith .keywords: PC, composite preconditioner, add 4124b9ad928SBarry Smith @*/ 4137087cfbeSBarry Smith PetscErrorCode PCCompositeAddPC(PC pc,PCType type) 4144b9ad928SBarry Smith { 4154ac538c5SBarry Smith PetscErrorCode ierr; 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith PetscFunctionBegin; 4180700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4194ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 4204b9ad928SBarry Smith PetscFunctionReturn(0); 4214b9ad928SBarry Smith } 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith #undef __FUNCT__ 4244b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC" 425f39d8e23SSatish Balay /*@ 4264b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 4274b9ad928SBarry Smith 4284b9ad928SBarry Smith Not Collective 4294b9ad928SBarry Smith 4304b9ad928SBarry Smith Input Parameter: 4312a6744ebSBarry Smith + pc - the preconditioner context 4322a6744ebSBarry Smith - n - the number of the pc requested 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Output Parameters: 4354b9ad928SBarry Smith . subpc - the PC requested 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith Level: Developer 4384b9ad928SBarry Smith 4394b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner 4404b9ad928SBarry Smith 4414b9ad928SBarry Smith .seealso: PCCompositeAddPC() 4424b9ad928SBarry Smith @*/ 4437087cfbeSBarry Smith PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 4444b9ad928SBarry Smith { 4454ac538c5SBarry Smith PetscErrorCode ierr; 4464b9ad928SBarry Smith 4474b9ad928SBarry Smith PetscFunctionBegin; 4480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4494482741eSBarry Smith PetscValidPointer(subpc,3); 4504ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC *),(pc,n,subpc));CHKERRQ(ierr); 4514b9ad928SBarry Smith PetscFunctionReturn(0); 4524b9ad928SBarry Smith } 4534b9ad928SBarry Smith 4544b9ad928SBarry Smith #undef __FUNCT__ 4554b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue" 4564b9ad928SBarry Smith /*@ 4574b9ad928SBarry Smith PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 4584b9ad928SBarry Smith the matrix used to define the preconditioner) is used to compute 4594b9ad928SBarry Smith the residual when the multiplicative scheme is used. 4604b9ad928SBarry Smith 461ad4df100SBarry Smith Logically Collective on PC 4624b9ad928SBarry Smith 4634b9ad928SBarry Smith Input Parameters: 4644b9ad928SBarry Smith . pc - the preconditioner context 4654b9ad928SBarry Smith 4664b9ad928SBarry Smith Options Database Key: 4674b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 4684b9ad928SBarry Smith 4694b9ad928SBarry Smith Note: 4704b9ad928SBarry Smith For the common case in which the preconditioning and linear 4714b9ad928SBarry Smith system matrices are identical, this routine is unnecessary. 4724b9ad928SBarry Smith 4734b9ad928SBarry Smith Level: Developer 4744b9ad928SBarry Smith 4754b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag 4764b9ad928SBarry Smith 4774b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 4784b9ad928SBarry Smith @*/ 4797087cfbeSBarry Smith PetscErrorCode PCCompositeSetUseTrue(PC pc) 4804b9ad928SBarry Smith { 4814ac538c5SBarry Smith PetscErrorCode ierr; 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith PetscFunctionBegin; 4840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4854ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeSetUseTrue_C",(PC),(pc));CHKERRQ(ierr); 4864b9ad928SBarry Smith PetscFunctionReturn(0); 4874b9ad928SBarry Smith } 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 4904b9ad928SBarry Smith 4914b9ad928SBarry Smith /*MC 4924b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 4934b9ad928SBarry Smith 4944b9ad928SBarry Smith Options Database Keys: 49551f519a2SBarry Smith + -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 4964b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 49751f519a2SBarry Smith - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 4984b9ad928SBarry Smith 4994b9ad928SBarry Smith Level: intermediate 5004b9ad928SBarry Smith 5014b9ad928SBarry Smith Concepts: composing solvers 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 5044b9ad928SBarry Smith inner PCs to be PCKSP. 5054b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 50679416396SBarry Smith the incorrect answer) unless you use KSPFGMRES as the outter Krylov method 5074b9ad928SBarry Smith 5084b9ad928SBarry Smith 5094b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 5104b9ad928SBarry Smith PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 5114b9ad928SBarry Smith PCCompositeGetPC(), PCCompositeSetUseTrue() 5124b9ad928SBarry Smith 5134b9ad928SBarry Smith M*/ 5144b9ad928SBarry Smith 5154b9ad928SBarry Smith EXTERN_C_BEGIN 5164b9ad928SBarry Smith #undef __FUNCT__ 5174b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite" 5187087cfbeSBarry Smith PetscErrorCode PCCreate_Composite(PC pc) 5194b9ad928SBarry Smith { 520dfbe8321SBarry Smith PetscErrorCode ierr; 5214b9ad928SBarry Smith PC_Composite *jac; 5224b9ad928SBarry Smith 5234b9ad928SBarry Smith PetscFunctionBegin; 52438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Composite,&jac);CHKERRQ(ierr); 5254b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 5264b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 527*69d2c0f9SBarry Smith pc->ops->reset = PCReset_Composite; 5284b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 5294b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 5304b9ad928SBarry Smith pc->ops->view = PCView_Composite; 5314b9ad928SBarry Smith pc->ops->applyrichardson = 0; 5324b9ad928SBarry Smith 5334b9ad928SBarry Smith pc->data = (void*)jac; 5344b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5354b9ad928SBarry Smith jac->work1 = 0; 5364b9ad928SBarry Smith jac->work2 = 0; 5374b9ad928SBarry Smith jac->head = 0; 5384b9ad928SBarry Smith 5394b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 5404b9ad928SBarry Smith PCCompositeSetType_Composite);CHKERRQ(ierr); 5414b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 5424b9ad928SBarry Smith PCCompositeAddPC_Composite);CHKERRQ(ierr); 5434b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 5444b9ad928SBarry Smith PCCompositeGetPC_Composite);CHKERRQ(ierr); 5454b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 5464b9ad928SBarry Smith PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 5474b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 5484b9ad928SBarry Smith PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith PetscFunctionReturn(0); 5514b9ad928SBarry Smith } 5524b9ad928SBarry Smith EXTERN_C_END 5534b9ad928SBarry Smith 554