1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 44b9ad928SBarry Smith */ 5af0996ceSBarry 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 236849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 244b9ad928SBarry Smith { 254b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 264b9ad928SBarry Smith PC_CompositeLink next = jac->head; 274b9ad928SBarry Smith Mat mat = pc->pmat; 284b9ad928SBarry Smith 294b9ad928SBarry Smith PetscFunctionBegin; 30450d59ebSPatrick Farrell 31*28b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 32450d59ebSPatrick Farrell 33450d59ebSPatrick Farrell /* Set the reuse flag on children PCs */ 34450d59ebSPatrick Farrell while (next) { 355f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetReusePreconditioner(next->pc,pc->reusepreconditioner)); 36450d59ebSPatrick Farrell next = next->next; 37450d59ebSPatrick Farrell } 38450d59ebSPatrick Farrell next = jac->head; 39450d59ebSPatrick Farrell 404b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(jac->work1,&jac->work2)); 424b9ad928SBarry Smith } 4349517cdeSBarry Smith if (pc->useAmat) mat = pc->mat; 445f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,x,y)); /* y <- B x */ 454b9ad928SBarry Smith while (next->next) { 464b9ad928SBarry Smith next = next->next; 475f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,y,jac->work1)); /* work1 <- A y */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(jac->work2,-1.0,jac->work1,x)); /* work2 <- x - work1 */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,jac->work2,jac->work1)); /* work1 <- C work2 */ 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); /* y <- y + work1 = B x + C (x - A B x) = (B + C (1 - A B)) x */ 514b9ad928SBarry Smith } 52421e10b8SBarry Smith if (jac->type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 53421e10b8SBarry Smith while (next->previous) { 54421e10b8SBarry Smith next = next->previous; 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(mat,y,jac->work1)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,jac->work2,jac->work1)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); 59421e10b8SBarry Smith } 60421e10b8SBarry Smith } 614b9ad928SBarry Smith PetscFunctionReturn(0); 624b9ad928SBarry Smith } 634b9ad928SBarry Smith 642533e041SBarry Smith static PetscErrorCode PCApplyTranspose_Composite_Multiplicative(PC pc,Vec x,Vec y) 652533e041SBarry Smith { 662533e041SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 672533e041SBarry Smith PC_CompositeLink next = jac->head; 682533e041SBarry Smith Mat mat = pc->pmat; 692533e041SBarry Smith 702533e041SBarry Smith PetscFunctionBegin; 71*28b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 722533e041SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 735f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(jac->work1,&jac->work2)); 742533e041SBarry Smith } 7549517cdeSBarry Smith if (pc->useAmat) mat = pc->mat; 762533e041SBarry Smith /* locate last PC */ 772533e041SBarry Smith while (next->next) { 782533e041SBarry Smith next = next->next; 792533e041SBarry Smith } 805f80ce2aSJacob Faibussowitsch CHKERRQ(PCApplyTranspose(next->pc,x,y)); 812533e041SBarry Smith while (next->previous) { 822533e041SBarry Smith next = next->previous; 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(mat,y,jac->work1)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PCApplyTranspose(next->pc,jac->work2,jac->work1)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); 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; 925f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTranspose(mat,y,jac->work1)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(PCApplyTranspose(next->pc,jac->work2,jac->work1)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); 962533e041SBarry Smith } 972533e041SBarry Smith } 982533e041SBarry Smith PetscFunctionReturn(0); 992533e041SBarry Smith } 1002533e041SBarry Smith 1014b9ad928SBarry Smith /* 1024b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 1034b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 1044b9ad928SBarry Smith alpha I + R 1054b9ad928SBarry Smith */ 1066849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y) 1074b9ad928SBarry Smith { 1084b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1094b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1104b9ad928SBarry Smith 1114b9ad928SBarry Smith PetscFunctionBegin; 112*28b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 1132c71b3e2SJacob Faibussowitsch PetscCheckFalse(!next->next || next->next->next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs"); 1144b9ad928SBarry Smith 115450d59ebSPatrick Farrell /* Set the reuse flag on children PCs */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetReusePreconditioner(next->pc,pc->reusepreconditioner)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetReusePreconditioner(next->next->pc,pc->reusepreconditioner)); 118450d59ebSPatrick Farrell 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,x,jac->work1)); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->next->pc,jac->work1,y)); 1214b9ad928SBarry Smith PetscFunctionReturn(0); 1224b9ad928SBarry Smith } 1234b9ad928SBarry Smith 1246849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y) 1254b9ad928SBarry Smith { 1264b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1274b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith PetscFunctionBegin; 130*28b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 131450d59ebSPatrick Farrell 132450d59ebSPatrick Farrell /* Set the reuse flag on children PCs */ 133450d59ebSPatrick Farrell while (next) { 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetReusePreconditioner(next->pc,pc->reusepreconditioner)); 135450d59ebSPatrick Farrell next = next->next; 136450d59ebSPatrick Farrell } 137450d59ebSPatrick Farrell next = jac->head; 138450d59ebSPatrick Farrell 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,x,y)); 1404b9ad928SBarry Smith while (next->next) { 1414b9ad928SBarry Smith next = next->next; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PCApply(next->pc,x,jac->work1)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); 1444b9ad928SBarry Smith } 1454b9ad928SBarry Smith PetscFunctionReturn(0); 1464b9ad928SBarry Smith } 1474b9ad928SBarry Smith 1482533e041SBarry Smith static PetscErrorCode PCApplyTranspose_Composite_Additive(PC pc,Vec x,Vec y) 1492533e041SBarry Smith { 1502533e041SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1512533e041SBarry Smith PC_CompositeLink next = jac->head; 1522533e041SBarry Smith 1532533e041SBarry Smith PetscFunctionBegin; 154*28b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PCApplyTranspose(next->pc,x,y)); 1562533e041SBarry Smith while (next->next) { 1572533e041SBarry Smith next = next->next; 1585f80ce2aSJacob Faibussowitsch CHKERRQ(PCApplyTranspose(next->pc,x,jac->work1)); 1595f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(y,1.0,jac->work1)); 1602533e041SBarry Smith } 1612533e041SBarry Smith PetscFunctionReturn(0); 1622533e041SBarry Smith } 1632533e041SBarry Smith 1646849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc) 1654b9ad928SBarry Smith { 1664b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1674b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1685a78d018SMatthew G. Knepley DM dm; 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith PetscFunctionBegin; 1714b9ad928SBarry Smith if (!jac->work1) { 1725f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(pc->pmat,&jac->work1,NULL)); 1734b9ad928SBarry Smith } 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetDM(pc,&dm)); 1754b9ad928SBarry Smith while (next) { 1762b1d202aSBarry Smith if (!next->pc->dm) { 1775f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetDM(next->pc,dm)); 1782b1d202aSBarry Smith } 1792b1d202aSBarry Smith if (!next->pc->mat) { 1805f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOperators(next->pc,pc->mat,pc->pmat)); 1812b1d202aSBarry Smith } 1824b9ad928SBarry Smith next = next->next; 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith PetscFunctionReturn(0); 1854b9ad928SBarry Smith } 1864b9ad928SBarry Smith 18769d2c0f9SBarry Smith static PetscErrorCode PCReset_Composite(PC pc) 18869d2c0f9SBarry Smith { 18969d2c0f9SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1905f48b12bSBarry Smith PC_CompositeLink next = jac->head; 19169d2c0f9SBarry Smith 19269d2c0f9SBarry Smith PetscFunctionBegin; 19369d2c0f9SBarry Smith while (next) { 1945f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset(next->pc)); 19569d2c0f9SBarry Smith next = next->next; 19669d2c0f9SBarry Smith } 1975f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->work1)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&jac->work2)); 19969d2c0f9SBarry Smith PetscFunctionReturn(0); 20069d2c0f9SBarry Smith } 20169d2c0f9SBarry Smith 2026849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc) 2034b9ad928SBarry Smith { 2044b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 205724c2c99SHong Zhang PC_CompositeLink next = jac->head,next_tmp; 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith PetscFunctionBegin; 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_Composite(pc)); 2094b9ad928SBarry Smith while (next) { 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PCDestroy(&next->pc)); 211724c2c99SHong Zhang next_tmp = next; 2124b9ad928SBarry Smith next = next->next; 2135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(next_tmp)); 2144b9ad928SBarry Smith } 2155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc->data)); 2164b9ad928SBarry Smith PetscFunctionReturn(0); 2174b9ad928SBarry Smith } 2184b9ad928SBarry Smith 2194416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,PC pc) 2204b9ad928SBarry Smith { 2214b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2229dcbbd2bSBarry Smith PetscInt nmax = 8,i; 2234b9ad928SBarry Smith PC_CompositeLink next; 224e5999256SBarry Smith char *pcs[8]; 225ace3abfcSBarry Smith PetscBool flg; 2264b9ad928SBarry Smith 2274b9ad928SBarry Smith PetscFunctionBegin; 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Composite preconditioner options")); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg)); 23051f519a2SBarry Smith if (flg) { 2315f80ce2aSJacob Faibussowitsch CHKERRQ(PCCompositeSetType(pc,jac->type)); 23251f519a2SBarry Smith } 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPCType",pcs,&nmax,&flg)); 2344b9ad928SBarry Smith if (flg) { 2354b9ad928SBarry Smith for (i=0; i<nmax; i++) { 2365f80ce2aSJacob Faibussowitsch CHKERRQ(PCCompositeAddPCType(pc,pcs[i])); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 2384b9ad928SBarry Smith } 2394b9ad928SBarry Smith } 2405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith next = jac->head; 2434b9ad928SBarry Smith while (next) { 2445f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetFromOptions(next->pc)); 2454b9ad928SBarry Smith next = next->next; 2464b9ad928SBarry Smith } 2474b9ad928SBarry Smith PetscFunctionReturn(0); 2484b9ad928SBarry Smith } 2494b9ad928SBarry Smith 2506849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 2514b9ad928SBarry Smith { 2524b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2534b9ad928SBarry Smith PC_CompositeLink next = jac->head; 254ace3abfcSBarry Smith PetscBool iascii; 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith PetscFunctionBegin; 2575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 25832077d6dSBarry Smith if (iascii) { 2595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type])); 2605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n")); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"---------------------------------\n")); 2624b9ad928SBarry Smith } 26332077d6dSBarry Smith if (iascii) { 2645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 2654b9ad928SBarry Smith } 2664b9ad928SBarry Smith while (next) { 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PCView(next->pc,viewer)); 2684b9ad928SBarry Smith next = next->next; 2694b9ad928SBarry Smith } 27032077d6dSBarry Smith if (iascii) { 2715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 2725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer,"---------------------------------\n")); 2734b9ad928SBarry Smith } 2744b9ad928SBarry Smith PetscFunctionReturn(0); 2754b9ad928SBarry Smith } 2764b9ad928SBarry Smith 2774b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2784b9ad928SBarry Smith 2791e6b0712SBarry Smith static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2804b9ad928SBarry Smith { 2814b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2825fd66863SKarl Rupp 2834b9ad928SBarry Smith PetscFunctionBegin; 2844b9ad928SBarry Smith jac->alpha = alpha; 2854b9ad928SBarry Smith PetscFunctionReturn(0); 2864b9ad928SBarry Smith } 2874b9ad928SBarry Smith 2881e6b0712SBarry Smith static PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 2894b9ad928SBarry Smith { 290fad69fbaSJed Brown PC_Composite *jac = (PC_Composite*)pc->data; 291fad69fbaSJed Brown 2924b9ad928SBarry Smith PetscFunctionBegin; 2934b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 2944b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 2952533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 296421e10b8SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 2974b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 2982533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 2994b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 3004b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 3010298fd71SBarry Smith pc->ops->applytranspose = NULL; 302a7261c6bSprj- } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown composite preconditioner type"); 303fad69fbaSJed Brown jac->type = type; 3044b9ad928SBarry Smith PetscFunctionReturn(0); 3054b9ad928SBarry Smith } 3064b9ad928SBarry Smith 307c60c7ad4SBarry Smith static PetscErrorCode PCCompositeGetType_Composite(PC pc,PCCompositeType *type) 308c60c7ad4SBarry Smith { 309c60c7ad4SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 310c60c7ad4SBarry Smith 311c60c7ad4SBarry Smith PetscFunctionBegin; 312c60c7ad4SBarry Smith *type = jac->type; 313c60c7ad4SBarry Smith PetscFunctionReturn(0); 314c60c7ad4SBarry Smith } 315c60c7ad4SBarry Smith 3168aa07aa6SMatthew G. Knepley static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc) 3174b9ad928SBarry Smith { 3184b9ad928SBarry Smith PC_Composite *jac; 3195a9f2f41SSatish Balay PC_CompositeLink next, ilink; 32079416396SBarry Smith PetscInt cnt = 0; 3212dcb1b2aSMatthew Knepley const char *prefix; 322d726e3a5SJed Brown char newprefix[20]; 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith PetscFunctionBegin; 3255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc, &ilink)); 3260a545947SLisandro Dalcin ilink->next = NULL; 3278aa07aa6SMatthew G. Knepley ilink->pc = subpc; 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith jac = (PC_Composite *) pc->data; 3304b9ad928SBarry Smith next = jac->head; 3314b9ad928SBarry Smith if (!next) { 3325a9f2f41SSatish Balay jac->head = ilink; 3330298fd71SBarry Smith ilink->previous = NULL; 3344b9ad928SBarry Smith } else { 3354b9ad928SBarry Smith cnt++; 3364b9ad928SBarry Smith while (next->next) { 3374b9ad928SBarry Smith next = next->next; 3384b9ad928SBarry Smith cnt++; 3394b9ad928SBarry Smith } 3405a9f2f41SSatish Balay next->next = ilink; 341421e10b8SBarry Smith ilink->previous = next; 3424b9ad928SBarry Smith } 3435f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc, &prefix)); 3445f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetOptionsPrefix(subpc, prefix)); 3455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(newprefix, 20, "sub_%d_", (int) cnt)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(PCAppendOptionsPrefix(subpc, newprefix)); 3475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject) subpc)); 3488aa07aa6SMatthew G. Knepley PetscFunctionReturn(0); 3498aa07aa6SMatthew G. Knepley } 3508aa07aa6SMatthew G. Knepley 3518aa07aa6SMatthew G. Knepley static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type) 3528aa07aa6SMatthew G. Knepley { 3538aa07aa6SMatthew G. Knepley PC subpc; 3548aa07aa6SMatthew G. Knepley 3558aa07aa6SMatthew G. Knepley PetscFunctionBegin; 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PCCreate(PetscObjectComm((PetscObject)pc), &subpc)); 3575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1)); 3585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject) pc, (PetscObject) subpc)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(PCCompositeAddPC_Composite(pc, subpc)); 3604b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 3615f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(subpc, type)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(PCDestroy(&subpc)); 3634b9ad928SBarry Smith PetscFunctionReturn(0); 3644b9ad928SBarry Smith } 3654b9ad928SBarry Smith 3668e6eba06SBarry Smith static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc,PetscInt *n) 3678e6eba06SBarry Smith { 3688e6eba06SBarry Smith PC_Composite *jac; 3698e6eba06SBarry Smith PC_CompositeLink next; 3708e6eba06SBarry Smith 3718e6eba06SBarry Smith PetscFunctionBegin; 3728e6eba06SBarry Smith jac = (PC_Composite*)pc->data; 3738e6eba06SBarry Smith next = jac->head; 3748e6eba06SBarry Smith *n = 0; 3758e6eba06SBarry Smith while (next) { 3768e6eba06SBarry Smith next = next->next; 3778e6eba06SBarry Smith (*n) ++; 3788e6eba06SBarry Smith } 3798e6eba06SBarry Smith PetscFunctionReturn(0); 3808e6eba06SBarry Smith } 3818e6eba06SBarry Smith 3821e6b0712SBarry Smith static PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 3834b9ad928SBarry Smith { 3844b9ad928SBarry Smith PC_Composite *jac; 3854b9ad928SBarry Smith PC_CompositeLink next; 38679416396SBarry Smith PetscInt i; 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith PetscFunctionBegin; 3894b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3904b9ad928SBarry Smith next = jac->head; 3914b9ad928SBarry Smith for (i=0; i<n; i++) { 392*28b400f6SJacob Faibussowitsch PetscCheck(next->next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 3934b9ad928SBarry Smith next = next->next; 3944b9ad928SBarry Smith } 3954b9ad928SBarry Smith *subpc = next->pc; 3964b9ad928SBarry Smith PetscFunctionReturn(0); 3974b9ad928SBarry Smith } 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 400f39d8e23SSatish Balay /*@ 4014b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 4024b9ad928SBarry Smith 403ad4df100SBarry Smith Logically Collective on PC 4044b9ad928SBarry Smith 405c60c7ad4SBarry Smith Input Parameters: 4062a6744ebSBarry Smith + pc - the preconditioner context 4072a6744ebSBarry Smith - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith Options Database Key: 4104b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Level: Developer 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith @*/ 4157087cfbeSBarry Smith PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 4164b9ad928SBarry Smith { 4174b9ad928SBarry Smith PetscFunctionBegin; 4180700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 419c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type))); 4214b9ad928SBarry Smith PetscFunctionReturn(0); 4224b9ad928SBarry Smith } 4234b9ad928SBarry Smith 424c60c7ad4SBarry Smith /*@ 425721f67b5SBarry Smith PCCompositeGetType - Gets the type of composite preconditioner. 426c60c7ad4SBarry Smith 427c60c7ad4SBarry Smith Logically Collective on PC 428c60c7ad4SBarry Smith 429c60c7ad4SBarry Smith Input Parameter: 430c60c7ad4SBarry Smith . pc - the preconditioner context 431c60c7ad4SBarry Smith 432c60c7ad4SBarry Smith Output Parameter: 433c60c7ad4SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 434c60c7ad4SBarry Smith 435c60c7ad4SBarry Smith Options Database Key: 436c60c7ad4SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 437c60c7ad4SBarry Smith 438c60c7ad4SBarry Smith Level: Developer 439c60c7ad4SBarry Smith 440c60c7ad4SBarry Smith @*/ 441c60c7ad4SBarry Smith PetscErrorCode PCCompositeGetType(PC pc,PCCompositeType *type) 442c60c7ad4SBarry Smith { 443c60c7ad4SBarry Smith PetscFunctionBegin; 444c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(pc,"PCCompositeGetType_C",(PC,PCCompositeType*),(pc,type))); 446c60c7ad4SBarry Smith PetscFunctionReturn(0); 447c60c7ad4SBarry Smith } 448c60c7ad4SBarry Smith 449f39d8e23SSatish Balay /*@ 4504b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 4514b9ad928SBarry Smith for alphaI + R + S 4524b9ad928SBarry Smith 453ad4df100SBarry Smith Logically Collective on PC 4544b9ad928SBarry Smith 455d8d19677SJose E. Roman Input Parameters: 4564b9ad928SBarry Smith + pc - the preconditioner context 4574b9ad928SBarry Smith - alpha - scale on identity 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith Level: Developer 4604b9ad928SBarry Smith 4614b9ad928SBarry Smith @*/ 4627087cfbeSBarry Smith PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 4634b9ad928SBarry Smith { 4644b9ad928SBarry Smith PetscFunctionBegin; 4650700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 466c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(pc,alpha,2); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha))); 4684b9ad928SBarry Smith PetscFunctionReturn(0); 4694b9ad928SBarry Smith } 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith /*@C 4728aa07aa6SMatthew G. Knepley PCCompositeAddPCType - Adds another PC of the given type to the composite PC. 4734b9ad928SBarry Smith 4744b9ad928SBarry Smith Collective on PC 4754b9ad928SBarry Smith 4764b9ad928SBarry Smith Input Parameters: 4772a6744ebSBarry Smith + pc - the preconditioner context 4782a6744ebSBarry Smith - type - the type of the new preconditioner 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith Level: Developer 4814b9ad928SBarry Smith 4828aa07aa6SMatthew G. Knepley .seealso: PCCompositeAddPC() 4834b9ad928SBarry Smith @*/ 4848aa07aa6SMatthew G. Knepley PetscErrorCode PCCompositeAddPCType(PC pc,PCType type) 4854b9ad928SBarry Smith { 4864b9ad928SBarry Smith PetscFunctionBegin; 4870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCCompositeAddPCType_C",(PC,PCType),(pc,type))); 4898aa07aa6SMatthew G. Knepley PetscFunctionReturn(0); 4908aa07aa6SMatthew G. Knepley } 4918aa07aa6SMatthew G. Knepley 4928aa07aa6SMatthew G. Knepley /*@ 4938aa07aa6SMatthew G. Knepley PCCompositeAddPC - Adds another PC to the composite PC. 4948aa07aa6SMatthew G. Knepley 4958aa07aa6SMatthew G. Knepley Collective on PC 4968aa07aa6SMatthew G. Knepley 4978aa07aa6SMatthew G. Knepley Input Parameters: 4988aa07aa6SMatthew G. Knepley + pc - the preconditioner context 4998aa07aa6SMatthew G. Knepley - subpc - the new preconditioner 5008aa07aa6SMatthew G. Knepley 5018aa07aa6SMatthew G. Knepley Level: Developer 5028aa07aa6SMatthew G. Knepley 5038aa07aa6SMatthew G. Knepley .seealso: PCCompositeAddPCType() 5048aa07aa6SMatthew G. Knepley @*/ 5058aa07aa6SMatthew G. Knepley PetscErrorCode PCCompositeAddPC(PC pc, PC subpc) 5068aa07aa6SMatthew G. Knepley { 5078aa07aa6SMatthew G. Knepley PetscFunctionBegin; 5088aa07aa6SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5098aa07aa6SMatthew G. Knepley PetscValidHeaderSpecific(subpc,PC_CLASSID,2); 5105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PC),(pc,subpc))); 5114b9ad928SBarry Smith PetscFunctionReturn(0); 5124b9ad928SBarry Smith } 5134b9ad928SBarry Smith 5148e6eba06SBarry Smith /*@ 5158e6eba06SBarry Smith PCCompositeGetNumberPC - Gets the number of PC objects in the composite PC. 5168e6eba06SBarry Smith 5178e6eba06SBarry Smith Not Collective 5188e6eba06SBarry Smith 5198e6eba06SBarry Smith Input Parameter: 5208e6eba06SBarry Smith . pc - the preconditioner context 5218e6eba06SBarry Smith 5228e6eba06SBarry Smith Output Parameter: 5238e6eba06SBarry Smith . num - the number of sub pcs 5248e6eba06SBarry Smith 5258e6eba06SBarry Smith Level: Developer 5268e6eba06SBarry Smith 5278e6eba06SBarry Smith .seealso: PCCompositeGetPC() 5288e6eba06SBarry Smith @*/ 5298e6eba06SBarry Smith PetscErrorCode PCCompositeGetNumberPC(PC pc,PetscInt *num) 5308e6eba06SBarry Smith { 5318e6eba06SBarry Smith PetscFunctionBegin; 5328e6eba06SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5338e6eba06SBarry Smith PetscValidIntPointer(num,2); 5345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(pc,"PCCompositeGetNumberPC_C",(PC,PetscInt*),(pc,num))); 5358e6eba06SBarry Smith PetscFunctionReturn(0); 5368e6eba06SBarry Smith } 5378e6eba06SBarry Smith 538f39d8e23SSatish Balay /*@ 5394b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 5404b9ad928SBarry Smith 5414b9ad928SBarry Smith Not Collective 5424b9ad928SBarry Smith 543d8d19677SJose E. Roman Input Parameters: 5442a6744ebSBarry Smith + pc - the preconditioner context 5452a6744ebSBarry Smith - n - the number of the pc requested 5464b9ad928SBarry Smith 5474b9ad928SBarry Smith Output Parameters: 5484b9ad928SBarry Smith . subpc - the PC requested 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith Level: Developer 5514b9ad928SBarry Smith 55295452b02SPatrick Sanan Notes: 55395452b02SPatrick Sanan To use a different operator to construct one of the inner preconditioners first call PCCompositeGetPC(), then 5542b1d202aSBarry Smith call PCSetOperators() on that PC. 5552b1d202aSBarry Smith 5568aa07aa6SMatthew G. Knepley .seealso: PCCompositeAddPCType(), PCCompositeGetNumberPC(), PCSetOperators() 5574b9ad928SBarry Smith @*/ 5587087cfbeSBarry Smith PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 5594b9ad928SBarry Smith { 5604b9ad928SBarry Smith PetscFunctionBegin; 5610700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5624482741eSBarry Smith PetscValidPointer(subpc,3); 5635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC*),(pc,n,subpc))); 5644b9ad928SBarry Smith PetscFunctionReturn(0); 5654b9ad928SBarry Smith } 5664b9ad928SBarry Smith 5674b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 5684b9ad928SBarry Smith 5694b9ad928SBarry Smith /*MC 5704b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 5714b9ad928SBarry Smith 5724b9ad928SBarry Smith Options Database Keys: 5732eab2d5bSJungho Lee + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 5742b1d202aSBarry Smith . -pc_use_amat - activates PCSetUseAmat() 57551f519a2SBarry Smith - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 5764b9ad928SBarry Smith 5774b9ad928SBarry Smith Level: intermediate 5784b9ad928SBarry Smith 57995452b02SPatrick Sanan Notes: 58095452b02SPatrick Sanan To use a Krylov method inside the composite preconditioner, set the PCType of one or more 5814b9ad928SBarry Smith inner PCs to be PCKSP. 5824b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 583b3ef52cdSBarry Smith the incorrect answer) unless you use KSPFGMRES as the outer Krylov method 5842b1d202aSBarry Smith To use a different operator to construct one of the inner preconditioners first call PCCompositeGetPC(), then 5852b1d202aSBarry Smith call PCSetOperators() on that PC. 5864b9ad928SBarry Smith 5874b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 5888aa07aa6SMatthew G. Knepley PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPCType(), 58949517cdeSBarry Smith PCCompositeGetPC(), PCSetUseAmat() 5904b9ad928SBarry Smith 5914b9ad928SBarry Smith M*/ 5924b9ad928SBarry Smith 5938cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 5944b9ad928SBarry Smith { 5954b9ad928SBarry Smith PC_Composite *jac; 5964b9ad928SBarry Smith 5974b9ad928SBarry Smith PetscFunctionBegin; 5985f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&jac)); 5992fa5cd67SKarl Rupp 6004b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 6012533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 6024b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 60369d2c0f9SBarry Smith pc->ops->reset = PCReset_Composite; 6044b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 6054b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 6064b9ad928SBarry Smith pc->ops->view = PCView_Composite; 6070a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 6084b9ad928SBarry Smith 6094b9ad928SBarry Smith pc->data = (void*)jac; 6104b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 6110a545947SLisandro Dalcin jac->work1 = NULL; 6120a545947SLisandro Dalcin jac->work2 = NULL; 6130a545947SLisandro Dalcin jac->head = NULL; 6144b9ad928SBarry Smith 6155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",PCCompositeSetType_Composite)); 6165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetType_C",PCCompositeGetType_Composite)); 6175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPCType_C",PCCompositeAddPCType_Composite)); 6185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",PCCompositeAddPC_Composite)); 6195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetNumberPC_C",PCCompositeGetNumberPC_Composite)); 6205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",PCCompositeGetPC_Composite)); 6215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",PCCompositeSpecialSetAlpha_Composite)); 6224b9ad928SBarry Smith PetscFunctionReturn(0); 6234b9ad928SBarry Smith } 624