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 3128b400f6SJacob 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) { 359566063dSJacob Faibussowitsch PetscCall(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 */ 419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(jac->work1,&jac->work2)); 424b9ad928SBarry Smith } 4349517cdeSBarry Smith if (pc->useAmat) mat = pc->mat; 449566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,x,y)); /* y <- B x */ 454b9ad928SBarry Smith while (next->next) { 464b9ad928SBarry Smith next = next->next; 479566063dSJacob Faibussowitsch PetscCall(MatMult(mat,y,jac->work1)); /* work1 <- A y */ 489566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->work2,-1.0,jac->work1,x)); /* work2 <- x - work1 */ 499566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,jac->work2,jac->work1)); /* work1 <- C work2 */ 509566063dSJacob Faibussowitsch PetscCall(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; 559566063dSJacob Faibussowitsch PetscCall(MatMult(mat,y,jac->work1)); 569566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 579566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,jac->work2,jac->work1)); 589566063dSJacob Faibussowitsch PetscCall(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; 7128b400f6SJacob 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 */ 739566063dSJacob Faibussowitsch PetscCall(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 } 809566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(next->pc,x,y)); 812533e041SBarry Smith while (next->previous) { 822533e041SBarry Smith next = next->previous; 839566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat,y,jac->work1)); 849566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 859566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(next->pc,jac->work2,jac->work1)); 869566063dSJacob Faibussowitsch PetscCall(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; 929566063dSJacob Faibussowitsch PetscCall(MatMultTranspose(mat,y,jac->work1)); 939566063dSJacob Faibussowitsch PetscCall(VecWAXPY(jac->work2,-1.0,jac->work1,x)); 949566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(next->pc,jac->work2,jac->work1)); 959566063dSJacob Faibussowitsch PetscCall(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; 11228b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 1137827d75bSBarry Smith PetscCheck(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 */ 1169566063dSJacob Faibussowitsch PetscCall(PCSetReusePreconditioner(next->pc,pc->reusepreconditioner)); 1179566063dSJacob Faibussowitsch PetscCall(PCSetReusePreconditioner(next->next->pc,pc->reusepreconditioner)); 118450d59ebSPatrick Farrell 1199566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,x,jac->work1)); 1209566063dSJacob Faibussowitsch PetscCall(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; 13028b400f6SJacob 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) { 1349566063dSJacob Faibussowitsch PetscCall(PCSetReusePreconditioner(next->pc,pc->reusepreconditioner)); 135450d59ebSPatrick Farrell next = next->next; 136450d59ebSPatrick Farrell } 137450d59ebSPatrick Farrell next = jac->head; 138450d59ebSPatrick Farrell 1399566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,x,y)); 1404b9ad928SBarry Smith while (next->next) { 1414b9ad928SBarry Smith next = next->next; 1429566063dSJacob Faibussowitsch PetscCall(PCApply(next->pc,x,jac->work1)); 1439566063dSJacob Faibussowitsch PetscCall(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; 15428b400f6SJacob Faibussowitsch PetscCheck(next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPCType() or -pc_composite_pcs"); 1559566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(next->pc,x,y)); 1562533e041SBarry Smith while (next->next) { 1572533e041SBarry Smith next = next->next; 1589566063dSJacob Faibussowitsch PetscCall(PCApplyTranspose(next->pc,x,jac->work1)); 1599566063dSJacob Faibussowitsch PetscCall(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) { 1729566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&jac->work1,NULL)); 1734b9ad928SBarry Smith } 1749566063dSJacob Faibussowitsch PetscCall(PCGetDM(pc,&dm)); 1754b9ad928SBarry Smith while (next) { 1762b1d202aSBarry Smith if (!next->pc->dm) { 1779566063dSJacob Faibussowitsch PetscCall(PCSetDM(next->pc,dm)); 1782b1d202aSBarry Smith } 1792b1d202aSBarry Smith if (!next->pc->mat) { 1809566063dSJacob Faibussowitsch PetscCall(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) { 1949566063dSJacob Faibussowitsch PetscCall(PCReset(next->pc)); 19569d2c0f9SBarry Smith next = next->next; 19669d2c0f9SBarry Smith } 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&jac->work1)); 1989566063dSJacob Faibussowitsch PetscCall(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; 2089566063dSJacob Faibussowitsch PetscCall(PCReset_Composite(pc)); 2094b9ad928SBarry Smith while (next) { 2109566063dSJacob Faibussowitsch PetscCall(PCDestroy(&next->pc)); 211724c2c99SHong Zhang next_tmp = next; 2124b9ad928SBarry Smith next = next->next; 2139566063dSJacob Faibussowitsch PetscCall(PetscFree(next_tmp)); 2144b9ad928SBarry Smith } 2152e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",NULL)); 2162e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetType_C",NULL)); 2172e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPCType_C",NULL)); 2182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",NULL)); 2192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetNumberPC_C",NULL)); 2202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",NULL)); 2212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 2234b9ad928SBarry Smith PetscFunctionReturn(0); 2244b9ad928SBarry Smith } 2254b9ad928SBarry Smith 2264416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PetscOptionItems *PetscOptionsObject,PC pc) 2274b9ad928SBarry Smith { 2284b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2299dcbbd2bSBarry Smith PetscInt nmax = 8,i; 2304b9ad928SBarry Smith PC_CompositeLink next; 231e5999256SBarry Smith char *pcs[8]; 232ace3abfcSBarry Smith PetscBool flg; 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith PetscFunctionBegin; 235d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Composite preconditioner options"); 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg)); 237*1baa6e33SBarry Smith if (flg) PetscCall(PCCompositeSetType(pc,jac->type)); 2389566063dSJacob Faibussowitsch PetscCall(PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPCType",pcs,&nmax,&flg)); 2394b9ad928SBarry Smith if (flg) { 2404b9ad928SBarry Smith for (i=0; i<nmax; i++) { 2419566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPCType(pc,pcs[i])); 2429566063dSJacob Faibussowitsch PetscCall(PetscFree(pcs[i])); /* deallocate string pcs[i], which is allocated in PetscOptionsStringArray() */ 2434b9ad928SBarry Smith } 2444b9ad928SBarry Smith } 245d0609cedSBarry Smith PetscOptionsHeadEnd(); 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith next = jac->head; 2484b9ad928SBarry Smith while (next) { 2499566063dSJacob Faibussowitsch PetscCall(PCSetFromOptions(next->pc)); 2504b9ad928SBarry Smith next = next->next; 2514b9ad928SBarry Smith } 2524b9ad928SBarry Smith PetscFunctionReturn(0); 2534b9ad928SBarry Smith } 2544b9ad928SBarry Smith 2556849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 2564b9ad928SBarry Smith { 2574b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2584b9ad928SBarry Smith PC_CompositeLink next = jac->head; 259ace3abfcSBarry Smith PetscBool iascii; 2604b9ad928SBarry Smith 2614b9ad928SBarry Smith PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 26332077d6dSBarry Smith if (iascii) { 2649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type])); 2659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n")); 2669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"---------------------------------\n")); 2674b9ad928SBarry Smith } 268*1baa6e33SBarry Smith if (iascii) PetscCall(PetscViewerASCIIPushTab(viewer)); 2694b9ad928SBarry Smith while (next) { 2709566063dSJacob Faibussowitsch PetscCall(PCView(next->pc,viewer)); 2714b9ad928SBarry Smith next = next->next; 2724b9ad928SBarry Smith } 27332077d6dSBarry Smith if (iascii) { 2749566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"---------------------------------\n")); 2764b9ad928SBarry Smith } 2774b9ad928SBarry Smith PetscFunctionReturn(0); 2784b9ad928SBarry Smith } 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2814b9ad928SBarry Smith 2821e6b0712SBarry Smith static PetscErrorCode PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2834b9ad928SBarry Smith { 2844b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2855fd66863SKarl Rupp 2864b9ad928SBarry Smith PetscFunctionBegin; 2874b9ad928SBarry Smith jac->alpha = alpha; 2884b9ad928SBarry Smith PetscFunctionReturn(0); 2894b9ad928SBarry Smith } 2904b9ad928SBarry Smith 2911e6b0712SBarry Smith static PetscErrorCode PCCompositeSetType_Composite(PC pc,PCCompositeType type) 2924b9ad928SBarry Smith { 293fad69fbaSJed Brown PC_Composite *jac = (PC_Composite*)pc->data; 294fad69fbaSJed Brown 2954b9ad928SBarry Smith PetscFunctionBegin; 2964b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 2974b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 2982533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 299421e10b8SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE || type == PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE) { 3004b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 3012533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Multiplicative; 3024b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 3034b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 3040298fd71SBarry Smith pc->ops->applytranspose = NULL; 305a7261c6bSprj- } else SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"Unknown composite preconditioner type"); 306fad69fbaSJed Brown jac->type = type; 3074b9ad928SBarry Smith PetscFunctionReturn(0); 3084b9ad928SBarry Smith } 3094b9ad928SBarry Smith 310c60c7ad4SBarry Smith static PetscErrorCode PCCompositeGetType_Composite(PC pc,PCCompositeType *type) 311c60c7ad4SBarry Smith { 312c60c7ad4SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 313c60c7ad4SBarry Smith 314c60c7ad4SBarry Smith PetscFunctionBegin; 315c60c7ad4SBarry Smith *type = jac->type; 316c60c7ad4SBarry Smith PetscFunctionReturn(0); 317c60c7ad4SBarry Smith } 318c60c7ad4SBarry Smith 3198aa07aa6SMatthew G. Knepley static PetscErrorCode PCCompositeAddPC_Composite(PC pc, PC subpc) 3204b9ad928SBarry Smith { 3214b9ad928SBarry Smith PC_Composite *jac; 3225a9f2f41SSatish Balay PC_CompositeLink next, ilink; 32379416396SBarry Smith PetscInt cnt = 0; 3242dcb1b2aSMatthew Knepley const char *prefix; 325d726e3a5SJed Brown char newprefix[20]; 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith PetscFunctionBegin; 3289566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &ilink)); 3290a545947SLisandro Dalcin ilink->next = NULL; 3308aa07aa6SMatthew G. Knepley ilink->pc = subpc; 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith jac = (PC_Composite *) pc->data; 3334b9ad928SBarry Smith next = jac->head; 3344b9ad928SBarry Smith if (!next) { 3355a9f2f41SSatish Balay jac->head = ilink; 3360298fd71SBarry Smith ilink->previous = NULL; 3374b9ad928SBarry Smith } else { 3384b9ad928SBarry Smith cnt++; 3394b9ad928SBarry Smith while (next->next) { 3404b9ad928SBarry Smith next = next->next; 3414b9ad928SBarry Smith cnt++; 3424b9ad928SBarry Smith } 3435a9f2f41SSatish Balay next->next = ilink; 344421e10b8SBarry Smith ilink->previous = next; 3454b9ad928SBarry Smith } 3469566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 3479566063dSJacob Faibussowitsch PetscCall(PCSetOptionsPrefix(subpc, prefix)); 3489566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(newprefix, 20, "sub_%d_", (int) cnt)); 3499566063dSJacob Faibussowitsch PetscCall(PCAppendOptionsPrefix(subpc, newprefix)); 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) subpc)); 3518aa07aa6SMatthew G. Knepley PetscFunctionReturn(0); 3528aa07aa6SMatthew G. Knepley } 3538aa07aa6SMatthew G. Knepley 3548aa07aa6SMatthew G. Knepley static PetscErrorCode PCCompositeAddPCType_Composite(PC pc, PCType type) 3558aa07aa6SMatthew G. Knepley { 3568aa07aa6SMatthew G. Knepley PC subpc; 3578aa07aa6SMatthew G. Knepley 3588aa07aa6SMatthew G. Knepley PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PCCreate(PetscObjectComm((PetscObject)pc), &subpc)); 3609566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject) subpc, (PetscObject) pc, 1)); 3619566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject) pc, (PetscObject) subpc)); 3629566063dSJacob Faibussowitsch PetscCall(PCCompositeAddPC_Composite(pc, subpc)); 3634b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 3649566063dSJacob Faibussowitsch PetscCall(PCSetType(subpc, type)); 3659566063dSJacob Faibussowitsch PetscCall(PCDestroy(&subpc)); 3664b9ad928SBarry Smith PetscFunctionReturn(0); 3674b9ad928SBarry Smith } 3684b9ad928SBarry Smith 3698e6eba06SBarry Smith static PetscErrorCode PCCompositeGetNumberPC_Composite(PC pc,PetscInt *n) 3708e6eba06SBarry Smith { 3718e6eba06SBarry Smith PC_Composite *jac; 3728e6eba06SBarry Smith PC_CompositeLink next; 3738e6eba06SBarry Smith 3748e6eba06SBarry Smith PetscFunctionBegin; 3758e6eba06SBarry Smith jac = (PC_Composite*)pc->data; 3768e6eba06SBarry Smith next = jac->head; 3778e6eba06SBarry Smith *n = 0; 3788e6eba06SBarry Smith while (next) { 3798e6eba06SBarry Smith next = next->next; 3808e6eba06SBarry Smith (*n) ++; 3818e6eba06SBarry Smith } 3828e6eba06SBarry Smith PetscFunctionReturn(0); 3838e6eba06SBarry Smith } 3848e6eba06SBarry Smith 3851e6b0712SBarry Smith static PetscErrorCode PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 3864b9ad928SBarry Smith { 3874b9ad928SBarry Smith PC_Composite *jac; 3884b9ad928SBarry Smith PC_CompositeLink next; 38979416396SBarry Smith PetscInt i; 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith PetscFunctionBegin; 3924b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3934b9ad928SBarry Smith next = jac->head; 3944b9ad928SBarry Smith for (i=0; i<n; i++) { 39528b400f6SJacob Faibussowitsch PetscCheck(next->next,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 3964b9ad928SBarry Smith next = next->next; 3974b9ad928SBarry Smith } 3984b9ad928SBarry Smith *subpc = next->pc; 3994b9ad928SBarry Smith PetscFunctionReturn(0); 4004b9ad928SBarry Smith } 4014b9ad928SBarry Smith 4024b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 403f39d8e23SSatish Balay /*@ 4044b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 4054b9ad928SBarry Smith 406ad4df100SBarry Smith Logically Collective on PC 4074b9ad928SBarry Smith 408c60c7ad4SBarry Smith Input Parameters: 4092a6744ebSBarry Smith + pc - the preconditioner context 4102a6744ebSBarry Smith - type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Options Database Key: 4134b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith Level: Developer 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith @*/ 4187087cfbeSBarry Smith PetscErrorCode PCCompositeSetType(PC pc,PCCompositeType type) 4194b9ad928SBarry Smith { 4204b9ad928SBarry Smith PetscFunctionBegin; 4210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 422c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 423cac4c232SBarry Smith PetscTryMethod(pc,"PCCompositeSetType_C",(PC,PCCompositeType),(pc,type)); 4244b9ad928SBarry Smith PetscFunctionReturn(0); 4254b9ad928SBarry Smith } 4264b9ad928SBarry Smith 427c60c7ad4SBarry Smith /*@ 428721f67b5SBarry Smith PCCompositeGetType - Gets the type of composite preconditioner. 429c60c7ad4SBarry Smith 430c60c7ad4SBarry Smith Logically Collective on PC 431c60c7ad4SBarry Smith 432c60c7ad4SBarry Smith Input Parameter: 433c60c7ad4SBarry Smith . pc - the preconditioner context 434c60c7ad4SBarry Smith 435c60c7ad4SBarry Smith Output Parameter: 436c60c7ad4SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 437c60c7ad4SBarry Smith 438c60c7ad4SBarry Smith Options Database Key: 439c60c7ad4SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 440c60c7ad4SBarry Smith 441c60c7ad4SBarry Smith Level: Developer 442c60c7ad4SBarry Smith 443c60c7ad4SBarry Smith @*/ 444c60c7ad4SBarry Smith PetscErrorCode PCCompositeGetType(PC pc,PCCompositeType *type) 445c60c7ad4SBarry Smith { 446c60c7ad4SBarry Smith PetscFunctionBegin; 447c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 448cac4c232SBarry Smith PetscUseMethod(pc,"PCCompositeGetType_C",(PC,PCCompositeType*),(pc,type)); 449c60c7ad4SBarry Smith PetscFunctionReturn(0); 450c60c7ad4SBarry Smith } 451c60c7ad4SBarry Smith 452f39d8e23SSatish Balay /*@ 4534b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 4544b9ad928SBarry Smith for alphaI + R + S 4554b9ad928SBarry Smith 456ad4df100SBarry Smith Logically Collective on PC 4574b9ad928SBarry Smith 458d8d19677SJose E. Roman Input Parameters: 4594b9ad928SBarry Smith + pc - the preconditioner context 4604b9ad928SBarry Smith - alpha - scale on identity 4614b9ad928SBarry Smith 4624b9ad928SBarry Smith Level: Developer 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith @*/ 4657087cfbeSBarry Smith PetscErrorCode PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 4664b9ad928SBarry Smith { 4674b9ad928SBarry Smith PetscFunctionBegin; 4680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 469c5eb9154SBarry Smith PetscValidLogicalCollectiveScalar(pc,alpha,2); 470cac4c232SBarry Smith PetscTryMethod(pc,"PCCompositeSpecialSetAlpha_C",(PC,PetscScalar),(pc,alpha)); 4714b9ad928SBarry Smith PetscFunctionReturn(0); 4724b9ad928SBarry Smith } 4734b9ad928SBarry Smith 4744b9ad928SBarry Smith /*@C 4758aa07aa6SMatthew G. Knepley PCCompositeAddPCType - Adds another PC of the given type to the composite PC. 4764b9ad928SBarry Smith 4774b9ad928SBarry Smith Collective on PC 4784b9ad928SBarry Smith 4794b9ad928SBarry Smith Input Parameters: 4802a6744ebSBarry Smith + pc - the preconditioner context 4812a6744ebSBarry Smith - type - the type of the new preconditioner 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith Level: Developer 4844b9ad928SBarry Smith 485db781477SPatrick Sanan .seealso: `PCCompositeAddPC()` 4864b9ad928SBarry Smith @*/ 4878aa07aa6SMatthew G. Knepley PetscErrorCode PCCompositeAddPCType(PC pc,PCType type) 4884b9ad928SBarry Smith { 4894b9ad928SBarry Smith PetscFunctionBegin; 4900700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 491cac4c232SBarry Smith PetscTryMethod(pc,"PCCompositeAddPCType_C",(PC,PCType),(pc,type)); 4928aa07aa6SMatthew G. Knepley PetscFunctionReturn(0); 4938aa07aa6SMatthew G. Knepley } 4948aa07aa6SMatthew G. Knepley 4958aa07aa6SMatthew G. Knepley /*@ 4968aa07aa6SMatthew G. Knepley PCCompositeAddPC - Adds another PC to the composite PC. 4978aa07aa6SMatthew G. Knepley 4988aa07aa6SMatthew G. Knepley Collective on PC 4998aa07aa6SMatthew G. Knepley 5008aa07aa6SMatthew G. Knepley Input Parameters: 5018aa07aa6SMatthew G. Knepley + pc - the preconditioner context 5028aa07aa6SMatthew G. Knepley - subpc - the new preconditioner 5038aa07aa6SMatthew G. Knepley 5048aa07aa6SMatthew G. Knepley Level: Developer 5058aa07aa6SMatthew G. Knepley 506db781477SPatrick Sanan .seealso: `PCCompositeAddPCType()` 5078aa07aa6SMatthew G. Knepley @*/ 5088aa07aa6SMatthew G. Knepley PetscErrorCode PCCompositeAddPC(PC pc, PC subpc) 5098aa07aa6SMatthew G. Knepley { 5108aa07aa6SMatthew G. Knepley PetscFunctionBegin; 5118aa07aa6SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5128aa07aa6SMatthew G. Knepley PetscValidHeaderSpecific(subpc,PC_CLASSID,2); 513cac4c232SBarry Smith PetscTryMethod(pc,"PCCompositeAddPC_C",(PC,PC),(pc,subpc)); 5144b9ad928SBarry Smith PetscFunctionReturn(0); 5154b9ad928SBarry Smith } 5164b9ad928SBarry Smith 5178e6eba06SBarry Smith /*@ 5188e6eba06SBarry Smith PCCompositeGetNumberPC - Gets the number of PC objects in the composite PC. 5198e6eba06SBarry Smith 5208e6eba06SBarry Smith Not Collective 5218e6eba06SBarry Smith 5228e6eba06SBarry Smith Input Parameter: 5238e6eba06SBarry Smith . pc - the preconditioner context 5248e6eba06SBarry Smith 5258e6eba06SBarry Smith Output Parameter: 5268e6eba06SBarry Smith . num - the number of sub pcs 5278e6eba06SBarry Smith 5288e6eba06SBarry Smith Level: Developer 5298e6eba06SBarry Smith 530db781477SPatrick Sanan .seealso: `PCCompositeGetPC()` 5318e6eba06SBarry Smith @*/ 5328e6eba06SBarry Smith PetscErrorCode PCCompositeGetNumberPC(PC pc,PetscInt *num) 5338e6eba06SBarry Smith { 5348e6eba06SBarry Smith PetscFunctionBegin; 5358e6eba06SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5368e6eba06SBarry Smith PetscValidIntPointer(num,2); 537cac4c232SBarry Smith PetscUseMethod(pc,"PCCompositeGetNumberPC_C",(PC,PetscInt*),(pc,num)); 5388e6eba06SBarry Smith PetscFunctionReturn(0); 5398e6eba06SBarry Smith } 5408e6eba06SBarry Smith 541f39d8e23SSatish Balay /*@ 5424b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 5434b9ad928SBarry Smith 5444b9ad928SBarry Smith Not Collective 5454b9ad928SBarry Smith 546d8d19677SJose E. Roman Input Parameters: 5472a6744ebSBarry Smith + pc - the preconditioner context 5482a6744ebSBarry Smith - n - the number of the pc requested 5494b9ad928SBarry Smith 5504b9ad928SBarry Smith Output Parameters: 5514b9ad928SBarry Smith . subpc - the PC requested 5524b9ad928SBarry Smith 5534b9ad928SBarry Smith Level: Developer 5544b9ad928SBarry Smith 55595452b02SPatrick Sanan Notes: 55695452b02SPatrick Sanan To use a different operator to construct one of the inner preconditioners first call PCCompositeGetPC(), then 5572b1d202aSBarry Smith call PCSetOperators() on that PC. 5582b1d202aSBarry Smith 559db781477SPatrick Sanan .seealso: `PCCompositeAddPCType()`, `PCCompositeGetNumberPC()`, `PCSetOperators()` 5604b9ad928SBarry Smith @*/ 5617087cfbeSBarry Smith PetscErrorCode PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 5624b9ad928SBarry Smith { 5634b9ad928SBarry Smith PetscFunctionBegin; 5640700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 5654482741eSBarry Smith PetscValidPointer(subpc,3); 566cac4c232SBarry Smith PetscUseMethod(pc,"PCCompositeGetPC_C",(PC,PetscInt,PC*),(pc,n,subpc)); 5674b9ad928SBarry Smith PetscFunctionReturn(0); 5684b9ad928SBarry Smith } 5694b9ad928SBarry Smith 5704b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 5714b9ad928SBarry Smith 5724b9ad928SBarry Smith /*MC 5734b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 5744b9ad928SBarry Smith 5754b9ad928SBarry Smith Options Database Keys: 5762eab2d5bSJungho Lee + -pc_composite_type <type: one of multiplicative, additive, symmetric_multiplicative, special> - Sets composite preconditioner type 5772b1d202aSBarry Smith . -pc_use_amat - activates PCSetUseAmat() 57851f519a2SBarry Smith - -pc_composite_pcs - <pc0,pc1,...> list of PCs to compose 5794b9ad928SBarry Smith 5804b9ad928SBarry Smith Level: intermediate 5814b9ad928SBarry Smith 58295452b02SPatrick Sanan Notes: 58395452b02SPatrick Sanan To use a Krylov method inside the composite preconditioner, set the PCType of one or more 5844b9ad928SBarry Smith inner PCs to be PCKSP. 5854b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 586b3ef52cdSBarry Smith the incorrect answer) unless you use KSPFGMRES as the outer Krylov method 5872b1d202aSBarry Smith To use a different operator to construct one of the inner preconditioners first call PCCompositeGetPC(), then 5882b1d202aSBarry Smith call PCSetOperators() on that PC. 5894b9ad928SBarry Smith 590db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 591db781477SPatrick Sanan `PCSHELL`, `PCKSP`, `PCCompositeSetType()`, `PCCompositeSpecialSetAlpha()`, `PCCompositeAddPCType()`, 592db781477SPatrick Sanan `PCCompositeGetPC()`, `PCSetUseAmat()` 5934b9ad928SBarry Smith 5944b9ad928SBarry Smith M*/ 5954b9ad928SBarry Smith 5968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Composite(PC pc) 5974b9ad928SBarry Smith { 5984b9ad928SBarry Smith PC_Composite *jac; 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith PetscFunctionBegin; 6019566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&jac)); 6022fa5cd67SKarl Rupp 6034b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 6042533e041SBarry Smith pc->ops->applytranspose = PCApplyTranspose_Composite_Additive; 6054b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 60669d2c0f9SBarry Smith pc->ops->reset = PCReset_Composite; 6074b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 6084b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 6094b9ad928SBarry Smith pc->ops->view = PCView_Composite; 6100a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 6114b9ad928SBarry Smith 6124b9ad928SBarry Smith pc->data = (void*)jac; 6134b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 6140a545947SLisandro Dalcin jac->work1 = NULL; 6150a545947SLisandro Dalcin jac->work2 = NULL; 6160a545947SLisandro Dalcin jac->head = NULL; 6174b9ad928SBarry Smith 6189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSetType_C",PCCompositeSetType_Composite)); 6199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetType_C",PCCompositeGetType_Composite)); 6209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPCType_C",PCCompositeAddPCType_Composite)); 6219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeAddPC_C",PCCompositeAddPC_Composite)); 6229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetNumberPC_C",PCCompositeGetNumberPC_Composite)); 6239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeGetPC_C",PCCompositeGetPC_Composite)); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",PCCompositeSpecialSetAlpha_Composite)); 6254b9ad928SBarry Smith PetscFunctionReturn(0); 6264b9ad928SBarry Smith } 627