1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 44b9ad928SBarry Smith */ 507475bc1SBarry Smith #include <petsc-private/pcimpl.h> 6c6db04a5SJed Brown #include <petscksp.h> /*I "petscksp.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink; 94b9ad928SBarry Smith struct _PC_CompositeLink { 104b9ad928SBarry Smith PC pc; 114b9ad928SBarry Smith PC_CompositeLink next; 12421e10b8SBarry Smith PC_CompositeLink previous; 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 } PC_Composite; 224b9ad928SBarry Smith 234b9ad928SBarry Smith #undef __FUNCT__ 244b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative" 256849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 264b9ad928SBarry Smith { 27dfbe8321SBarry Smith PetscErrorCode ierr; 284b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 294b9ad928SBarry Smith PC_CompositeLink next = jac->head; 304b9ad928SBarry Smith Mat mat = pc->pmat; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 33ce94432eSBarry Smith if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 344b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 354b9ad928SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 364b9ad928SBarry Smith } 3749517cdeSBarry Smith if (pc->useAmat) mat = pc->mat; 382533e041SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 394b9ad928SBarry Smith while (next->next) { 404b9ad928SBarry Smith next = next->next; 414b9ad928SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 42efb30889SBarry Smith ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 4351f519a2SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 44d8fd42c4SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 45efb30889SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 464b9ad928SBarry Smith } 47421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 48421e10b8SBarry Smith while (next->previous) { 49421e10b8SBarry Smith next = next->previous; 50421e10b8SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 51421e10b8SBarry Smith ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 52421e10b8SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 53421e10b8SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 54421e10b8SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 55421e10b8SBarry Smith } 56421e10b8SBarry Smith } 574b9ad928SBarry Smith PetscFunctionReturn(0); 584b9ad928SBarry Smith } 594b9ad928SBarry Smith 602533e041SBarry Smith #undef __FUNCT__ 612533e041SBarry Smith #define __FUNCT__ "PCApplyTranspose_Composite_Multiplicative" 622533e041SBarry Smith static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc,Vec x,Vec y) 632533e041SBarry Smith { 642533e041SBarry Smith PetscErrorCode ierr; 652533e041SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 662533e041SBarry Smith PC_CompositeLink next = jac->head; 672533e041SBarry Smith Mat mat = pc->pmat; 682533e041SBarry Smith 692533e041SBarry Smith PetscFunctionBegin; 70ce94432eSBarry Smith if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 712533e041SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 722533e041SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 732533e041SBarry Smith } 7449517cdeSBarry Smith if (pc->useAmat) mat = pc->mat; 752533e041SBarry Smith /* locate last PC */ 762533e041SBarry Smith while (next->next) { 772533e041SBarry Smith next = next->next; 782533e041SBarry Smith } 792533e041SBarry Smith ierr = PCApplyTranspose(next->pc,x,y);CHKERRQ(ierr); 802533e041SBarry Smith while (next->previous) { 812533e041SBarry Smith next = next->previous; 822533e041SBarry Smith ierr = MatMultTranspose(mat,y,jac->work1);CHKERRQ(ierr); 832533e041SBarry Smith ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 842533e041SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 852533e041SBarry Smith ierr = PCApplyTranspose(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 862533e041SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 872533e041SBarry Smith } 882533e041SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 892533e041SBarry Smith next = jac->head; 902533e041SBarry Smith while (next->next) { 912533e041SBarry Smith next = next->next; 922533e041SBarry Smith ierr = MatMultTranspose(mat,y,jac->work1);CHKERRQ(ierr); 932533e041SBarry Smith ierr = VecWAXPY(jac->work2,-1.0,jac->work1,x);CHKERRQ(ierr); 942533e041SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 952533e041SBarry Smith ierr = PCApplyTranspose(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 962533e041SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 972533e041SBarry Smith } 982533e041SBarry Smith } 992533e041SBarry Smith PetscFunctionReturn(0); 1002533e041SBarry Smith } 1012533e041SBarry Smith 1024b9ad928SBarry Smith /* 1034b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 1044b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 1054b9ad928SBarry Smith alpha I + R 1064b9ad928SBarry Smith */ 1074b9ad928SBarry Smith #undef __FUNCT__ 1084b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special" 1096849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y) 1104b9ad928SBarry Smith { 111dfbe8321SBarry Smith PetscErrorCode ierr; 1124b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1134b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1144b9ad928SBarry Smith 1154b9ad928SBarry Smith PetscFunctionBegin; 116ce94432eSBarry Smith if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 117ce94432eSBarry Smith if (!next->next || next->next->next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs"); 1184b9ad928SBarry Smith 119d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 120d8fd42c4SBarry Smith ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr); 1214b9ad928SBarry Smith PetscFunctionReturn(0); 1224b9ad928SBarry Smith } 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith #undef __FUNCT__ 1254b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive" 1266849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y) 1274b9ad928SBarry Smith { 128dfbe8321SBarry Smith PetscErrorCode ierr; 1294b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1304b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1314b9ad928SBarry Smith 1324b9ad928SBarry Smith PetscFunctionBegin; 133ce94432eSBarry Smith if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 134d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 1354b9ad928SBarry Smith while (next->next) { 1364b9ad928SBarry Smith next = next->next; 13751f519a2SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 138d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 139efb30889SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 1404b9ad928SBarry Smith } 1414b9ad928SBarry Smith PetscFunctionReturn(0); 1424b9ad928SBarry Smith } 1434b9ad928SBarry Smith 1444b9ad928SBarry Smith #undef __FUNCT__ 1452533e041SBarry Smith #define __FUNCT__ "PCApplyTranspose_Composite_Additive" 1462533e041SBarry Smith static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc,Vec x,Vec y) 1472533e041SBarry Smith { 1482533e041SBarry Smith PetscErrorCode ierr; 1492533e041SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1502533e041SBarry Smith PC_CompositeLink next = jac->head; 1512533e041SBarry Smith 1522533e041SBarry Smith PetscFunctionBegin; 153ce94432eSBarry Smith if (!next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC() or -pc_composite_pcs"); 1542533e041SBarry Smith ierr = PCApplyTranspose(next->pc,x,y);CHKERRQ(ierr); 1552533e041SBarry Smith while (next->next) { 1562533e041SBarry Smith next = next->next; 1572533e041SBarry Smith ierr = VecSet(jac->work1,0.0);CHKERRQ(ierr); /* zero since some PC's may not set all entries in the result */ 1582533e041SBarry Smith ierr = PCApplyTranspose(next->pc,x,jac->work1);CHKERRQ(ierr); 1592533e041SBarry Smith ierr = VecAXPY(y,1.0,jac->work1);CHKERRQ(ierr); 1602533e041SBarry Smith } 1612533e041SBarry Smith PetscFunctionReturn(0); 1622533e041SBarry Smith } 1632533e041SBarry Smith 1642533e041SBarry Smith #undef __FUNCT__ 1654b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite" 1666849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc) 1674b9ad928SBarry Smith { 168dfbe8321SBarry Smith PetscErrorCode ierr; 1694b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1704b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1714b9ad928SBarry Smith 1724b9ad928SBarry Smith PetscFunctionBegin; 1734b9ad928SBarry Smith if (!jac->work1) { 17423ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 1754b9ad928SBarry Smith } 1764b9ad928SBarry Smith while (next) { 1774b9ad928SBarry Smith ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1784b9ad928SBarry Smith next = next->next; 1794b9ad928SBarry Smith } 1804b9ad928SBarry Smith PetscFunctionReturn(0); 1814b9ad928SBarry Smith } 1824b9ad928SBarry Smith 1834b9ad928SBarry Smith #undef __FUNCT__ 18469d2c0f9SBarry Smith #define __FUNCT__ "PCReset_Composite" 18569d2c0f9SBarry Smith static PetscErrorCode PCReset_Composite(PC pc) 18669d2c0f9SBarry Smith { 18769d2c0f9SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 18869d2c0f9SBarry Smith PetscErrorCode ierr; 1895f48b12bSBarry Smith PC_CompositeLink next = jac->head; 19069d2c0f9SBarry Smith 19169d2c0f9SBarry Smith PetscFunctionBegin; 19269d2c0f9SBarry Smith while (next) { 19369d2c0f9SBarry Smith ierr = PCReset(next->pc);CHKERRQ(ierr); 19469d2c0f9SBarry Smith next = next->next; 19569d2c0f9SBarry Smith } 1966bf464f9SBarry Smith ierr = VecDestroy(&jac->work1);CHKERRQ(ierr); 1976bf464f9SBarry Smith ierr = VecDestroy(&jac->work2);CHKERRQ(ierr); 19869d2c0f9SBarry Smith PetscFunctionReturn(0); 19969d2c0f9SBarry Smith } 20069d2c0f9SBarry Smith 20169d2c0f9SBarry Smith #undef __FUNCT__ 2024b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 2036849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc) 2044b9ad928SBarry Smith { 2054b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 206dfbe8321SBarry Smith PetscErrorCode ierr; 207724c2c99SHong Zhang PC_CompositeLink next = jac->head,next_tmp; 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith PetscFunctionBegin; 21069d2c0f9SBarry Smith ierr = PCReset_Composite(pc);CHKERRQ(ierr); 2114b9ad928SBarry Smith while (next) { 2126bf464f9SBarry Smith ierr = PCDestroy(&next->pc);CHKERRQ(ierr); 213724c2c99SHong Zhang next_tmp = next; 2144b9ad928SBarry Smith next = next->next; 215724c2c99SHong Zhang ierr = PetscFree(next_tmp);CHKERRQ(ierr); 2164b9ad928SBarry Smith } 217c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2184b9ad928SBarry Smith PetscFunctionReturn(0); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith 2214b9ad928SBarry Smith #undef __FUNCT__ 2224b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 2236849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc) 2244b9ad928SBarry Smith { 2254b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 226dfbe8321SBarry Smith PetscErrorCode ierr; 2279dcbbd2bSBarry Smith PetscInt nmax = 8,i; 2284b9ad928SBarry Smith PC_CompositeLink next; 229e5999256SBarry Smith char *pcs[8]; 230ace3abfcSBarry Smith PetscBool flg; 2314b9ad928SBarry Smith 2324b9ad928SBarry Smith PetscFunctionBegin; 2334b9ad928SBarry Smith ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 2349dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 23551f519a2SBarry Smith if (flg) { 23651f519a2SBarry Smith ierr = PCCompositeSetType(pc,jac->type);CHKERRQ(ierr); 23751f519a2SBarry Smith } 2384b9ad928SBarry Smith ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 2394b9ad928SBarry Smith if (flg) { 2404b9ad928SBarry Smith for (i=0; i<nmax; i++) { 2414b9ad928SBarry Smith ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 242724c2c99SHong Zhang ierr = PetscFree(pcs[i]);CHKERRQ(ierr); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 2434b9ad928SBarry Smith } 2444b9ad928SBarry Smith } 2454b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith next = jac->head; 2484b9ad928SBarry Smith while (next) { 2494b9ad928SBarry Smith ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 2504b9ad928SBarry Smith next = next->next; 2514b9ad928SBarry Smith } 2524b9ad928SBarry Smith PetscFunctionReturn(0); 2534b9ad928SBarry Smith } 2544b9ad928SBarry Smith 2554b9ad928SBarry Smith #undef __FUNCT__ 2564b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite" 2576849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 2584b9ad928SBarry Smith { 2594b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 260dfbe8321SBarry Smith PetscErrorCode ierr; 2614b9ad928SBarry Smith PC_CompositeLink next = jac->head; 262ace3abfcSBarry Smith PetscBool iascii; 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith PetscFunctionBegin; 265251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 26632077d6dSBarry Smith if (iascii) { 2679dcbbd2bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 2684b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 2694b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2704b9ad928SBarry Smith } 27132077d6dSBarry Smith if (iascii) { 2724b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 2734b9ad928SBarry Smith } 2744b9ad928SBarry Smith while (next) { 2754b9ad928SBarry Smith ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 2764b9ad928SBarry Smith next = next->next; 2774b9ad928SBarry Smith } 27832077d6dSBarry Smith if (iascii) { 2794b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2804b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2814b9ad928SBarry Smith } 2824b9ad928SBarry Smith PetscFunctionReturn(0); 2834b9ad928SBarry Smith } 2844b9ad928SBarry Smith 2854b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith #undef __FUNCT__ 2884b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 2891e6b0712SBarry Smith static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2904b9ad928SBarry Smith { 2914b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2925fd66863SKarl Rupp 2934b9ad928SBarry Smith PetscFunctionBegin; 2944b9ad928SBarry Smith jac->alpha = alpha; 2954b9ad928SBarry Smith PetscFunctionReturn(0); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith 2984b9ad928SBarry Smith #undef __FUNCT__ 2994b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite" 3001e6b0712SBarry Smith static PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 3014b9ad928SBarry Smith { 302fad69fbaSJed Brown PC_Composite *jac = (PC_Composite*)pc->data; 303fad69fbaSJed Brown 3044b9ad928SBarry Smith PetscFunctionBegin; 3054b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 3064b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 3072533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 308421e10b8SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 3094b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 3102533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 3114b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 3124b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 3130298fd71SBarry Smith pc->ops->applytranspose = NULL; 314ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 315fad69fbaSJed Brown jac->type = type; 3164b9ad928SBarry Smith PetscFunctionReturn(0); 3174b9ad928SBarry Smith } 3184b9ad928SBarry Smith 3194b9ad928SBarry Smith #undef __FUNCT__ 3204b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite" 3211e6b0712SBarry Smith static PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type) 3224b9ad928SBarry Smith { 3234b9ad928SBarry Smith PC_Composite *jac; 3245a9f2f41SSatish Balay PC_CompositeLink next,ilink; 325dfbe8321SBarry Smith PetscErrorCode ierr; 32679416396SBarry Smith PetscInt cnt = 0; 3272dcb1b2aSMatthew Knepley const char *prefix; 3282dcb1b2aSMatthew Knepley char newprefix[8]; 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith PetscFunctionBegin; 331*b00a9115SJed Brown ierr = PetscNewLog(pc,&ilink);CHKERRQ(ierr); 3325a9f2f41SSatish Balay ilink->next = 0; 333ce94432eSBarry Smith ierr = PCCreate(PetscObjectComm((PetscObject)pc),&ilink->pc);CHKERRQ(ierr); 334ace3abfcSBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ilink->pc);CHKERRQ(ierr); 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3374b9ad928SBarry Smith next = jac->head; 3384b9ad928SBarry Smith if (!next) { 3395a9f2f41SSatish Balay jac->head = ilink; 3400298fd71SBarry Smith ilink->previous = NULL; 3414b9ad928SBarry Smith } else { 3424b9ad928SBarry Smith cnt++; 3434b9ad928SBarry Smith while (next->next) { 3444b9ad928SBarry Smith next = next->next; 3454b9ad928SBarry Smith cnt++; 3464b9ad928SBarry Smith } 3475a9f2f41SSatish Balay next->next = ilink; 348421e10b8SBarry Smith ilink->previous = next; 3494b9ad928SBarry Smith } 3504b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3515a9f2f41SSatish Balay ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 35213f74950SBarry Smith sprintf(newprefix,"sub_%d_",(int)cnt); 3535a9f2f41SSatish Balay ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 3544b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 3555a9f2f41SSatish Balay ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 3564b9ad928SBarry Smith PetscFunctionReturn(0); 3574b9ad928SBarry Smith } 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith #undef __FUNCT__ 3604b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite" 3611e6b0712SBarry Smith static PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 3624b9ad928SBarry Smith { 3634b9ad928SBarry Smith PC_Composite *jac; 3644b9ad928SBarry Smith PC_CompositeLink next; 36579416396SBarry Smith PetscInt i; 3664b9ad928SBarry Smith 3674b9ad928SBarry Smith PetscFunctionBegin; 3684b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3694b9ad928SBarry Smith next = jac->head; 3704b9ad928SBarry Smith for (i=0; i<n; i++) { 371ce94432eSBarry Smith if (!next->next) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 3724b9ad928SBarry Smith next = next->next; 3734b9ad928SBarry Smith } 3744b9ad928SBarry Smith *subpc = next->pc; 3754b9ad928SBarry Smith PetscFunctionReturn(0); 3764b9ad928SBarry Smith } 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 3794b9ad928SBarry Smith #undef __FUNCT__ 3804b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType" 381f39d8e23SSatish Balay /*@ 3824b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 3834b9ad928SBarry Smith 384ad4df100SBarry Smith Logically Collective on PC 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith Input Parameter: 3872a6744ebSBarry Smith + pc - the preconditioner context 3882a6744ebSBarry Smith - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith Options Database Key: 3914b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith Level: Developer 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3964b9ad928SBarry Smith @*/ 3977087cfbeSBarry Smith PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 3984b9ad928SBarry Smith { 3997bb14e67SBarry Smith PetscErrorCode ierr; 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith PetscFunctionBegin; 4020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 403c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 4047bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type));CHKERRQ(ierr); 4054b9ad928SBarry Smith PetscFunctionReturn(0); 4064b9ad928SBarry Smith } 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith #undef __FUNCT__ 4094b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha" 410f39d8e23SSatish Balay /*@ 4114b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 4124b9ad928SBarry Smith for alphaI + R + S 4134b9ad928SBarry Smith 414ad4df100SBarry Smith Logically Collective on PC 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith Input Parameter: 4174b9ad928SBarry Smith + pc - the preconditioner context 4184b9ad928SBarry Smith - alpha - scale on identity 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith Level: Developer 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 4234b9ad928SBarry Smith @*/ 4247087cfbeSBarry Smith PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 4254b9ad928SBarry Smith { 4264ac538c5SBarry Smith PetscErrorCode ierr; 4274b9ad928SBarry Smith 4284b9ad928SBarry Smith PetscFunctionBegin; 4290700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 430c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(pc,alpha,2); 4314ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha));CHKERRQ(ierr); 4324b9ad928SBarry Smith PetscFunctionReturn(0); 4334b9ad928SBarry Smith } 4344b9ad928SBarry Smith 4354b9ad928SBarry Smith #undef __FUNCT__ 4364b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC" 4374b9ad928SBarry Smith /*@C 4384b9ad928SBarry Smith PCCompositeAddPC - Adds another PC to the composite PC. 4394b9ad928SBarry Smith 4404b9ad928SBarry Smith Collective on PC 4414b9ad928SBarry Smith 4424b9ad928SBarry Smith Input Parameters: 4432a6744ebSBarry Smith + pc - the preconditioner context 4442a6744ebSBarry Smith - type - the type of the new preconditioner 4454b9ad928SBarry Smith 4464b9ad928SBarry Smith Level: Developer 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith .keywords: PC, composite preconditioner, add 4494b9ad928SBarry Smith @*/ 45019fd82e9SBarry Smith PetscErrorCode PCCompositeAddPC(PC pc,PCType type) 4514b9ad928SBarry Smith { 4524ac538c5SBarry Smith PetscErrorCode ierr; 4534b9ad928SBarry Smith 4544b9ad928SBarry Smith PetscFunctionBegin; 4550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 45619fd82e9SBarry Smith ierr = PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PCType),(pc,type));CHKERRQ(ierr); 4574b9ad928SBarry Smith PetscFunctionReturn(0); 4584b9ad928SBarry Smith } 4594b9ad928SBarry Smith 4604b9ad928SBarry Smith #undef __FUNCT__ 4614b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC" 462f39d8e23SSatish Balay /*@ 4634b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 4644b9ad928SBarry Smith 4654b9ad928SBarry Smith Not Collective 4664b9ad928SBarry Smith 4674b9ad928SBarry Smith Input Parameter: 4682a6744ebSBarry Smith + pc - the preconditioner context 4692a6744ebSBarry Smith - n - the number of the pc requested 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith Output Parameters: 4724b9ad928SBarry Smith . subpc - the PC requested 4734b9ad928SBarry Smith 4744b9ad928SBarry Smith Level: Developer 4754b9ad928SBarry Smith 4764b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner 4774b9ad928SBarry Smith 4784b9ad928SBarry Smith .seealso: PCCompositeAddPC() 4794b9ad928SBarry Smith @*/ 4807087cfbeSBarry Smith PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 4814b9ad928SBarry Smith { 4824ac538c5SBarry Smith PetscErrorCode ierr; 4834b9ad928SBarry Smith 4844b9ad928SBarry Smith PetscFunctionBegin; 4850700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4864482741eSBarry Smith PetscValidPointer(subpc,3); 4874ac538c5SBarry Smith ierr = PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC*),(pc,n,subpc));CHKERRQ(ierr); 4884b9ad928SBarry Smith PetscFunctionReturn(0); 4894b9ad928SBarry Smith } 4904b9ad928SBarry Smith 4914b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 4924b9ad928SBarry Smith 4934b9ad928SBarry Smith /*MC 4944b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith Options Database Keys: 4972eab2d5bSJungho Lee + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 49849517cdeSBarry Smith . -pc_use_amat - Activates PCSetUseAmat() 49951f519a2SBarry Smith - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 5004b9ad928SBarry Smith 5014b9ad928SBarry Smith Level: intermediate 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith Concepts: composing solvers 5044b9ad928SBarry Smith 5054b9ad928SBarry Smith Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 5064b9ad928SBarry Smith inner PCs to be PCKSP. 5074b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 508b3ef52cdSBarry Smith the incorrect answer) unless you use KSPFGMRES as the outer Krylov method 5094b9ad928SBarry Smith 5104b9ad928SBarry Smith 5114b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 5124b9ad928SBarry Smith PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 51349517cdeSBarry Smith PCCompositeGetPC(), PCSetUseAmat() 5144b9ad928SBarry Smith 5154b9ad928SBarry Smith M*/ 5164b9ad928SBarry Smith 5174b9ad928SBarry Smith #undef __FUNCT__ 5184b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite" 5198cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 5204b9ad928SBarry Smith { 521dfbe8321SBarry Smith PetscErrorCode ierr; 5224b9ad928SBarry Smith PC_Composite *jac; 5234b9ad928SBarry Smith 5244b9ad928SBarry Smith PetscFunctionBegin; 525*b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 5262fa5cd67SKarl Rupp 5274b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 5282533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 5294b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 53069d2c0f9SBarry Smith pc->ops->reset = PCReset_Composite; 5314b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 5324b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 5334b9ad928SBarry Smith pc->ops->view = PCView_Composite; 5344b9ad928SBarry Smith pc->ops->applyrichardson = 0; 5354b9ad928SBarry Smith 5364b9ad928SBarry Smith pc->data = (void*)jac; 5374b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5384b9ad928SBarry Smith jac->work1 = 0; 5394b9ad928SBarry Smith jac->work2 = 0; 5404b9ad928SBarry Smith jac->head = 0; 5414b9ad928SBarry Smith 542bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",PCCompositeSetType_Composite);CHKERRQ(ierr); 543bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",PCCompositeAddPC_Composite);CHKERRQ(ierr); 544bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",PCCompositeGetPC_Composite);CHKERRQ(ierr); 545bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 5464b9ad928SBarry Smith PetscFunctionReturn(0); 5474b9ad928SBarry Smith } 5484b9ad928SBarry Smith 549